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