CGEM BOSS 6.6.5.h
BESIII Offline Software System
Loading...
Searching...
No Matches
incllambda.cxx
Go to the documentation of this file.
1#include "GaudiKernel/MsgStream.h"
2#include "GaudiKernel/AlgFactory.h"
3#include "GaudiKernel/ISvcLocator.h"
4#include "GaudiKernel/SmartDataPtr.h"
5#include "GaudiKernel/IDataProviderSvc.h"
6#include "GaudiKernel/PropertyMgr.h"
7
9#include "EventModel/Event.h"
10
15
16#include "TMath.h"
17#include "GaudiKernel/INTupleSvc.h"
18#include "GaudiKernel/NTuple.h"
19#include "GaudiKernel/Bootstrap.h"
20//#include "GaudiKernel/IHistogramSvc.h"
21#include "GaudiKernel/ITHistSvc.h"
22
23#include "CLHEP/Vector/ThreeVector.h"
24#include "CLHEP/Vector/LorentzVector.h"
25#include "CLHEP/Vector/TwoVector.h"
26using CLHEP::Hep3Vector;
27using CLHEP::Hep2Vector;
28using CLHEP::HepLorentzVector;
29#include "CLHEP/Geometry/Point3D.h"
30#ifndef ENABLE_BACKWARDS_COMPATIBILITY
32#endif
33
35#include "VertexFit/VertexFit.h"
38
40
41#include <vector>
42
43const double ecms = 3.097;
44const double mpi = 0.13957;
45const double mp = 0.93827203;
46const double mlam = 1.115683;
47const double xmass[5] = {0.000511, 0.105658, 0.139570,0.493677, 0.938272};
48// e mu pi K p
49
50typedef std::vector<int> Vint;
51typedef std::vector<HepLorentzVector> Vp4;
52
53/////////////////////////////////////////////////////////////////////////////
54// //
55// Lambda + X //
56// |-->p+- pi-+ //
57// //
58// xugf 2009.06 //
59/////////////////////////////////////////////////////////////////////////////
60
61incllambda::incllambda(const std::string& name, ISvcLocator* pSvcLocator) :
62 Algorithm(name, pSvcLocator) {
63
64 m_tuple1 = 0;
65 for(int i = 0; i < 10; i++) m_pass[i] = 0;
66
67//Declare the properties
68 declareProperty("Vr0cut", m_vr0cut=20.0);
69 declareProperty("Vz0cut", m_vz0cut=50.0);
70 declareProperty("CheckDedx", m_checkDedx = 1);
71 declareProperty("CheckTof", m_checkTof = 1);
72
73}
74
75// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
77 MsgStream log(msgSvc(), name());
78
79 log << MSG::INFO << "in initialize()" << endmsg;
80
81 StatusCode status;
82
83 if(service("THistSvc", m_thsvc).isFailure()) {
84 log << MSG::ERROR << "Couldn't get THistSvc" << endreq;
85 return StatusCode::FAILURE;
86 }
87 // "DQAHist" is fixed
88 TH1F* incllambda_mass = new TH1F("InclLam_mass","INCLUSIVE_LAMBDA_MASS",150,1.11,1.125);
89 incllambda_mass->GetXaxis()->SetTitle("M_{p#pi} (GeV)");
90 incllambda_mass->GetYaxis()->SetTitle("Nentries/0.1MeV/c^{2}");
91 if(m_thsvc->regHist("/DQAHist/InclLam/InclLam_mass", incllambda_mass).isFailure()) {
92 log << MSG::ERROR << "Couldn't register incllambda_mass" << endreq;
93 }
94
95//*****************************************
96
97 NTuplePtr nt1(ntupleSvc(), "DQAFILE/InclLam");
98 if ( nt1 ) m_tuple1 = nt1;
99 else {
100 m_tuple1 = ntupleSvc()->book ("DQAFILE/InclLam", CLID_ColumnWiseTuple, "incllam Ntuple");
101 if ( m_tuple1 ) {
102 status = m_tuple1->addItem ("npi", m_npi);
103 status = m_tuple1->addItem ("np", m_np);
104 status = m_tuple1->addItem ("lxyz", m_lxyz);
105 status = m_tuple1->addItem ("ctau", m_ctau);
106 status = m_tuple1->addItem ("elxyz", m_elxyz);
107 status = m_tuple1->addItem ("lamchi", m_chis);
108 status = m_tuple1->addItem ("mlamraw", m_lamRawMass);
109 status = m_tuple1->addItem ("plam", m_plam);
110 }
111 else {
112 log << MSG::ERROR << " Cannot book N-tuple:" << long(m_tuple1) << endmsg;
113 return StatusCode::FAILURE;
114 }
115 }
116 //
117 //--------end of book--------
118 //
119
120 log << MSG::INFO << "successfully return from initialize()" <<endmsg;
121 return StatusCode::SUCCESS;
122
123}
124
125// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
127 StatusCode sc = StatusCode::SUCCESS;
128
129 MsgStream log(msgSvc(), name());
130 log << MSG::INFO << "in execute()" << endreq;
131
132 // DQA
133 // Add the line below at the beginning of execute()
134 //
135 setFilterPassed(false);
136
137 m_pass[0] += 1;
138
139 SmartDataPtr<Event::EventHeader> eventHeader(eventSvc(),"/Event/EventHeader");
140 SmartDataPtr<EvtRecEvent> evtRecEvent(eventSvc(), EventModel::EvtRec::EvtRecEvent);
141 SmartDataPtr<EvtRecTrackCol> evtRecTrkCol(eventSvc(), EventModel::EvtRec::EvtRecTrackCol);
142
143 Vint ilam, ip, ipi, iGood;
144 iGood.clear();
145 ilam.clear();
146 ip.clear();
147 ipi.clear();
148
149 Vp4 ppi, pp;
150 ppi.clear();
151 pp.clear();
152
153 int TotCharge = 0;
154 for(int i = 0; i < evtRecEvent->totalCharged(); i++){
155 EvtRecTrackIterator itTrk=evtRecTrkCol->begin() + i;
156 if(!(*itTrk)->isMdcTrackValid()) continue;
157 RecMdcTrack *mdcTrk = (*itTrk)->mdcTrack();
158 if(fabs(mdcTrk->z()) >= m_vz0cut) continue;
159 if(mdcTrk->r() >= m_vr0cut) continue;
160 iGood.push_back(i);
161 TotCharge += mdcTrk->charge();
162 }
163 //
164 // Finish Good Charged Track Selection
165 //
166 int nGood = iGood.size();
167
168 //
169 // Charge track number cut
170 //
171
172 if((nGood < 2) || (TotCharge!=0)) return sc;
173
174 m_pass[1] += 1;
175 //
176 // Assign 4-momentum to each charged track
177 //
179 for(int i = 0; i < nGood; i++) {
180 EvtRecTrackIterator itTrk = evtRecTrkCol->begin() + iGood[i];
181 // if(pid) delete pid;
182 pid->init();
183 pid->setMethod(pid->methodProbability());
184 pid->setChiMinCut(4);
185 pid->setRecTrack(*itTrk);
186 pid->usePidSys(pid->useDedx() | pid->useTof1() | pid->useTof2()); // use PID sub-system
187 pid->identify(pid->onlyPion() | pid->onlyProton()); // seperater Pion/Kaon
188 pid->calculate();
189 if(!(pid->IsPidInfoValid())) continue;
190// RecMdcTrack* mdcTrk = (*itTrk)->mdcTrack();
191 if(!(*itTrk)->isMdcKalTrackValid()) continue ;
192 RecMdcKalTrack* mdcKalTrk = (*itTrk)->mdcKalTrack();
193
194 if(mdcKalTrk->charge() == 0) continue;
195
196 if(pid->probPion() > 0.001 && (pid->probPion() > pid->probProton())){
198 ipi.push_back(iGood[i]);
199 HepLorentzVector ptrk;
200 ptrk.setPx(mdcKalTrk->px());
201 ptrk.setPy(mdcKalTrk->py());
202 ptrk.setPz(mdcKalTrk->pz());
203 double p3 = ptrk.mag();
204 ptrk.setE(sqrt(p3*p3+mpi*mpi));
205 ppi.push_back(ptrk);
206 }
207 if(pid->probProton() > 0.001 && (pid->probProton() > pid->probPion())){
209 ip.push_back(iGood[i]);
210 HepLorentzVector ptrk;
211 ptrk.setPx(mdcKalTrk->px());
212 ptrk.setPy(mdcKalTrk->py());
213 ptrk.setPz(mdcKalTrk->pz());
214 double p3 = ptrk.mag();
215 ptrk.setE(sqrt(p3*p3+mp*mp));
216 pp.push_back(ptrk);
217 }
218 }
219
220 m_pass[2] += 1;
221
222 int npi = ipi.size();
223 int np = ip.size();
224 m_npi = npi;
225 m_np = np;
226 if(npi < 1 || np <1) return sc;
227
228 m_pass[3] += 1;
229//
230//****** Second Vertex Check************
231//
232 double chi=999.;
233 double delm;
234
235 HepPoint3D vx(0., 0., 0.);
236 HepSymMatrix Evx(3, 0);
237 double bx = 1E+6;
238 double by = 1E+6;
239 double bz = 1E+6;
240 Evx[0][0] = bx*bx;
241 Evx[1][1] = by*by;
242 Evx[2][2] = bz*bz;
243 VertexParameter vxpar;
244 vxpar.setVx(vx);
245 vxpar.setEvx(Evx);
246
247 int ipion=-1, iproton=-1;
248
249 VertexFit *vtxfit0 = VertexFit::instance();
251
252 for(unsigned int i1 = 0; i1 < ipi.size(); i1++) {
253 RecMdcKalTrack *piTrk = (*(evtRecTrkCol->begin()+ipi[i1]))->mdcKalTrack();
255 WTrackParameter wpitrk(mpi, piTrk->helix(), piTrk->err());
256
257 for(unsigned int i2 = 0; i2 < ip.size(); i2++) {
258 RecMdcKalTrack *protonTrk = (*(evtRecTrkCol->begin()+ip[i2]))->mdcKalTrack();
260 WTrackParameter wptrk(mp, protonTrk->helix(), protonTrk->err());
261
262 int NCharge = 0;
263 NCharge =piTrk->charge() + protonTrk->charge();
264
265// cout << "TotNcharge = " << NCharge << endl;
266
267 if(NCharge != 0) continue;
268
269 vtxfit0->init();
270 vtxfit0->AddTrack(0, wpitrk);
271 vtxfit0->AddTrack(1, wptrk);
272 vtxfit0->AddVertex(0, vxpar, 0, 1);
273 if(!(vtxfit0->Fit(0))) continue;
274 vtxfit0->BuildVirtualParticle(0);
275
276 vtxfit->init();
277// vtxfit->setPrimaryVertex(bs);
278 vtxfit->AddTrack(0, vtxfit0->wVirtualTrack(0));
279 vtxfit->setVpar(vtxfit0->vpar(0));
280 if(!vtxfit->Fit()) continue;
281
282 //double mass = vtxfit->p4par().m();
283 //if(mass < m_massRangeLower) continue;
284 //if(mass > m_massRangeHigher) continue;
285 if(vtxfit->chisq() > chi) continue;
286// if(fabs((vtxfit->p4par()).m()-mlam) > delm) continue;
287 delm = fabs((vtxfit->p4par()).m()-mlam);
288 chi = vtxfit->chisq();
289 ipion=ipi[i1];
290 iproton=ip[i2];
291 }
292 }
293//
294// Lammda check
295//
296 if(ipion==-1 || iproton==-1) return sc;
297 m_pass[4] += 1;
298
299 RecMdcKalTrack *pionTrk = (*(evtRecTrkCol->begin()+ipion))->mdcKalTrack();
301 WTrackParameter wpiontrk(mpi, pionTrk->helix(), pionTrk->err());
302
303 RecMdcKalTrack *protonTrk = (*(evtRecTrkCol->begin()+iproton))->mdcKalTrack();
305 WTrackParameter wprotontrk(mp, protonTrk->helix(), protonTrk->err());
306 vtxfit0->init();
307 vtxfit0->AddTrack(0, wpiontrk);
308 vtxfit0->AddTrack(1, wprotontrk);
309 vtxfit0->AddVertex(0, vxpar, 0, 1);
310 if(!(vtxfit0->Fit(0))) return sc;
311 vtxfit0->BuildVirtualParticle(0);
312
313 vtxfit->init();
314// vtxfit->setPrimaryVertex(bs);
315 vtxfit->AddTrack(0, vtxfit0->wVirtualTrack(0));
316 vtxfit->setVpar(vtxfit0->vpar(0));
317 if(!vtxfit->Fit()) return sc;
318
319 m_lamRawMass = vtxfit->p4par().m();
320 m_ctau = vtxfit->ctau();
321 m_lxyz = vtxfit->decayLength();
322 m_elxyz = vtxfit->decayLengthError();
323 m_chis = vtxfit->chisq();
324 m_plam = vtxfit->p4par().rho();
325
326 if(m_lxyz>3.5){
327
328 TH1 *h(0);
329 if (m_thsvc->getHist("/DQAHist/InclLam/InclLam_mass", h).isSuccess()) {
330 h->Fill(m_lamRawMass);
331 } else {
332 log << MSG::ERROR << "Couldn't retrieve incllam_mass" << endreq;
333 }
334
335 }
336 m_tuple1->write();
337////////////////////////////////////
338//DQA
339// set tag and quality
340 // Pid: 1 - electron, 2 - muon, 3 - pion, 4 - kaon, 5 - proton
341// (*(evtRecTrkCol->begin()+ipion))->setPartId(3);
342// (*(evtRecTrkCol->begin()+iproton))->setPartId(5);
343 (*(evtRecTrkCol->begin()+ipion))->tagPion();
344 (*(evtRecTrkCol->begin()+iproton))->tagProton();
345 // Quality: defined by whether dE/dx or TOF is used to identify particle
346 // 0 - no dE/dx, no TOF (can be used for dE/dx and TOF calibration)
347 // 1 - only dE/dx (can be used for TOF calibration)
348 // 2 - only TOF (can be used for dE/dx calibration)
349 // 3 - Both dE/dx and TOF
350 (*(evtRecTrkCol->begin()+ipion))->setQuality(1);
351 (*(evtRecTrkCol->begin()+iproton))->setQuality(1);
352
353//--------------------------------------------------
354 // Add the line below at the end of execute(), (before return)
355 //
356 setFilterPassed(true);
357
358 return StatusCode::SUCCESS;
359}
360
361// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
363
364 MsgStream log(msgSvc(), name());
365 log << MSG::INFO << "in finalize()" << endmsg;
366 log << MSG::INFO << "Total Entries : " << m_pass[0] << endreq;
367 log << MSG::INFO << "Qtot Cut : " << m_pass[1] << endreq;
368 log << MSG::INFO << "Nks : " << m_pass[2] << endreq;
369 log << MSG::INFO << "Good Ks cut : " << m_pass[3] << endreq;
370 return StatusCode::SUCCESS;
371}
372
EvtRecTrackCol::iterator EvtRecTrackIterator
std::vector< HepLorentzVector > Vp4
Definition Gam4pikp.cxx:53
const double mpi
Definition Gam4pikp.cxx:47
std::vector< int > Vint
Definition Gam4pikp.cxx:52
INTupleSvc * ntupleSvc()
IMessageSvc * msgSvc()
const HepVector & helix() const
const double px() const
static void setPidType(PidType pidType)
const double pz() const
const double py() const
const HepSymMatrix & err() const
const int charge() const
const double r() const
Definition DstMdcTrack.h:64
const int charge() const
Definition DstMdcTrack.h:53
const double z() const
Definition DstMdcTrack.h:63
int useTof2() const
int onlyProton() const
int methodProbability() const
int useDedx() const
int onlyPion() const
int useTof1() const
void setChiMinCut(const double chi=4)
void setRecTrack(EvtRecTrack *trk)
void setMethod(const int method)
Definition ParticleID.h:101
void identify(const int pidcase)
Definition ParticleID.h:110
void usePidSys(const int pidsys)
Definition ParticleID.h:104
static ParticleID * instance()
bool IsPidInfoValid() const
double probPion() const
Definition ParticleID.h:130
void calculate()
void init()
double probProton() const
Definition ParticleID.h:132
HepLorentzVector p4par() const
double ctau() const
double decayLength() const
double decayLengthError() const
static SecondVertexFit * instance()
void setVpar(const VertexParameter vpar)
double chisq() const
void AddTrack(const int number, const double mass, const RecMdcTrack *trk)
Definition TrackPool.cxx:22
WTrackParameter wVirtualTrack(int n) const
Definition VertexFit.h:91
void init()
Definition VertexFit.cxx:29
void AddVertex(int number, VertexParameter vpar, std::vector< int > lis)
Definition VertexFit.cxx:89
static VertexFit * instance()
Definition VertexFit.cxx:15
VertexParameter vpar(int n) const
Definition VertexFit.h:88
void BuildVirtualParticle(int number)
bool Fit()
void setEvx(const HepSymMatrix &eVx)
void setVx(const HepPoint3D &vx)
StatusCode initialize()
StatusCode finalize()
incllambda(const std::string &name, ISvcLocator *pSvcLocator)
StatusCode execute()
HepGeom::Point3D< double > HepPoint3D
std::vector< HepLorentzVector > Vp4
const double mp
const double xmass[5]
const double mpi
const double mlam
const double ecms
std::vector< int > Vint
_EXTERN_ std::string EvtRecEvent
Definition EventModel.h:134
_EXTERN_ std::string EvtRecTrackCol
Definition EventModel.h:135