BOSS 8.0.0
BESIII Offline Software System
Loading...
Searching...
No Matches
incllambda Class Reference

#include <incllambda.h>

Inheritance diagram for incllambda:

Public Member Functions

 incllambda (const std::string &name, ISvcLocator *pSvcLocator)
StatusCode initialize ()
StatusCode execute ()
StatusCode finalize ()

Detailed Description

Definition at line 15 of file incllambda.h.

Constructor & Destructor Documentation

◆ incllambda()

incllambda::incllambda ( const std::string & name,
ISvcLocator * pSvcLocator )

Definition at line 61 of file incllambda.cxx.

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}

Referenced by incllambda().

Member Function Documentation

◆ execute()

StatusCode incllambda::execute ( )

Definition at line 128 of file incllambda.cxx.

128 {
129 StatusCode sc = StatusCode::SUCCESS;
130
131 MsgStream log( msgSvc(), name() );
132 log << MSG::INFO << "in execute()" << endmsg;
133
134 // DQA
135 // Add the line below at the beginning of execute()
136 //
137 setFilterPassed( false );
138
139 m_pass[0] += 1;
140
141 SmartDataPtr<Event::EventHeader> eventHeader( eventSvc(), "/Event/EventHeader" );
142 SmartDataPtr<EvtRecEvent> evtRecEvent( eventSvc(), EventModel::EvtRec::EvtRecEvent );
143 SmartDataPtr<EvtRecTrackCol> evtRecTrkCol( eventSvc(), EventModel::EvtRec::EvtRecTrackCol );
144
145 Vint ilam, ip, ipi, iGood;
146 iGood.clear();
147 ilam.clear();
148 ip.clear();
149 ipi.clear();
150
151 Vp4 ppi, pp;
152 ppi.clear();
153 pp.clear();
154
155 int TotCharge = 0;
156 for ( int i = 0; i < evtRecEvent->totalCharged(); i++ )
157 {
158 EvtRecTrackIterator itTrk = evtRecTrkCol->begin() + i;
159 if ( !( *itTrk )->isMdcTrackValid() ) continue;
160 RecMdcTrack* mdcTrk = ( *itTrk )->mdcTrack();
161 if ( fabs( mdcTrk->z() ) >= m_vz0cut ) continue;
162 if ( mdcTrk->r() >= m_vr0cut ) continue;
163 iGood.push_back( i );
164 TotCharge += mdcTrk->charge();
165 }
166 //
167 // Finish Good Charged Track Selection
168 //
169 int nGood = iGood.size();
170
171 //
172 // Charge track number cut
173 //
174
175 if ( ( nGood < 2 ) || ( TotCharge != 0 ) ) return sc;
176
177 m_pass[1] += 1;
178 //
179 // Assign 4-momentum to each charged track
180 //
181 ParticleID* pid = ParticleID::instance();
182 for ( int i = 0; i < nGood; i++ )
183 {
184 EvtRecTrackIterator itTrk = evtRecTrkCol->begin() + iGood[i];
185 // if(pid) delete pid;
186 pid->init();
187 pid->setMethod( pid->methodProbability() );
188 pid->setChiMinCut( 4 );
189 pid->setRecTrack( *itTrk );
190 pid->usePidSys( pid->useDedx() | pid->useTof1() | pid->useTof2() ); // use PID sub-system
191 pid->identify( pid->onlyPion() | pid->onlyProton() ); // seperater Pion/Kaon
192 pid->calculate();
193 if ( !( pid->IsPidInfoValid() ) ) continue;
194 // RecMdcTrack* mdcTrk = (*itTrk)->mdcTrack();
195 if ( !( *itTrk )->isMdcKalTrackValid() ) continue;
196 RecMdcKalTrack* mdcKalTrk = ( *itTrk )->mdcKalTrack();
197
198 if ( mdcKalTrk->charge() == 0 ) continue;
199
200 if ( pid->probPion() > 0.001 && ( pid->probPion() > pid->probProton() ) )
201 {
202 mdcKalTrk->setPidType( RecMdcKalTrack::pion );
203 ipi.push_back( iGood[i] );
204 HepLorentzVector ptrk;
205 ptrk.setPx( mdcKalTrk->px() );
206 ptrk.setPy( mdcKalTrk->py() );
207 ptrk.setPz( mdcKalTrk->pz() );
208 double p3 = ptrk.mag();
209 ptrk.setE( sqrt( p3 * p3 + mpi * mpi ) );
210 ppi.push_back( ptrk );
211 }
212 if ( pid->probProton() > 0.001 && ( pid->probProton() > pid->probPion() ) )
213 {
215 ip.push_back( iGood[i] );
216 HepLorentzVector ptrk;
217 ptrk.setPx( mdcKalTrk->px() );
218 ptrk.setPy( mdcKalTrk->py() );
219 ptrk.setPz( mdcKalTrk->pz() );
220 double p3 = ptrk.mag();
221 ptrk.setE( sqrt( p3 * p3 + mp * mp ) );
222 pp.push_back( ptrk );
223 }
224 }
225
226 m_pass[2] += 1;
227
228 int npi = ipi.size();
229 int np = ip.size();
230 m_npi = npi;
231 m_np = np;
232 if ( npi < 1 || np < 1 ) return sc;
233
234 m_pass[3] += 1;
235 //
236 //****** Second Vertex Check************
237 //
238 double chi = 999.;
239 double delm;
240
241 HepPoint3D vx( 0., 0., 0. );
242 HepSymMatrix Evx( 3, 0 );
243 double bx = 1E+6;
244 double by = 1E+6;
245 double bz = 1E+6;
246 Evx[0][0] = bx * bx;
247 Evx[1][1] = by * by;
248 Evx[2][2] = bz * bz;
249 VertexParameter vxpar;
250 vxpar.setVx( vx );
251 vxpar.setEvx( Evx );
252
253 int ipion = -1, iproton = -1;
254
255 VertexFit* vtxfit0 = VertexFit::instance();
256 SecondVertexFit* vtxfit = SecondVertexFit::instance();
257
258 for ( unsigned int i1 = 0; i1 < ipi.size(); i1++ )
259 {
260 RecMdcKalTrack* piTrk = ( *( evtRecTrkCol->begin() + ipi[i1] ) )->mdcKalTrack();
262 WTrackParameter wpitrk( mpi, piTrk->helix(), piTrk->err() );
263
264 for ( unsigned int i2 = 0; i2 < ip.size(); i2++ )
265 {
266 RecMdcKalTrack* protonTrk = ( *( evtRecTrkCol->begin() + ip[i2] ) )->mdcKalTrack();
268 WTrackParameter wptrk( mp, protonTrk->helix(), protonTrk->err() );
269
270 int NCharge = 0;
271 NCharge = piTrk->charge() + protonTrk->charge();
272
273 // cout << "TotNcharge = " << NCharge << endl;
274
275 if ( NCharge != 0 ) continue;
276
277 vtxfit0->init();
278 vtxfit0->AddTrack( 0, wpitrk );
279 vtxfit0->AddTrack( 1, wptrk );
280 vtxfit0->AddVertex( 0, vxpar, 0, 1 );
281 if ( !( vtxfit0->Fit( 0 ) ) ) continue;
282 vtxfit0->BuildVirtualParticle( 0 );
283
284 vtxfit->init();
285 // vtxfit->setPrimaryVertex(bs);
286 vtxfit->AddTrack( 0, vtxfit0->wVirtualTrack( 0 ) );
287 vtxfit->setVpar( vtxfit0->vpar( 0 ) );
288 if ( !vtxfit->Fit() ) continue;
289
290 // double mass = vtxfit->p4par().m();
291 // if(mass < m_massRangeLower) continue;
292 // if(mass > m_massRangeHigher) continue;
293 if ( vtxfit->chisq() > chi ) continue;
294 // if(fabs((vtxfit->p4par()).m()-mlam) > delm) continue;
295 delm = fabs( ( vtxfit->p4par() ).m() - mlam );
296 chi = vtxfit->chisq();
297 ipion = ipi[i1];
298 iproton = ip[i2];
299 }
300 }
301 //
302 // Lammda check
303 //
304 if ( ipion == -1 || iproton == -1 ) return sc;
305 m_pass[4] += 1;
306
307 RecMdcKalTrack* pionTrk = ( *( evtRecTrkCol->begin() + ipion ) )->mdcKalTrack();
309 WTrackParameter wpiontrk( mpi, pionTrk->helix(), pionTrk->err() );
310
311 RecMdcKalTrack* protonTrk = ( *( evtRecTrkCol->begin() + iproton ) )->mdcKalTrack();
313 WTrackParameter wprotontrk( mp, protonTrk->helix(), protonTrk->err() );
314 vtxfit0->init();
315 vtxfit0->AddTrack( 0, wpiontrk );
316 vtxfit0->AddTrack( 1, wprotontrk );
317 vtxfit0->AddVertex( 0, vxpar, 0, 1 );
318 if ( !( vtxfit0->Fit( 0 ) ) ) return sc;
319 vtxfit0->BuildVirtualParticle( 0 );
320
321 vtxfit->init();
322 // vtxfit->setPrimaryVertex(bs);
323 vtxfit->AddTrack( 0, vtxfit0->wVirtualTrack( 0 ) );
324 vtxfit->setVpar( vtxfit0->vpar( 0 ) );
325 if ( !vtxfit->Fit() ) return sc;
326
327 m_lamRawMass = vtxfit->p4par().m();
328 m_ctau = vtxfit->ctau();
329 m_lxyz = vtxfit->decayLength();
330 m_elxyz = vtxfit->decayLengthError();
331 m_chis = vtxfit->chisq();
332 m_plam = vtxfit->p4par().rho();
333
334 if ( m_lxyz > 3.5 )
335 {
336
337 TH1* h( 0 );
338 if ( m_thsvc->getHist( "/DQAHist/InclLam/InclLam_mass", h ).isSuccess() )
339 { h->Fill( m_lamRawMass ); }
340 else { log << MSG::ERROR << "Couldn't retrieve incllam_mass" << endmsg; }
341 }
342 m_tuple1->write();
343 ////////////////////////////////////
344 // DQA
345 // set tag and quality
346 // Pid: 1 - electron, 2 - muon, 3 - pion, 4 - kaon, 5 - proton
347 // (*(evtRecTrkCol->begin()+ipion))->setPartId(3);
348 // (*(evtRecTrkCol->begin()+iproton))->setPartId(5);
349 ( *( evtRecTrkCol->begin() + ipion ) )->tagPion();
350 ( *( evtRecTrkCol->begin() + iproton ) )->tagProton();
351 // Quality: defined by whether dE/dx or TOF is used to identify particle
352 // 0 - no dE/dx, no TOF (can be used for dE/dx and TOF calibration)
353 // 1 - only dE/dx (can be used for TOF calibration)
354 // 2 - only TOF (can be used for dE/dx calibration)
355 // 3 - Both dE/dx and TOF
356 ( *( evtRecTrkCol->begin() + ipion ) )->setQuality( 1 );
357 ( *( evtRecTrkCol->begin() + iproton ) )->setQuality( 1 );
358
359 //--------------------------------------------------
360 // Add the line below at the end of execute(), (before return)
361 //
362 setFilterPassed( true );
363
364 return StatusCode::SUCCESS;
365}
HepGeom::Point3D< double > HepPoint3D
EvtRecTrackCol::iterator EvtRecTrackIterator
double mpi
double mp
std::vector< HepLorentzVector > Vp4
Definition Gam4pikp.cxx:38
std::vector< int > Vint
Definition Gam4pikp.cxx:37
IMessageSvc * msgSvc()
static ParticleID * instance()
bool IsPidInfoValid() const
void calculate()
void init()
static SecondVertexFit * instance()
void AddTrack(const int number, const double mass, const RecMdcTrack *trk)
Definition TrackPool.cxx:21
WTrackParameter wVirtualTrack(int n) const
void init()
Definition VertexFit.cxx:27
void AddVertex(int number, VertexParameter vpar, std::vector< int > lis)
Definition VertexFit.cxx:85
static VertexFit * instance()
Definition VertexFit.cxx:15
void BuildVirtualParticle(int number)
bool Fit()
const double mlam

