BOSS 8.0.0
BESIII Offline Software System
Loading...
Searching...
No Matches
DQAKsKpi.cxx
Go to the documentation of this file.
1#include "CLHEP/Vector/LorentzVector.h"
2#include "CLHEP/Vector/ThreeVector.h"
3#include "CLHEP/Vector/TwoVector.h"
4#include "GaudiKernel/Bootstrap.h"
5#include "GaudiKernel/IDataProviderSvc.h"
6#include "GaudiKernel/INTupleSvc.h"
7#include "GaudiKernel/ISvcLocator.h"
8#include "GaudiKernel/ITHistSvc.h"
9#include "GaudiKernel/MsgStream.h"
10#include "GaudiKernel/NTuple.h"
11#include "GaudiKernel/PropertyMgr.h"
12#include "GaudiKernel/SmartDataPtr.h"
13#include "TH1F.h"
14#include "TMath.h"
15#include <vector>
16
17#include "DstEvent/TofHitStatus.h"
18#include "EventModel/Event.h"
19#include "EventModel/EventHeader.h"
20#include "EventModel/EventModel.h"
21#include "EvtRecEvent/EvtRecEvent.h"
22#include "EvtRecEvent/EvtRecTrack.h"
23#include "ParticleID/ParticleID.h"
24#include "VertexDbSvc/IVertexDbSvc.h"
25#include "VertexFit/Helix.h"
26#include "VertexFit/KinematicFit.h"
27#include "VertexFit/SecondVertexFit.h"
28#include "VertexFit/VertexFit.h"
29
30using CLHEP::Hep2Vector;
31using CLHEP::Hep3Vector;
32using CLHEP::HepLorentzVector;
33#include "CLHEP/Geometry/Point3D.h"
34#ifndef ENABLE_BACKWARDS_COMPATIBILITY
35typedef HepGeom::Point3D<double> HepPoint3D;
36#endif
37
38#include "DQAKsKpi.h"
39
40const double mpi = 0.13957;
41const double mk = 0.493677;
42const double mks0 = 0.497648;
43const double xmass[5] = { 0.000511, 0.105658, 0.139570, 0.493677, 0.938272 };
44const double velc = 299.792458; // tof path unit in mm
45typedef std::vector<int> Vint;
46typedef std::vector<HepLorentzVector> Vp4;
47
48// boost
49const HepLorentzVector p_cms( 0.034067, 0.0, 0.0, 3.097 );
50const Hep3Vector u_cms = -p_cms.boostVector();
51// declare one counter
52static int counter[10] = { 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 };
53
54/////////////////////////////////////////////////////////////////////////////
56DQAKsKpi::DQAKsKpi( const std::string& name, ISvcLocator* pSvcLocator )
57 : Algorithm( name, pSvcLocator ) {
58
59 // Declare the properties
60 declareProperty( "Vr0cut", m_vr0cut = 5.0 );
61 declareProperty( "Vz0cut", m_vz0cut = 20.0 );
62 declareProperty( "Vr1cut", m_vr1cut = 1.0 );
63 declareProperty( "Vz1cut", m_vz1cut = 5.0 );
64 declareProperty( "Vctcut", m_cthcut = 0.93 );
65 declareProperty( "EnergyThreshold", m_energyThreshold = 0.04 );
66 declareProperty( "GammaAngCut", m_gammaAngCut = 20.0 );
67}
68
69// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
71 MsgStream log( msgSvc(), name() );
72
73 log << MSG::INFO << "in initialize()" << endmsg;
74 StatusCode status;
75
76 if ( service( "THistSvc", m_thsvc ).isFailure() )
77 {
78 log << MSG::ERROR << "Couldn't get THistSvc" << endmsg;
79 return StatusCode::FAILURE;
80 }
81
82 // "DQAHist" is fixed
83 // "DQAKsKpi" is control sample name, just as ntuple case.
84 TH1F* hks_dl = new TH1F( "ks_dl", "ks_dl", 300, -5.0, 25.0 );
85 if ( m_thsvc->regHist( "/DQAHist/DQAKsKpi/hks_dl", hks_dl ).isFailure() )
86 { log << MSG::ERROR << "Couldn't register ks_dl" << endmsg; }
87
88 TH1F* hks_m = new TH1F( "ks_m", "ks_m", 200, 0.4, 0.6 );
89 if ( m_thsvc->regHist( "/DQAHist/DQAKsKpi/hks_m", hks_m ).isFailure() )
90 { log << MSG::ERROR << "Couldn't register ks_m" << endmsg; }
91
92 TH1F* hkspi_m = new TH1F( "kspi_m", "kspi_m", 200, 0.6, 2.6 );
93 if ( m_thsvc->regHist( "/DQAHist/DQAKsKpi/hkspi_m", hkspi_m ).isFailure() )
94 { log << MSG::ERROR << "Couldn't register kspi_m" << endmsg; }
95
96 TH1F* hks_p = new TH1F( "ks_p", "ks_p", 100, 0.0, 1.5 );
97 if ( m_thsvc->regHist( "/DQAHist/DQAKsKpi/hks_p", hks_p ).isFailure() )
98 { log << MSG::ERROR << "Couldn't register ks_p" << endmsg; }
99
100 TH1F* hkpi_m = new TH1F( "kpi_m", "kpi_m", 200, 0.6, 2.6 );
101 if ( m_thsvc->regHist( "/DQAHist/DQAKsKpi/hkpi_m", hkpi_m ).isFailure() )
102 { log << MSG::ERROR << "Couldn't register kpi_m" << endmsg; }
103
104 // DQA
105 // The first directory specifier must be "DQAFILE"
106 // The second is the control sample name, the first letter is in upper format. eg. "Rhopi"
107 NTuplePtr nt( ntupleSvc(), "DQAFILE/KsKpi" );
108 if ( nt ) m_tuple = nt;
109 else
110 {
111 m_tuple = ntupleSvc()->book( "DQAFILE/KsKpi", CLID_ColumnWiseTuple, "KsKpi ntuple" );
112 if ( m_tuple )
113 {
114 status = m_tuple->addItem( "runNo", m_runNo );
115 status = m_tuple->addItem( "event", m_event );
116 // status = m_tuple->addItem("Nchrg", m_nchrg);
117 // status = m_tuple->addItem("Nneu", m_nneu);
118
119 status = m_tuple->addItem( "npip", m_npip );
120 status = m_tuple->addItem( "npim", m_npim );
121 status = m_tuple->addItem( "nkp", m_nkp );
122 status = m_tuple->addItem( "nkm", m_nkm );
123 status = m_tuple->addItem( "np", m_np );
124 status = m_tuple->addItem( "npb", m_npb );
125
126 status = m_tuple->addItem( "vfits_chi", m_vfits_chi );
127 status = m_tuple->addItem( "vfits_vx", m_vfits_vx );
128 status = m_tuple->addItem( "vfits_vy", m_vfits_vy );
129 status = m_tuple->addItem( "vfits_vz", m_vfits_vz );
130 status = m_tuple->addItem( "vfits_vr", m_vfits_vr );
131
132 status = m_tuple->addItem( "vfitp_chi", m_vfitp_chi );
133 status = m_tuple->addItem( "vfitp_vx", m_vfitp_vx );
134 status = m_tuple->addItem( "vfitp_vy", m_vfitp_vy );
135 status = m_tuple->addItem( "vfitp_vz", m_vfitp_vz );
136 status = m_tuple->addItem( "vfitp_vr", m_vfitp_vr );
137
138 status = m_tuple->addItem( "vfit2_chi", m_vfit2_chi );
139 status = m_tuple->addItem( "vfit2_mks", m_vfit2_mks );
140 status = m_tuple->addItem( "vfit2_ct", m_vfit2_ct );
141 status = m_tuple->addItem( "vfit2_dl", m_vfit2_dl );
142 status = m_tuple->addItem( "vfit2_dle", m_vfit2_dle );
143
144 status = m_tuple->addItem( "chi2_fs4c", m_chi2_fs4c );
145 status = m_tuple->addItem( "mks_fs4c", m_mks_fs4c );
146 status = m_tuple->addItem( "mkspi_fs4c", m_mkspi_fs4c );
147 status = m_tuple->addItem( "mksk_fs4c", m_mksk_fs4c );
148 status = m_tuple->addItem( "mkpi_fs4c", m_mkpi_fs4c );
149
150 status = m_tuple->addItem( "4c_chi2", m_4c_chi2 );
151 status = m_tuple->addItem( "4c_mks", m_4c_mks );
152 status = m_tuple->addItem( "4c_mkspi", m_4c_mkspi );
153 status = m_tuple->addItem( "4c_mksk", m_4c_mksk );
154 status = m_tuple->addItem( "4c_mkpi", m_4c_mkpi );
155 status = m_tuple->addItem( "4c_ks_px", m_4c_ks_px );
156 status = m_tuple->addItem( "4c_ks_py", m_4c_ks_py );
157 status = m_tuple->addItem( "4c_ks_pz", m_4c_ks_pz );
158 status = m_tuple->addItem( "4c_ks_p", m_4c_ks_p );
159 status = m_tuple->addItem( "4c_ks_cos", m_4c_ks_cos );
160
161 status = m_tuple->addItem( "NGch", m_ngch, 0, 10 );
162 status = m_tuple->addIndexedItem( "pidcode", m_ngch, m_pidcode );
163 status = m_tuple->addIndexedItem( "pidprob", m_ngch, m_pidprob );
164 status = m_tuple->addIndexedItem( "pidchiDedx", m_ngch, m_pidchiDedx );
165 status = m_tuple->addIndexedItem( "pidchiTof1", m_ngch, m_pidchiTof1 );
166 status = m_tuple->addIndexedItem( "pidchiTof2", m_ngch, m_pidchiTof2 );
167
168 status = m_tuple->addIndexedItem( "charge", m_ngch, m_charge );
169 status = m_tuple->addIndexedItem( "vx0", m_ngch, m_vx0 );
170 status = m_tuple->addIndexedItem( "vy0", m_ngch, m_vy0 );
171 status = m_tuple->addIndexedItem( "vz0", m_ngch, m_vz0 );
172 status = m_tuple->addIndexedItem( "vr0", m_ngch, m_vr0 );
173
174 status = m_tuple->addIndexedItem( "vx", m_ngch, m_vx );
175 status = m_tuple->addIndexedItem( "vy", m_ngch, m_vy );
176 status = m_tuple->addIndexedItem( "vz", m_ngch, m_vz );
177 status = m_tuple->addIndexedItem( "vr", m_ngch, m_vr );
178
179 status = m_tuple->addIndexedItem( "px", m_ngch, m_px );
180 status = m_tuple->addIndexedItem( "py", m_ngch, m_py );
181 status = m_tuple->addIndexedItem( "pz", m_ngch, m_pz );
182 status = m_tuple->addIndexedItem( "p", m_ngch, m_p );
183 status = m_tuple->addIndexedItem( "cost", m_ngch, m_cost );
184
185 status = m_tuple->addIndexedItem( "probPH", m_ngch, m_probPH );
186 status = m_tuple->addIndexedItem( "normPH", m_ngch, m_normPH );
187 status = m_tuple->addIndexedItem( "chie", m_ngch, m_chie );
188 status = m_tuple->addIndexedItem( "chimu", m_ngch, m_chimu );
189 status = m_tuple->addIndexedItem( "chipi", m_ngch, m_chipi );
190 status = m_tuple->addIndexedItem( "chik", m_ngch, m_chik );
191 status = m_tuple->addIndexedItem( "chip", m_ngch, m_chip );
192 status = m_tuple->addIndexedItem( "ghit", m_ngch, m_ghit );
193 status = m_tuple->addIndexedItem( "thit", m_ngch, m_thit );
194
195 status = m_tuple->addIndexedItem( "e_emc", m_ngch, m_e_emc );
196
197 status = m_tuple->addIndexedItem( "qual_etof", m_ngch, m_qual_etof );
198 status = m_tuple->addIndexedItem( "tof_etof", m_ngch, m_tof_etof );
199 status = m_tuple->addIndexedItem( "te_etof", m_ngch, m_te_etof );
200 status = m_tuple->addIndexedItem( "tmu_etof", m_ngch, m_tmu_etof );
201 status = m_tuple->addIndexedItem( "tpi_etof", m_ngch, m_tpi_etof );
202 status = m_tuple->addIndexedItem( "tk_etof", m_ngch, m_tk_etof );
203 status = m_tuple->addIndexedItem( "tp_etof", m_ngch, m_tp_etof );
204
205 status = m_tuple->addIndexedItem( "qual_btof1", m_ngch, m_qual_btof1 );
206 status = m_tuple->addIndexedItem( "tof_btof1", m_ngch, m_tof_btof1 );
207 status = m_tuple->addIndexedItem( "te_btof1", m_ngch, m_te_btof1 );
208 status = m_tuple->addIndexedItem( "tmu_btof1", m_ngch, m_tmu_btof1 );
209 status = m_tuple->addIndexedItem( "tpi_btof1", m_ngch, m_tpi_btof1 );
210 status = m_tuple->addIndexedItem( "tk_btof1", m_ngch, m_tk_btof1 );
211 status = m_tuple->addIndexedItem( "tp_btof1", m_ngch, m_tp_btof1 );
212
213 status = m_tuple->addIndexedItem( "qual_btof2", m_ngch, m_qual_btof2 );
214 status = m_tuple->addIndexedItem( "tof_btof2", m_ngch, m_tof_btof2 );
215 status = m_tuple->addIndexedItem( "te_btof2", m_ngch, m_te_btof2 );
216 status = m_tuple->addIndexedItem( "tmu_btof2", m_ngch, m_tmu_btof2 );
217 status = m_tuple->addIndexedItem( "tpi_btof2", m_ngch, m_tpi_btof2 );
218 status = m_tuple->addIndexedItem( "tk_btof2", m_ngch, m_tk_btof2 );
219 status = m_tuple->addIndexedItem( "tp_btof2", m_ngch, m_tp_btof2 );
220 }
221 else { log << MSG::ERROR << "Can not book N-tuple:" << long( m_tuple ) << endmsg; }
222 }
223
224 //--------end of book--------
225
226 log << MSG::INFO << "successfully return from initialize()" << endmsg;
227 return StatusCode::SUCCESS;
228}
229
230// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
231StatusCode DQAKsKpi::execute() {
232
233 MsgStream log( msgSvc(), name() );
234 log << MSG::INFO << "in execute()" << endmsg;
235
236 // DQA
237 // Add the line below at the beginning of execute()
238 // setFilterPassed(false);
239
240 SmartDataPtr<Event::EventHeader> eventHeader( eventSvc(), "/Event/EventHeader" );
241 m_runNo = eventHeader->runNumber();
242 m_event = eventHeader->eventNumber();
243 log << MSG::DEBUG << "run, evtnum = " << m_runNo << " , " << m_event << endmsg;
244
245 SmartDataPtr<EvtRecEvent> evtRecEvent( eventSvc(), EventModel::EvtRec::EvtRecEvent );
246 // m_nchrg = evtRecEvent->totalCharged();
247 // m_nneu = evtRecEvent->totalNeutral();
248 log << MSG::DEBUG << "ncharg, nneu, tottks = " << evtRecEvent->totalCharged() << " , "
249 << evtRecEvent->totalNeutral() << " , " << evtRecEvent->totalTracks() << endmsg;
250
251 SmartDataPtr<EvtRecTrackCol> evtRecTrkCol( eventSvc(), EventModel::EvtRec::EvtRecTrackCol );
252
253 if ( evtRecEvent->totalNeutral() > 100 ) { return StatusCode::SUCCESS; }
254
255 Vint iGood;
256 iGood.clear();
257
258 Hep3Vector xorigin( 0, 0, 0 );
259
260 IVertexDbSvc* vtxsvc;
261 Gaudi::svcLocator()->service( "VertexDbSvc", vtxsvc ).ignore();
262 if ( vtxsvc->isVertexValid() )
263 {
264 double* dbv = vtxsvc->PrimaryVertex();
265 double* vv = vtxsvc->SigmaPrimaryVertex();
266 xorigin.setX( dbv[0] );
267 xorigin.setY( dbv[1] );
268 xorigin.setZ( dbv[2] );
269 log << MSG::INFO << "xorigin.x=" << xorigin.x() << ", "
270 << "xorigin.y=" << xorigin.y() << ", "
271 << "xorigin.z=" << xorigin.z() << ", " << endmsg;
272 }
273
274 int nCharge = 0;
275 for ( int i = 0; i < evtRecEvent->totalCharged(); i++ )
276 {
277 EvtRecTrackIterator itTrk = evtRecTrkCol->begin() + i;
278 if ( !( *itTrk )->isMdcTrackValid() ) continue;
279 if ( !( *itTrk )->isMdcKalTrackValid() ) continue;
280 RecMdcTrack* mdcTrk = ( *itTrk )->mdcTrack();
281
282 double pch = mdcTrk->p();
283 double x0 = mdcTrk->x();
284 double y0 = mdcTrk->y();
285 double z0 = mdcTrk->z();
286 double phi0 = mdcTrk->helix( 1 );
287 double xv = xorigin.x();
288 double yv = xorigin.y();
289 double Rxy = fabs( ( x0 - xv ) * cos( phi0 ) + ( y0 - yv ) * sin( phi0 ) );
290 double vx0 = x0;
291 double vy0 = y0;
292 double vz0 = z0 - xorigin.z();
293 double vr0 = Rxy;
294 double Vct = cos( mdcTrk->theta() );
295
296 HepVector a = mdcTrk->helix();
297 HepSymMatrix Ea = mdcTrk->err();
298 HepPoint3D point0( 0., 0., 0. ); // the initial point for MDC recosntruction
299 HepPoint3D IP( xorigin[0], xorigin[1], xorigin[2] );
300 VFHelix helixip( point0, a, Ea );
301 helixip.pivot( IP );
302 HepVector vecipa = helixip.a();
303 double Rvxy0 = fabs( vecipa[0] ); // the nearest distance to IP in xy plane
304 double Rvz0 = vecipa[3]; // the nearest distance to IP in z direction
305 double Rvphi0 = vecipa[1];
306 // m_rvxy0=Rvxy0;
307 // m_rvz0=Rvz0;
308 // m_rvphi0=Rvphi0;
309
310 if ( fabs( Rvz0 ) >= m_vz0cut ) continue;
311 if ( Rvxy0 >= m_vr0cut ) continue;
312 // if(fabs(Vct)>=m_cthcut) continue;
313 // iGood.push_back((*itTrk)->trackId());
314 iGood.push_back( i );
315 nCharge += mdcTrk->charge();
316 }
317
318 //
319 // Finish Good Charged Track Selection
320 //
321 int m_ngch = iGood.size();
322 log << MSG::DEBUG << "ngood, totcharge = " << m_ngch << " , " << nCharge << endmsg;
323 if ( ( m_ngch != 4 ) || ( nCharge != 0 ) ) { return StatusCode::SUCCESS; }
324
325 // Particle ID
326 //
327 Vint ipip, ipim, ikp, ikm, ipp, ipm;
328 ipip.clear();
329 ipim.clear();
330 ikp.clear();
331 ikm.clear();
332 ipp.clear();
333 ipm.clear();
334
335 Vp4 p_pip, p_pim, p_kp, p_km, p_pp, p_pm;
336 p_pip.clear();
337 p_pim.clear();
338 p_kp.clear();
339 p_km.clear();
340 p_pp.clear();
341 p_pm.clear();
342
344 for ( int i = 0; i < m_ngch; i++ )
345 {
346 EvtRecTrackIterator itTrk = evtRecTrkCol->begin() + iGood[i];
347 // if(pid) delete pid;
348 pid->init();
349 pid->setMethod( pid->methodProbability() );
350 pid->setChiMinCut( 4 );
351 pid->setRecTrack( *itTrk );
352 pid->usePidSys( pid->useDedx() ); // just use dedx PID
353 // pid->usePidSys(pid->useDedx() | pid->useTof1() | pid->useTof2() | pid->useTofE() |
354 // pid->useTofQ()); // use PID sub-system
355 pid->identify( pid->onlyPionKaonProton() ); // seperater Pion/Kaon/Proton
356 // pid->identify(pid->onlyPion() | pid->onlyKaon()); // seperater Pion/Kaon
357 // pid->identify(pid->onlyPion());
358 // pid->identify(pid->onlyKaon());
359 pid->calculate();
360 if ( !( pid->IsPidInfoValid() ) ) continue;
361
362 RecMdcTrack* mdcTrk = ( *itTrk )->mdcTrack();
363 RecMdcKalTrack* mdcKalTrk = 0;
364 if ( ( *itTrk )->isMdcKalTrackValid() ) mdcKalTrk = ( *itTrk )->mdcKalTrack();
365
366 double prob_pi = pid->probPion();
367 double prob_K = pid->probKaon();
368 double prob_p = pid->probProton();
369 // std::cout << "prob "<< prob_pi << ", "<< prob_K << ", "<< prob_p << std::endl;
370 HepLorentzVector ptrk;
371 ptrk.setPx( mdcTrk->px() );
372 ptrk.setPy( mdcTrk->py() );
373 ptrk.setPz( mdcTrk->pz() );
374 double p3 = ptrk.mag();
375
376 if ( prob_pi > prob_K && prob_pi > prob_p )
377 {
378 m_pidcode[i] = 2;
379 m_pidprob[i] = pid->prob( 2 );
380 m_pidchiDedx[i] = pid->chiDedx( 2 );
381 m_pidchiTof1[i] = pid->chiTof1( 2 );
382 m_pidchiTof2[i] = pid->chiTof2( 2 );
383 ptrk.setE( sqrt( p3 * p3 + xmass[2] * xmass[2] ) );
384 if ( mdcTrk->charge() > 0 )
385 {
386 ipip.push_back( iGood[i] );
387 p_pip.push_back( ptrk );
388 }
389 if ( mdcTrk->charge() < 0 )
390 {
391 ipim.push_back( iGood[i] );
392 p_pim.push_back( ptrk );
393 }
394 }
395
396 if ( prob_K > prob_pi && prob_K > prob_p )
397 {
398 m_pidcode[i] = 3;
399 m_pidprob[i] = pid->prob( 3 );
400 m_pidchiDedx[i] = pid->chiDedx( 3 );
401 m_pidchiTof1[i] = pid->chiTof1( 3 );
402 m_pidchiTof2[i] = pid->chiTof2( 3 );
403 ptrk.setE( sqrt( p3 * p3 + xmass[3] * xmass[3] ) );
404 if ( mdcTrk->charge() > 0 )
405 {
406 ikp.push_back( iGood[i] );
407 p_kp.push_back( ptrk );
408 }
409 if ( mdcTrk->charge() < 0 )
410 {
411 ikm.push_back( iGood[i] );
412 p_km.push_back( ptrk );
413 }
414 }
415
416 if ( prob_p > prob_pi && prob_p > prob_K )
417 {
418 m_pidcode[i] = 4;
419 m_pidprob[i] = pid->prob( 4 );
420 m_pidchiDedx[i] = pid->chiDedx( 4 );
421 m_pidchiTof1[i] = pid->chiTof1( 4 );
422 m_pidchiTof2[i] = pid->chiTof2( 4 );
423 ptrk.setE( sqrt( p3 * p3 + xmass[4] * xmass[4] ) );
424 if ( mdcTrk->charge() > 0 )
425 {
426 ipp.push_back( iGood[i] );
427 p_pp.push_back( ptrk );
428 }
429 if ( mdcTrk->charge() < 0 )
430 {
431 ipm.push_back( iGood[i] );
432 p_pm.push_back( ptrk );
433 }
434 }
435 }
436
437 m_npip = ipip.size();
438 m_npim = ipim.size();
439 m_nkp = ikp.size();
440 m_nkm = ikm.size();
441 m_np = ipp.size();
442 m_npb = ipm.size();
443
444 // cout<<"m_npip*m_npim: "<<m_npip*m_npim<<endl;
445 if ( m_npip * m_npim != 2 ) { return StatusCode::SUCCESS; }
446 // cout<<"m_nkp+m_nkm: "<<m_nkp+m_nkm<<endl;
447 if ( m_nkp + m_nkm != 1 ) { return StatusCode::SUCCESS; }
448 // cout<<"haha"<<endl;
449
450 // Vertex Fit
451 HepPoint3D vx( 0., 0., 0. );
452 HepSymMatrix Evx( 3, 0 );
453 double bx = 1E+6;
454 double by = 1E+6;
455 double bz = 1E+6;
456 Evx[0][0] = bx * bx;
457 Evx[1][1] = by * by;
458 Evx[2][2] = bz * bz;
459
460 VertexParameter vxpar;
461 vxpar.setVx( vx );
462 vxpar.setEvx( Evx );
463
464 VertexFit* vtxfit_s = VertexFit::instance(); // S second vertex
465 VertexFit* vtxfit_p = VertexFit::instance(); // P primary vertex
467
468 RecMdcKalTrack *pi1KalTrk{ nullptr }, *pi2KalTrk{ nullptr }, *pi3KalTrk{ nullptr },
469 *k1KalTrk{ nullptr };
470 RecMdcKalTrack *pipKalTrk{ nullptr }, *pimKalTrk{ nullptr }, *piKalTrk{ nullptr },
471 *kKalTrk{ nullptr };
472 WTrackParameter wks, vwks;
473
474 double chi_temp = 999.0;
475 double mks_temp = 10.0;
476 bool okloop = false;
477 for ( unsigned int i1 = 0; i1 < m_npip; i1++ )
478 { // pion plus
479 RecMdcKalTrack* pi1KalTrk = ( *( evtRecTrkCol->begin() + ipip[i1] ) )->mdcKalTrack();
480 pi1KalTrk->setPidType( RecMdcKalTrack::pion );
481 WTrackParameter wpi1trk( xmass[2], pi1KalTrk->getZHelix(), pi1KalTrk->getZError() );
482
483 for ( unsigned int i2 = 0; i2 < m_npim; i2++ )
484 { // pion minus
485 RecMdcKalTrack* pi2KalTrk = ( *( evtRecTrkCol->begin() + ipim[i2] ) )->mdcKalTrack();
486 pi2KalTrk->setPidType( RecMdcKalTrack::pion );
487 WTrackParameter wpi2trk( xmass[2], pi2KalTrk->getZHelix(), pi2KalTrk->getZError() );
488
489 vtxfit_s->init();
490 vtxfit_s->AddTrack( 0, wpi1trk );
491 vtxfit_s->AddTrack( 1, wpi2trk );
492 vtxfit_s->AddVertex( 0, vxpar, 0, 1 );
493
494 if ( !( vtxfit_s->Fit( 0 ) ) ) continue;
495 vtxfit_s->BuildVirtualParticle( 0 );
496 m_vfits_chi = vtxfit_s->chisq( 0 );
497 WTrackParameter wkshort = vtxfit_s->wVirtualTrack( 0 );
498 VertexParameter vparks = vtxfit_s->vpar( 0 );
499
500 m_vfits_vx = ( vparks.Vx() )[0];
501 m_vfits_vy = ( vparks.Vx() )[1];
502 m_vfits_vz = ( vparks.Vx() )[2];
503 m_vfits_vr = sqrt( m_vfits_vx * m_vfits_vx + m_vfits_vy * m_vfits_vy );
504
505 if ( m_npip == 2 )
506 {
507 int j = i1;
508 int jj = ( i1 == 1 ) ? 0 : 1;
509 pi3KalTrk = ( *( evtRecTrkCol->begin() + ipip[jj] ) )->mdcKalTrack();
510 k1KalTrk = ( *( evtRecTrkCol->begin() + ikm[0] ) )->mdcKalTrack();
511 }
512 if ( m_npim == 2 )
513 {
514 int j = i2;
515 int jj = ( i2 == 1 ) ? 0 : 1;
516 pi3KalTrk = ( *( evtRecTrkCol->begin() + ipim[jj] ) )->mdcKalTrack();
517 k1KalTrk = ( *( evtRecTrkCol->begin() + ikp[0] ) )->mdcKalTrack();
518 }
519
520 pi3KalTrk->setPidType( RecMdcKalTrack::pion );
521 WTrackParameter wpi3trk( xmass[2], pi3KalTrk->getZHelix(), pi3KalTrk->getZError() );
522 k1KalTrk->setPidType( RecMdcKalTrack::kaon );
523 WTrackParameter wk1trk( xmass[3], k1KalTrk->getZHelixK(), k1KalTrk->getZErrorK() );
524
525 vtxfit_p->init();
526 // vtxfit_p->AddTrack(0, wkshort);
527 vtxfit_p->AddTrack( 0, wpi3trk );
528 vtxfit_p->AddTrack( 1, wk1trk );
529 vtxfit_p->AddVertex( 0, vxpar, 0, 1 );
530 if ( !( vtxfit_p->Fit( 0 ) ) ) continue;
531
532 m_vfitp_chi = vtxfit_p->chisq( 0 );
533
534 VertexParameter primaryVpar = vtxfit_p->vpar( 0 );
535 m_vfitp_vx = ( primaryVpar.Vx() )[0];
536 m_vfitp_vy = ( primaryVpar.Vx() )[1];
537 m_vfitp_vz = ( primaryVpar.Vx() )[2];
538 m_vfitp_vr = sqrt( m_vfitp_vx * m_vfitp_vx + m_vfitp_vy * m_vfitp_vy );
539
540 vtxfit2->init();
541 vtxfit2->setPrimaryVertex( vtxfit_p->vpar( 0 ) );
542 vtxfit2->AddTrack( 0, wkshort );
543 vtxfit2->setVpar( vtxfit_s->vpar( 0 ) );
544 if ( !vtxfit2->Fit() ) continue;
545
546 if ( fabs( ( ( vtxfit2->wpar() ).p() ).m() - mks0 ) < mks_temp )
547 {
548 mks_temp = fabs( ( ( vtxfit2->wpar() ).p() ).m() - mks0 );
549
550 okloop = true;
551
552 wks = vtxfit2->wpar();
553 m_vfit2_mks = ( wks.p() ).m();
554 m_vfit2_chi = vtxfit2->chisq();
555 m_vfit2_ct = vtxfit2->ctau();
556 m_vfit2_dl = vtxfit2->decayLength();
557 m_vfit2_dle = vtxfit2->decayLengthError();
558
559 pipKalTrk = pi1KalTrk;
560 pimKalTrk = pi2KalTrk;
561 piKalTrk = pi3KalTrk;
562 kKalTrk = k1KalTrk;
563 }
564 }
565 }
566
567 if ( !okloop ) { return StatusCode::SUCCESS; }
568
569 pipKalTrk->setPidType( RecMdcKalTrack::pion );
570 pimKalTrk->setPidType( RecMdcKalTrack::pion );
571 piKalTrk->setPidType( RecMdcKalTrack::pion );
572 kKalTrk->setPidType( RecMdcKalTrack::kaon );
573
574 WTrackParameter wpip( xmass[2], pipKalTrk->getZHelix(),
575 pipKalTrk->getZError() ); // pi+ from Ks
576 WTrackParameter wpim( xmass[2], pimKalTrk->getZHelix(),
577 pimKalTrk->getZError() ); // pi- from Ks
578
579 WTrackParameter wpi( xmass[2], piKalTrk->getZHelix(), piKalTrk->getZError() );
580 WTrackParameter wk( xmass[3], kKalTrk->getZHelixK(), kKalTrk->getZErrorK() );
581
582 //
583 // check good charged track's infomation
584 int ii{ 0 };
585 for ( int j = 0; j < m_ngch; j++ )
586 {
587 m_charge[j] = 9999.0;
588 m_vx0[j] = 9999.0;
589 m_vy0[j] = 9999.0;
590 m_vz0[j] = 9999.0;
591 m_vr0[j] = 9999.0;
592
593 m_vx[j] = 9999.0;
594 m_vy[j] = 9999.0;
595 m_vz[j] = 9999.0;
596 m_vr[j] = 9999.0;
597
598 m_px[j] = 9999.0;
599 m_py[j] = 9999.0;
600 m_pz[j] = 9999.0;
601 m_p[j] = 9999.0;
602 m_cost[j] = 9999.0;
603
604 m_probPH[j] = 9999.0;
605 m_normPH[j] = 9999.0;
606 m_chie[j] = 9999.0;
607 m_chimu[j] = 9999.0;
608 m_chipi[j] = 9999.0;
609 m_chik[j] = 9999.0;
610 m_chip[j] = 9999.0;
611 m_ghit[j] = 9999.0;
612 m_thit[j] = 9999.0;
613
614 m_e_emc[j] = 9999.0;
615
616 m_qual_etof[j] = 9999.0;
617 m_tof_etof[j] = 9999.0;
618 m_te_etof[j] = 9999.0;
619 m_tmu_etof[j] = 9999.0;
620 m_tpi_etof[j] = 9999.0;
621 m_tk_etof[j] = 9999.0;
622 m_tp_etof[j] = 9999.0;
623
624 m_qual_btof1[j] = 9999.0;
625 m_tof_btof1[j] = 9999.0;
626 m_te_btof1[j] = 9999.0;
627 m_tmu_btof1[j] = 9999.0;
628 m_tpi_btof1[j] = 9999.0;
629 m_tk_btof1[j] = 9999.0;
630 m_tp_btof1[j] = 9999.0;
631
632 m_qual_btof2[j] = 9999.0;
633 m_tof_btof2[j] = 9999.0;
634 m_te_btof2[j] = 9999.0;
635 m_tmu_btof2[j] = 9999.0;
636 m_tpi_btof2[j] = 9999.0;
637 m_tk_btof2[j] = 9999.0;
638 m_tp_btof2[j] = 9999.0;
639
640 m_pidcode[j] = 9999.0;
641 m_pidprob[j] = 9999.0;
642 m_pidchiDedx[j] = 9999.0;
643 m_pidchiTof1[j] = 9999.0;
644 m_pidchiTof2[j] = 99999.0;
645 }
646
647 for ( int i = 0; i < m_ngch; i++ )
648 {
649
650 EvtRecTrackIterator itTrk = evtRecTrkCol->begin() + iGood[i];
651 if ( !( *itTrk )->isMdcTrackValid() ) continue; // MDC information
652 RecMdcTrack* mdcTrk = ( *itTrk )->mdcTrack();
653 RecMdcKalTrack* mdcKalTrk = ( *itTrk )->mdcKalTrack();
654
655 if ( mdcKalTrk == pipKalTrk )
656 {
657 ii = 0;
658 mdcKalTrk->setPidType( RecMdcKalTrack::pion );
659 }
660 if ( mdcKalTrk == pimKalTrk )
661 {
662 ii = 1;
663 mdcKalTrk->setPidType( RecMdcKalTrack::pion );
664 }
665 if ( mdcKalTrk == piKalTrk )
666 {
667 ii = 2;
668 mdcKalTrk->setPidType( RecMdcKalTrack::pion );
669 }
670 if ( mdcKalTrk == kKalTrk )
671 {
672 ii = 3;
673 mdcKalTrk->setPidType( RecMdcKalTrack::kaon );
674 }
675
676 m_charge[ii] = mdcTrk->charge();
677 double x0 = mdcTrk->x();
678 double y0 = mdcTrk->y();
679 double z0 = mdcTrk->z();
680 double phi0 = mdcTrk->helix( 1 );
681 double xv = xorigin.x();
682 double yv = xorigin.y();
683 double zv = xorigin.z();
684 double rv = ( x0 - xv ) * cos( phi0 ) + ( y0 - yv ) * sin( phi0 );
685
686 m_vx0[ii] = x0 - xv;
687 m_vy0[ii] = y0 - yv;
688 m_vz0[ii] = z0 - zv;
689 m_vr0[ii] = rv;
690
691 x0 = mdcKalTrk->x();
692 y0 = mdcKalTrk->y();
693 z0 = mdcKalTrk->z();
694 rv = ( x0 - xv ) * cos( phi0 ) + ( y0 - yv ) * sin( phi0 );
695
696 m_vx[ii] = x0 - xv;
697 m_vy[ii] = y0 - yv;
698 m_vz[ii] = z0 - zv;
699 m_vr[ii] = rv;
700
701 m_px[ii] = mdcKalTrk->px();
702 m_py[ii] = mdcKalTrk->py();
703 m_pz[ii] = mdcKalTrk->pz();
704 m_p[ii] = sqrt( m_px[ii] * m_px[ii] + m_py[ii] * m_py[ii] + m_pz[ii] * m_pz[ii] );
705 m_cost[ii] = m_pz[ii] / m_p[ii];
706
707 double ptrk = m_p[ii];
708 if ( ( *itTrk )->isMdcDedxValid() )
709 { // DEDX information
710 RecMdcDedx* dedxTrk = ( *itTrk )->mdcDedx();
711 m_probPH[ii] = dedxTrk->probPH();
712 m_normPH[ii] = dedxTrk->normPH();
713
714 m_chie[ii] = dedxTrk->chiE();
715 m_chimu[ii] = dedxTrk->chiMu();
716 m_chipi[ii] = dedxTrk->chiPi();
717 m_chik[ii] = dedxTrk->chiK();
718 m_chip[ii] = dedxTrk->chiP();
719 m_ghit[ii] = dedxTrk->numGoodHits();
720 m_thit[ii] = dedxTrk->numTotalHits();
721 }
722
723 if ( ( *itTrk )->isEmcShowerValid() )
724 {
725 RecEmcShower* emcTrk = ( *itTrk )->emcShower();
726 m_e_emc[ii] = emcTrk->energy();
727 }
728
729 if ( ( *itTrk )->isTofTrackValid() )
730 { // TOF information
731 SmartRefVector<RecTofTrack> tofTrkCol = ( *itTrk )->tofTrack();
732
733 SmartRefVector<RecTofTrack>::iterator iter_tof = tofTrkCol.begin();
734
735 for ( ; iter_tof != tofTrkCol.end(); iter_tof++ )
736 {
737 TofHitStatus* status = new TofHitStatus;
738 status->setStatus( ( *iter_tof )->status() );
739
740 if ( !( status->is_barrel() ) )
741 { // endcap
742 if ( !( status->is_counter() ) ) continue; // ?
743 if ( status->layer() != 0 ) continue; // layer1
744 double path = ( *iter_tof )->path(); // ?
745 double tof = ( *iter_tof )->tof();
746 double ph = ( *iter_tof )->ph();
747 double rhit = ( *iter_tof )->zrhit();
748 double qual = 0.0 + ( *iter_tof )->quality();
749 double cntr = 0.0 + ( *iter_tof )->tofID();
750 double texp[5];
751 for ( int j = 0; j < 5; j++ )
752 {
753 double gb = ptrk / xmass[j];
754 double beta = gb / sqrt( 1 + gb * gb );
755 texp[j] = path / beta / velc;
756 }
757
758 m_qual_etof[ii] = qual;
759 m_tof_etof[ii] = tof;
760 }
761 else
762 { // barrel
763 if ( !( status->is_counter() ) ) continue; // ?
764 if ( status->layer() == 1 )
765 { // layer1
766 double path = ( *iter_tof )->path(); // ?
767 double tof = ( *iter_tof )->tof();
768 double ph = ( *iter_tof )->ph();
769 double rhit = ( *iter_tof )->zrhit();
770 double qual = 0.0 + ( *iter_tof )->quality();
771 double cntr = 0.0 + ( *iter_tof )->tofID();
772 double texp[5];
773 for ( int j = 0; j < 5; j++ )
774 {
775 double gb = ptrk / xmass[j];
776 double beta = gb / sqrt( 1 + gb * gb );
777 texp[j] = path / beta / velc;
778 }
779
780 m_qual_btof1[ii] = qual;
781 m_tof_btof1[ii] = tof;
782 }
783
784 if ( status->layer() == 2 )
785 { // layer2
786 double path = ( *iter_tof )->path(); // ?
787 double tof = ( *iter_tof )->tof();
788 double ph = ( *iter_tof )->ph();
789 double rhit = ( *iter_tof )->zrhit();
790 double qual = 0.0 + ( *iter_tof )->quality();
791 double cntr = 0.0 + ( *iter_tof )->tofID();
792 double texp[5];
793 for ( int j = 0; j < 5; j++ )
794 {
795 double gb = ptrk / xmass[j];
796 double beta = gb / sqrt( 1 + gb * gb );
797 texp[j] = path / beta / velc;
798 }
799
800 m_qual_btof2[ii] = qual;
801 m_tof_btof2[ii] = tof;
802 }
803 }
804 }
805 }
806 }
807
808 // Kinamatic Fit
810
811 double ecms = 3.097;
812 double chisq = 9999.;
813 m_4c_chi2 = 9999.;
814 m_4c_mks = 10.0;
815 m_4c_mkspi = 10.0;
816 m_4c_mksk = 10.0;
817 m_4c_mkpi = 10.0;
818 m_4c_ks_px = 10.0;
819 m_4c_ks_py = 10.0;
820 m_4c_ks_pz = 10.0;
821 m_4c_ks_p = 10.0;
822 m_4c_ks_cos = 10.0;
823
824 kmfit->init();
825 kmfit->AddTrack( 0, wpi );
826 kmfit->AddTrack( 1, wk );
827 kmfit->AddTrack( 2, wks );
828 kmfit->AddFourMomentum( 0, p_cms );
829 bool oksq = kmfit->Fit();
830 if ( oksq )
831 {
832 chisq = kmfit->chisq();
833
834 HepLorentzVector pk = kmfit->pfit( 1 );
835 HepLorentzVector pks = kmfit->pfit( 2 );
836 HepLorentzVector pkspi = kmfit->pfit( 0 ) + kmfit->pfit( 2 );
837 HepLorentzVector pksk = kmfit->pfit( 1 ) + kmfit->pfit( 2 );
838 HepLorentzVector pkpi = kmfit->pfit( 0 ) + kmfit->pfit( 1 );
839
840 pks.boost( u_cms );
841 pkspi.boost( u_cms );
842 pksk.boost( u_cms );
843 pkpi.boost( u_cms );
844
845 m_4c_chi2 = chisq;
846 m_4c_mks = pks.m();
847 m_4c_mkspi = pkspi.m();
848 m_4c_mksk = pksk.m();
849 m_4c_mkpi = pkpi.m();
850
851 m_4c_ks_px = pks.px();
852 m_4c_ks_py = pks.py();
853 m_4c_ks_pz = pks.pz();
854 m_4c_ks_p = ( pks.vect() ).mag();
855 m_4c_ks_cos = m_4c_ks_pz / m_4c_ks_p;
856 }
857
858 chisq = 9999.;
859 m_chi2_fs4c = 9999.;
860 m_mks_fs4c = 10.0;
861 m_mkspi_fs4c = 10.0;
862 m_mksk_fs4c = 10.0;
863 m_mkpi_fs4c = 10.0;
864
865 kmfit->init();
866 kmfit->AddTrack( 0, wpip );
867 kmfit->AddTrack( 1, wpim );
868 kmfit->AddTrack( 2, wpi );
869 kmfit->AddTrack( 3, wk );
870 kmfit->AddFourMomentum( 0, p_cms );
871 oksq = kmfit->Fit();
872 if ( oksq )
873 {
874 chisq = kmfit->chisq();
875
876 HepLorentzVector pks = kmfit->pfit( 0 ) + kmfit->pfit( 1 );
877 HepLorentzVector pkspi = pks + kmfit->pfit( 2 );
878 HepLorentzVector pksk = pks + kmfit->pfit( 3 );
879 HepLorentzVector pkpi = kmfit->pfit( 2 ) + kmfit->pfit( 3 );
880
881 pks.boost( u_cms );
882 pkspi.boost( u_cms );
883 pksk.boost( u_cms );
884 pkpi.boost( u_cms );
885
886 m_chi2_fs4c = chisq;
887 m_mks_fs4c = pks.m();
888 m_mkspi_fs4c = pkspi.m();
889 m_mksk_fs4c = pksk.m();
890 m_mkpi_fs4c = pkpi.m();
891 }
892
893 // finale selection
894 if ( chisq > 20 ) { return StatusCode::SUCCESS; } // 4C chi2
895 if ( m_vfit2_dl < 0.5 ) { return StatusCode::SUCCESS; } // Ks decay length
896 if ( fabs( m_4c_mks - mks0 ) > 0.01 ) { return StatusCode::SUCCESS; } // Ks mass
897 if ( m_4c_mkspi < 1.25 ) { return StatusCode::SUCCESS; } // Kspi mass
898
899 // DQA
900 TH1* h( 0 );
901 if ( m_thsvc->getHist( "/DQAHist/DQAKsKpi/hks_dl", h ).isSuccess() )
902 { h->Fill( m_vfit2_dl ); }
903 else { log << MSG::ERROR << "Couldn't retrieve hks_dl" << endmsg; }
904
905 if ( m_thsvc->getHist( "/DQAHist/DQAKsKpi/hks_m", h ).isSuccess() ) { h->Fill( m_4c_mks ); }
906 else { log << MSG::ERROR << "Couldn't retrieve hks_m" << endmsg; }
907
908 if ( m_thsvc->getHist( "/DQAHist/DQAKsKpi/hkspi_m", h ).isSuccess() )
909 { h->Fill( m_4c_mkspi ); }
910 else { log << MSG::ERROR << "Couldn't retrieve hkspi_m" << endmsg; }
911
912 if ( m_thsvc->getHist( "/DQAHist/DQAKsKpi/hks_p", h ).isSuccess() ) { h->Fill( m_4c_ks_p ); }
913 else { log << MSG::ERROR << "Couldn't retrieve hks_p" << endmsg; }
914
915 if ( m_thsvc->getHist( "/DQAHist/DQAKsKpi/hkpi_m", h ).isSuccess() )
916 { h->Fill( m_4c_mkpi ); }
917 else { log << MSG::ERROR << "Couldn't retrieve hkpi_m" << endmsg; }
918
919 ////////////////////////////////////////////////////////////
920 // DQA
921 // set tag and quality
922
923 // Pid: 1 - electron, 2 - muon, 3 - pion, 4 - kaon, 5 - proton
924 if ( m_npip == 2 && m_npim == 1 )
925 {
926 ( *( evtRecTrkCol->begin() + ipip[0] ) )->setPartId( 2 );
927 ( *( evtRecTrkCol->begin() + ipip[1] ) )->setPartId( 2 );
928 ( *( evtRecTrkCol->begin() + ipim[0] ) )->setPartId( 2 );
929 ( *( evtRecTrkCol->begin() + ikm[0] ) )->setPartId( 4 );
930 }
931 if ( m_npip == 1 && m_npim == 2 )
932 {
933 ( *( evtRecTrkCol->begin() + ipip[0] ) )->setPartId( 2 );
934 ( *( evtRecTrkCol->begin() + ipim[0] ) )->setPartId( 2 );
935 ( *( evtRecTrkCol->begin() + ipim[1] ) )->setPartId( 2 );
936 ( *( evtRecTrkCol->begin() + ikp[0] ) )->setPartId( 4 );
937 }
938 // Quality: defined by whether dE/dx or TOF is used to identify particle
939 // 0 - no dE/dx, no TOF (can be used for dE/dx and TOF calibration)
940 // 1 - only dE/dx (can be used for TOF calibration)
941 // 2 - only TOF (can be used for dE/dx calibration)
942 // 3 - Both dE/dx and TOF
943 if ( m_npip == 2 && m_npim == 1 )
944 {
945 ( *( evtRecTrkCol->begin() + ipip[0] ) )->setQuality( 1 );
946 ( *( evtRecTrkCol->begin() + ipip[1] ) )->setQuality( 1 );
947 ( *( evtRecTrkCol->begin() + ipim[0] ) )->setQuality( 1 );
948 ( *( evtRecTrkCol->begin() + ikm[0] ) )->setQuality( 1 );
949 }
950 if ( m_npip == 1 && m_npim == 2 )
951 {
952 ( *( evtRecTrkCol->begin() + ipip[0] ) )->setQuality( 1 );
953 ( *( evtRecTrkCol->begin() + ipim[0] ) )->setQuality( 1 );
954 ( *( evtRecTrkCol->begin() + ipim[1] ) )->setQuality( 1 );
955 ( *( evtRecTrkCol->begin() + ikp[0] ) )->setQuality( 1 );
956 }
957 // DQA
958 // Add the line below at the end of execute(), (before return)
959 //
960 setFilterPassed( true );
961 ////////////////////////////////////////////////////////////
962
963 m_tuple->write().ignore();
964
965 return StatusCode::SUCCESS;
966}
967
968// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
969StatusCode DQAKsKpi::finalize() {
970
971 MsgStream log( msgSvc(), name() );
972 log << MSG::INFO << "in finalize()" << endmsg;
973 return StatusCode::SUCCESS;
974}
DECLARE_COMPONENT(BesBdkRc)
HepGeom::Point3D< double > HepPoint3D
const Hep3Vector u_cms
const double mks0
const HepLorentzVector p_cms(0.034067, 0.0, 0.0, 3.097)
const HepLorentzVector p_cms(0.034067, 0.0, 0.0, 3.097)
EvtRecTrackCol::iterator EvtRecTrackIterator
double mpi
std::vector< HepLorentzVector > Vp4
Definition Gam4pikp.cxx:38
const double xmass[5]
Definition Gam4pikp.cxx:35
const double velc
Definition Gam4pikp.cxx:36
const double mk
Definition Gam4pikp.cxx:33
std::vector< int > Vint
Definition Gam4pikp.cxx:37
INTupleSvc * ntupleSvc()
IMessageSvc * msgSvc()
DQAKsKpi(const std::string &name, ISvcLocator *pSvcLocator)
Definition DQAKsKpi.cxx:56
StatusCode finalize()
Definition DQAKsKpi.cxx:969
StatusCode execute()
Definition DQAKsKpi.cxx:231
StatusCode initialize()
Definition DQAKsKpi.cxx:70
const HepSymMatrix err() const
const HepVector helix() const
......
virtual bool isVertexValid()=0
virtual double * SigmaPrimaryVertex()=0
virtual double * PrimaryVertex()=0
static KinematicFit * instance()
void AddFourMomentum(int number, HepLorentzVector p4)
double chiTof2(int n) const
static ParticleID * instance()
bool IsPidInfoValid() const
double chiTof1(int n) const
void calculate()
void init()
double chiDedx(int n) const
void setPrimaryVertex(const VertexParameter vpar)
static SecondVertexFit * instance()
void setStatus(unsigned int status)
void AddTrack(const int number, const double mass, const RecMdcTrack *trk)
Definition TrackPool.cxx:21
const HepPoint3D & pivot(void) const
returns pivot position.
const HepVector & a(void) const
returns helix parameters.
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 ecms
Definition inclkstar.cxx:26