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

#include <inclks.h>

Inheritance diagram for inclks:

Public Member Functions

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

Detailed Description

Definition at line 8 of file inclks.h.

Constructor & Destructor Documentation

◆ inclks()

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

Definition at line 48 of file inclks.cxx.

49 : Algorithm( name, pSvcLocator ) {
50
51 m_tuple1 = 0;
52 for ( int i = 0; i < 10; i++ ) m_pass[i] = 0;
53
54 // Declare the properties
55 declareProperty( "Vr0cut", m_vr0cut = 10.0 );
56 declareProperty( "Vz0cut", m_vz0cut = 50.0 );
57 declareProperty( "CheckDedx", m_checkDedx = 1 );
58 declareProperty( "CheckTof", m_checkTof = 1 );
59}

Referenced by inclks().

Member Function Documentation

◆ execute()

StatusCode inclks::execute ( )

Definition at line 114 of file inclks.cxx.

114 {
115 StatusCode sc = StatusCode::SUCCESS;
116
117 MsgStream log( msgSvc(), name() );
118 log << MSG::INFO << "in execute()" << endmsg;
119
120 // DQA
121 // Add the line below at the beginning of execute()
122 //
123 setFilterPassed( false );
124
125 m_pass[0] += 1;
126
127 SmartDataPtr<Event::EventHeader> eventHeader( eventSvc(), "/Event/EventHeader" );
128 SmartDataPtr<EvtRecEvent> evtRecEvent( eventSvc(), EventModel::EvtRec::EvtRecEvent );
129 SmartDataPtr<EvtRecTrackCol> evtRecTrkCol( eventSvc(), EventModel::EvtRec::EvtRecTrackCol );
130
131 Vint iks, ipip, ipim, iGood;
132 iGood.clear();
133 iks.clear();
134 ipip.clear();
135 ipim.clear();
136
137 Vp4 ppip, ppim;
138 ppip.clear();
139 ppim.clear();
140
141 int TotCharge = 0;
142 for ( int i = 0; i < evtRecEvent->totalCharged(); i++ )
143 {
144 EvtRecTrackIterator itTrk = evtRecTrkCol->begin() + i;
145 if ( !( *itTrk )->isMdcTrackValid() ) continue;
146 RecMdcTrack* mdcTrk = ( *itTrk )->mdcTrack();
147 if ( fabs( mdcTrk->z() ) >= m_vz0cut ) continue;
148 if ( mdcTrk->r() >= m_vr0cut ) continue;
149 iGood.push_back( i );
150 TotCharge += mdcTrk->charge();
151 }
152 //
153 // Finish Good Charged Track Selection
154 //
155 int nGood = iGood.size();
156
157 //
158 // Charge track number cut
159 //
160
161 if ( ( nGood < 2 ) || ( TotCharge != 0 ) ) return sc;
162
163 m_pass[1] += 1;
164 //
165 // Assign 4-momentum to each charged track
166 //
167 ParticleID* pid = ParticleID::instance();
168 for ( int i = 0; i < nGood; i++ )
169 {
170 EvtRecTrackIterator itTrk = evtRecTrkCol->begin() + iGood[i];
171 // if(pid) delete pid;
172 pid->init();
173 pid->setMethod( pid->methodProbability() );
174 pid->setChiMinCut( 4 );
175 pid->setRecTrack( *itTrk );
176 pid->usePidSys( pid->useDedx() | pid->useTof1() | pid->useTof2() ); // use PID sub-system
177 pid->identify( pid->onlyPion() | pid->onlyKaon() ); // seperater Pion/Kaon
178 // pid->identify(pid->onlyPion());
179 // pid->identify(pid->onlyKaon());
180 pid->calculate();
181 if ( !( pid->IsPidInfoValid() ) ) continue;
182 if ( !( *itTrk )->isMdcKalTrackValid() ) continue;
183 RecMdcKalTrack* mdcKalTrk = ( *itTrk )->mdcKalTrack();
184 // RecMdcTrack* mdcTrk = (*itTrk)->mdcTrack();
185
186 if ( pid->probPion() < 0.001 || ( pid->probPion() < pid->probKaon() ) ) continue;
187 // if(pid->probPion() < 0.001) continue;
188 mdcKalTrk->setPidType( RecMdcKalTrack::pion );
189 HepLorentzVector ptrk;
190 ptrk.setPx( mdcKalTrk->px() );
191 ptrk.setPy( mdcKalTrk->py() );
192 ptrk.setPz( mdcKalTrk->pz() );
193 double p3 = ptrk.mag();
194 ptrk.setE( sqrt( p3 * p3 + xmass[2] * xmass[2] ) );
195
196 if ( mdcKalTrk->charge() > 0 )
197 {
198 ipip.push_back( iGood[i] );
199 ppip.push_back( ptrk );
200 }
201 else
202 {
203 ipim.push_back( iGood[i] );
204 ppim.push_back( ptrk );
205 }
206 }
207
208 m_pass[2] += 1;
209 int npip = ipip.size();
210 int npim = ipim.size();
211 m_npip = npip;
212 m_npim = npim;
213
214 if ( npip < 1 || npim < 1 ) return sc;
215
216 m_pass[3] += 1;
217
218 //
219 //****** Second Vertex Check************
220 //
221 double chi, delm;
222 double chisq = 999.;
223 // double delm=100.;
224
225 HepPoint3D vx( 0., 0., 0. );
226 HepSymMatrix Evx( 3, 0 );
227 double bx = 1E+6;
228 double by = 1E+6;
229 double bz = 1E+6;
230 Evx[0][0] = bx * bx;
231 Evx[1][1] = by * by;
232 Evx[2][2] = bz * bz;
233 VertexParameter vxpar;
234 vxpar.setVx( vx );
235 vxpar.setEvx( Evx );
236
237 // HepPoint3D bvx(0., 0., 0.);
238 // HepSymMatrix bEvx(3, 0);
239 // bEvx[0][0] = 0.038*0.038;
240 // bEvx[1][1] = 0.00057*0.00057;
241 // bEvx[2][2] = 1.5*1.5;
242 // VertexParameter bs;
243 // bs.setVx(bvx);
244 // bs.setEvx(bEvx);
245
246 int ip1 = -1, ip2 = -1;
247
248 VertexFit* vtxfit0 = VertexFit::instance();
249 SecondVertexFit* vtxfit = SecondVertexFit::instance();
250
251 for ( int i1 = 0; i1 < ipip.size(); i1++ )
252 {
253 RecMdcKalTrack* pipTrk = ( *( evtRecTrkCol->begin() + ipip[i1] ) )->mdcKalTrack();
255 WTrackParameter wpiptrk( xmass[2], pipTrk->helix(), pipTrk->err() );
256
257 for ( int i2 = 0; i2 < ipim.size(); i2++ )
258 {
259 RecMdcKalTrack* pimTrk = ( *( evtRecTrkCol->begin() + ipim[i2] ) )->mdcKalTrack();
261 WTrackParameter wpimtrk( xmass[2], pimTrk->helix(), pimTrk->err() );
262
263 vtxfit0->init();
264 vtxfit0->AddTrack( 0, wpiptrk );
265 vtxfit0->AddTrack( 1, wpimtrk );
266 vtxfit0->AddVertex( 0, vxpar, 0, 1 );
267 if ( !( vtxfit0->Fit( 0 ) ) ) continue;
268 vtxfit0->BuildVirtualParticle( 0 );
269
270 vtxfit->init();
271 vtxfit->AddTrack( 0, vtxfit0->wVirtualTrack( 0 ) );
272 vtxfit->setVpar( vtxfit0->vpar( 0 ) );
273 if ( !vtxfit->Fit() ) continue;
274 chi = vtxfit->chisq();
275
276 // if(fabs((vtxfit->p4par().m())-mk0) > delm) continue;
277 if ( chi > chisq ) continue;
278 delm = fabs( ( vtxfit->p4par().m() ) - mk0 );
279 chisq = chi;
280 ip1 = ipip[i1];
281 ip2 = ipim[i2];
282 }
283 }
284 //
285 // Kshort check
286 //
287 if ( ip1 == -1 || ip2 == -1 ) return sc;
288 m_pass[4] += 1;
289
290 RecMdcKalTrack* pi1Trk = ( *( evtRecTrkCol->begin() + ip1 ) )->mdcKalTrack();
292 WTrackParameter wpi1trk( xmass[2], pi1Trk->helix(), pi1Trk->err() );
293
294 RecMdcKalTrack* pi2Trk = ( *( evtRecTrkCol->begin() + ip2 ) )->mdcKalTrack();
296 WTrackParameter wpi2trk( xmass[2], pi2Trk->helix(), pi2Trk->err() );
297 vtxfit0->init();
298 vtxfit0->AddTrack( 0, wpi1trk );
299 vtxfit0->AddTrack( 1, wpi2trk );
300 vtxfit0->AddVertex( 0, vxpar, 0, 1 );
301 if ( !( vtxfit0->Fit( 0 ) ) ) return sc;
302 vtxfit0->BuildVirtualParticle( 0 );
303
304 vtxfit->init();
305 vtxfit->AddTrack( 0, vtxfit0->wVirtualTrack( 0 ) );
306 vtxfit->setVpar( vtxfit0->vpar( 0 ) );
307 if ( !vtxfit->Fit() ) return sc;
308
309 m_ksRawMass = vtxfit->p4par().m();
310 m_ctau = vtxfit->ctau();
311 m_lxyz = vtxfit->decayLength();
312 m_elxyz = vtxfit->decayLengthError();
313 m_chis = vtxfit->chisq();
314 m_pk0 = vtxfit->p4par().rho();
315
316 // DQA
317 if ( m_lxyz > 0.4 && m_chis < 20. )
318 {
319
320 TH1* h( 0 );
321 if ( m_thsvc->getHist( "/DQAHist/InclKs/InclKs_mass", h ).isSuccess() )
322 { h->Fill( m_ksRawMass ); }
323 else { log << MSG::ERROR << "Couldn't retrieve inclks_mass" << endmsg; }
324 }
325
326 m_tuple1->write().ignore();
327 ////////////////////////////////////
328 // DQA
329 // set tag and quality
330 // Pid: 1 - electron, 2 - muon, 3 - pion, 4 - kaon, 5 - proton
331 // (*(evtRecTrkCol->begin()+ip1))->setPartId(3);
332 // (*(evtRecTrkCol->begin()+ip2))->setPartId(3);
333 ( *( evtRecTrkCol->begin() + ip1 ) )->tagPion();
334 ( *( evtRecTrkCol->begin() + ip2 ) )->tagPion();
335 // Quality: defined by whether dE/dx or TOF is used to identify particle
336 // 0 - no dE/dx, no TOF (can be used for dE/dx and TOF calibration)
337 // 1 - only dE/dx (can be used for TOF calibration)
338 // 2 - only TOF (can be used for dE/dx calibration)
339 // 3 - Both dE/dx and TOF
340 ( *( evtRecTrkCol->begin() + ip1 ) )->setQuality( 2 );
341 ( *( evtRecTrkCol->begin() + ip2 ) )->setQuality( 2 );
342 //--------------------------------------------------
343 // Add the line below at the end of execute(), (before return)
344 //
345 setFilterPassed( true );
346
347 return StatusCode::SUCCESS;
348}
HepGeom::Point3D< double > HepPoint3D
EvtRecTrackCol::iterator EvtRecTrackIterator
std::vector< HepLorentzVector > Vp4
Definition Gam4pikp.cxx:38
const double xmass[5]
Definition Gam4pikp.cxx:35
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 mk0
Definition inclks.cxx:33

◆ finalize()

StatusCode inclks::finalize ( )

Definition at line 351 of file inclks.cxx.

351 {
352
353 MsgStream log( msgSvc(), name() );
354 log << MSG::INFO << "in finalize()" << endmsg;
355 log << MSG::INFO << "Total Entries : " << m_pass[0] << endmsg;
356 log << MSG::INFO << "Qtot Cut : " << m_pass[1] << endmsg;
357 log << MSG::INFO << "PID : " << m_pass[2] << endmsg;
358 log << MSG::INFO << "NPI : " << m_pass[3] << endmsg;
359 log << MSG::INFO << "Second Vertex Cut : " << m_pass[4] << endmsg;
360 return StatusCode::SUCCESS;
361}

◆ initialize()

StatusCode inclks::initialize ( )

Definition at line 61 of file inclks.cxx.

61 {
62 MsgStream log( msgSvc(), name() );
63
64 log << MSG::INFO << "in initialize()" << endmsg;
65
66 StatusCode status;
67
68 if ( service( "THistSvc", m_thsvc ).isFailure() )
69 {
70 log << MSG::ERROR << "Couldn't get THistSvc" << endmsg;
71 return StatusCode::FAILURE;
72 }
73 // "DQAHist" is fixed
74 TH1F* inclks_mass = new TH1F( "InclKs_mass", "INCLUSIVE_Ks_MASS", 70, 0.46, 0.53 );
75 inclks_mass->GetXaxis()->SetTitle( "M_{#pi#pi} (GeV)" );
76 inclks_mass->GetYaxis()->SetTitle( "Nentries/0.1MeV/c^{2}" );
77 if ( m_thsvc->regHist( "/DQAHist/InclKs/InclKs_mass", inclks_mass ).isFailure() )
78 { log << MSG::ERROR << "Couldn't register inclks_mass" << endmsg; }
79
80 //*****************************************
81
82 NTuplePtr nt1( ntupleSvc(), "DQAFILE/InclKs" );
83 if ( nt1 ) m_tuple1 = nt1;
84 else
85 {
86 m_tuple1 = ntupleSvc()->book( "DQAFILE/InclKs", CLID_ColumnWiseTuple, "ksklx Ntuple" );
87 if ( m_tuple1 )
88 {
89 status = m_tuple1->addItem( "npip", m_npip );
90 status = m_tuple1->addItem( "npim", m_npim );
91 status = m_tuple1->addItem( "ctau", m_ctau );
92 status = m_tuple1->addItem( "lxyz", m_lxyz );
93 status = m_tuple1->addItem( "elxyz", m_elxyz );
94 status = m_tuple1->addItem( "kschi", m_chis );
95 status = m_tuple1->addItem( "mksraw", m_ksRawMass );
96 status = m_tuple1->addItem( "pk0", m_pk0 );
97 }
98 else
99 {
100 log << MSG::ERROR << " Cannot book N-tuple:" << long( m_tuple1 ) << endmsg;
101 return StatusCode::FAILURE;
102 }
103 }
104
105 //
106 //--------end of book--------
107 //
108
109 log << MSG::INFO << "successfully return from initialize()" << endmsg;
110 return StatusCode::SUCCESS;
111}
INTupleSvc * ntupleSvc()

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