◆ finalize()

StatusCode incllambda::finalize ( )

Definition at line 368 of file incllambda.cxx.

368 {
369
370 MsgStream log( msgSvc(), name() );
371 log << MSG::INFO << "in finalize()" << endmsg;
372 log << MSG::INFO << "Total Entries : " << m_pass[0] << endmsg;
373 log << MSG::INFO << "Qtot Cut : " << m_pass[1] << endmsg;
374 log << MSG::INFO << "Nks : " << m_pass[2] << endmsg;
375 log << MSG::INFO << "Good Ks cut : " << m_pass[3] << endmsg;
376 return StatusCode::SUCCESS;
377}

◆ initialize()

StatusCode incllambda::initialize ( )

Definition at line 75 of file incllambda.cxx.

75 {
76 MsgStream log( msgSvc(), name() );
77
78 log << MSG::INFO << "in initialize()" << endmsg;
79
80 StatusCode status;
81
82 if ( service( "THistSvc", m_thsvc ).isFailure() )
83 {
84 log << MSG::ERROR << "Couldn't get THistSvc" << endmsg;
85 return StatusCode::FAILURE;
86 }
87 // "DQAHist" is fixed
88 TH1F* incllambda_mass =
89 new TH1F( "InclLam_mass", "INCLUSIVE_LAMBDA_MASS", 150, 1.11, 1.125 );
90 incllambda_mass->GetXaxis()->SetTitle( "M_{p#pi} (GeV)" );
91 incllambda_mass->GetYaxis()->SetTitle( "Nentries/0.1MeV/c^{2}" );
92 if ( m_thsvc->regHist( "/DQAHist/InclLam/InclLam_mass", incllambda_mass ).isFailure() )
93 { log << MSG::ERROR << "Couldn't register incllambda_mass" << endmsg; }
94
95 //*****************************************
96
97 NTuplePtr nt1( ntupleSvc(), "DQAFILE/InclLam" );
98 if ( nt1 ) m_tuple1 = nt1;
99 else
100 {
101 m_tuple1 = ntupleSvc()->book( "DQAFILE/InclLam", CLID_ColumnWiseTuple, "incllam Ntuple" );
102 if ( m_tuple1 )
103 {
104 status = m_tuple1->addItem( "npi", m_npi );
105 status = m_tuple1->addItem( "np", m_np );
106 status = m_tuple1->addItem( "lxyz", m_lxyz );
107 status = m_tuple1->addItem( "ctau", m_ctau );
108 status = m_tuple1->addItem( "elxyz", m_elxyz );
109 status = m_tuple1->addItem( "lamchi", m_chis );
110 status = m_tuple1->addItem( "mlamraw", m_lamRawMass );
111 status = m_tuple1->addItem( "plam", m_plam );
112 }
113 else
114 {
115 log << MSG::ERROR << " Cannot book N-tuple:" << long( m_tuple1 ) << endmsg;
116 return StatusCode::FAILURE;
117 }
118 }
119 //
120 //--------end of book--------
121 //
122
123 log << MSG::INFO << "successfully return from initialize()" << endmsg;
124 return StatusCode::SUCCESS;
125}
INTupleSvc * ntupleSvc()

The documentation for this class was generated from the following files: