BOSS 7.1.0
BESIII Offline Software System
Loading...
Searching...
No Matches
Rhopi.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"
8#include "GaudiKernel/Bootstrap.h"
9#include "GaudiKernel/ISvcLocator.h"
10
12#include "EventModel/Event.h"
13
18
19
20
21#include "TMath.h"
22#include "GaudiKernel/INTupleSvc.h"
23#include "GaudiKernel/NTuple.h"
24#include "GaudiKernel/Bootstrap.h"
25#include "GaudiKernel/IHistogramSvc.h"
26#include "CLHEP/Vector/ThreeVector.h"
27#include "CLHEP/Vector/LorentzVector.h"
28#include "CLHEP/Vector/TwoVector.h"
29using CLHEP::Hep3Vector;
30using CLHEP::Hep2Vector;
31using CLHEP::HepLorentzVector;
32#include "CLHEP/Geometry/Point3D.h"
33#ifndef ENABLE_BACKWARDS_COMPATIBILITY
35#endif
36#include "RhopiAlg/Rhopi.h"
37
38//#include "VertexFit/KinematicFit.h"
40#include "VertexFit/VertexFit.h"
41#include "VertexFit/Helix.h"
43
44#include <vector>
45//const double twopi = 6.2831853;
46//const double pi = 3.1415927;
47const double mpi = 0.13957;
48const double xmass[5] = {0.000511, 0.105658, 0.139570,0.493677, 0.938272};
49//const double velc = 29.9792458; tof_path unit in cm.
50const double velc = 299.792458; // tof path unit in mm
51typedef std::vector<int> Vint;
52typedef std::vector<HepLorentzVector> Vp4;
53
55
56
57DECLARE_COMPONENT(Rhopi)
58
59Rhopi::Rhopi(const std::string& name, ISvcLocator* pSvcLocator) :
60 Algorithm(name, pSvcLocator) {
61
62 //Declare the properties
63 declareProperty("Vr0cut", m_vr0cut=1.0);
64 declareProperty("Vz0cut", m_vz0cut=5.0);
65 declareProperty("EnergyThreshold", m_energyThreshold=0.05);
66 declareProperty("GammaPhiCut", m_gammaPhiCut=20.0);
67 declareProperty("GammaThetaCut", m_gammaThetaCut=20.0);
68 declareProperty("GammaAngleCut", m_gammaAngleCut=10.0);
69 declareProperty("Test4C", m_test4C = 1);
70 declareProperty("Test5C", m_test5C = 1);
71 declareProperty("CheckDedx", m_checkDedx = 1);
72 declareProperty("CheckTof", m_checkTof = 1);
73}
74
75// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
76StatusCode Rhopi::initialize(){
77 MsgStream log(msgSvc(), name());
78
79 log << MSG::INFO << "in initialize()" << endmsg;
80
81 StatusCode status;
82 NTuplePtr nt1(ntupleSvc(), "FILE1/vxyz");
83 if ( nt1 ) m_tuple1 = nt1;
84 else {
85 m_tuple1 = ntupleSvc()->book ("FILE1/vxyz", CLID_ColumnWiseTuple, "ks N-Tuple example");
86 if ( m_tuple1 ) {
87 status = m_tuple1->addItem ("vx0", m_vx0);
88 status = m_tuple1->addItem ("vy0", m_vy0);
89 status = m_tuple1->addItem ("vz0", m_vz0);
90 status = m_tuple1->addItem ("vr0", m_vr0);
91 status = m_tuple1->addItem ("rvxy0", m_rvxy0);
92 status = m_tuple1->addItem ("rvz0", m_rvz0);
93 status = m_tuple1->addItem ("rvphi0", m_rvphi0);
94 }
95 else {
96 log << MSG::ERROR << " Cannot book N-tuple:" << long(m_tuple1) << endmsg;
97 return StatusCode::FAILURE;
98 }
99 }
100
101 NTuplePtr nt2(ntupleSvc(), "FILE1/photon");
102 if ( nt2 ) m_tuple2 = nt2;
103 else {
104 m_tuple2 = ntupleSvc()->book ("FILE1/photon", CLID_ColumnWiseTuple, "ks N-Tuple example");
105 if ( m_tuple2 ) {
106 status = m_tuple2->addItem ("dthe", m_dthe);
107 status = m_tuple2->addItem ("dphi", m_dphi);
108 status = m_tuple2->addItem ("dang", m_dang);
109 status = m_tuple2->addItem ("eraw", m_eraw);
110 }
111 else {
112 log << MSG::ERROR << " Cannot book N-tuple:" << long(m_tuple2) << endmsg;
113 return StatusCode::FAILURE;
114 }
115 }
116
117
118 NTuplePtr nt3(ntupleSvc(), "FILE1/etot");
119 if ( nt3 ) m_tuple3 = nt3;
120 else {
121 m_tuple3 = ntupleSvc()->book ("FILE1/etot", CLID_ColumnWiseTuple, "ks N-Tuple example");
122 if ( m_tuple3 ) {
123 status = m_tuple3->addItem ("m2gg", m_m2gg);
124 status = m_tuple3->addItem ("etot", m_etot);
125 }
126 else {
127 log << MSG::ERROR << " Cannot book N-tuple:" << long(m_tuple3) << endmsg;
128 return StatusCode::FAILURE;
129 }
130 }
131 if(m_test4C==1) {
132 NTuplePtr nt4(ntupleSvc(), "FILE1/fit4c");
133 if ( nt4 ) m_tuple4 = nt4;
134 else {
135 m_tuple4 = ntupleSvc()->book ("FILE1/fit4c", CLID_ColumnWiseTuple, "ks N-Tuple example");
136 if ( m_tuple4 ) {
137 status = m_tuple4->addItem ("chi2", m_chi1);
138 status = m_tuple4->addItem ("mpi0", m_mpi0);
139 }
140 else {
141 log << MSG::ERROR << " Cannot book N-tuple:" << long(m_tuple4) << endmsg;
142 return StatusCode::FAILURE;
143 }
144 }
145 } // test 4C
146
147
148 if(m_test5C==1) {
149 NTuplePtr nt5(ntupleSvc(), "FILE1/fit5c");
150 if ( nt5 ) m_tuple5 = nt5;
151 else {
152 m_tuple5 = ntupleSvc()->book ("FILE1/fit5c", CLID_ColumnWiseTuple, "ks N-Tuple example");
153 if ( m_tuple5 ) {
154 status = m_tuple5->addItem ("chi2", m_chi2);
155 status = m_tuple5->addItem ("mrh0", m_mrh0);
156 status = m_tuple5->addItem ("mrhp", m_mrhp);
157 status = m_tuple5->addItem ("mrhm", m_mrhm);
158 }
159 else {
160 log << MSG::ERROR << " Cannot book N-tuple:" << long(m_tuple5) << endmsg;
161 return StatusCode::FAILURE;
162 }
163 }
164
165 NTuplePtr nt6(ntupleSvc(), "FILE1/geff");
166 if ( nt6 ) m_tuple6 = nt6;
167 else {
168 m_tuple6 = ntupleSvc()->book ("FILE1/geff", CLID_ColumnWiseTuple, "ks N-Tuple example");
169 if ( m_tuple6 ) {
170 status = m_tuple6->addItem ("fcos", m_fcos);
171 status = m_tuple6->addItem ("elow", m_elow);
172 }
173 else {
174 log << MSG::ERROR << " Cannot book N-tuple:" << long(m_tuple6) << endmsg;
175 return StatusCode::FAILURE;
176 }
177 }
178 } // test 5c
179
180 if(m_checkDedx == 1) {
181 NTuplePtr nt7(ntupleSvc(), "FILE1/dedx");
182 if ( nt7 ) m_tuple7 = nt7;
183 else {
184 m_tuple7 = ntupleSvc()->book ("FILE1/dedx", CLID_ColumnWiseTuple, "ks N-Tuple example");
185 if ( m_tuple7 ) {
186 status = m_tuple7->addItem ("ptrk", m_ptrk);
187 status = m_tuple7->addItem ("chie", m_chie);
188 status = m_tuple7->addItem ("chimu", m_chimu);
189 status = m_tuple7->addItem ("chipi", m_chipi);
190 status = m_tuple7->addItem ("chik", m_chik);
191 status = m_tuple7->addItem ("chip", m_chip);
192 status = m_tuple7->addItem ("probPH", m_probPH);
193 status = m_tuple7->addItem ("normPH", m_normPH);
194 status = m_tuple7->addItem ("ghit", m_ghit);
195 status = m_tuple7->addItem ("thit", m_thit);
196 }
197 else {
198 log << MSG::ERROR << " Cannot book N-tuple:" << long(m_tuple7) << endmsg;
199 return StatusCode::FAILURE;
200 }
201 }
202 } // check dE/dx
203
204 if(m_checkTof == 1) {
205 NTuplePtr nt8(ntupleSvc(), "FILE1/tofe");
206 if ( nt8 ) m_tuple8 = nt8;
207 else {
208 m_tuple8 = ntupleSvc()->book ("FILE1/tofe",CLID_ColumnWiseTuple, "ks N-Tuple example");
209 if ( m_tuple8 ) {
210 status = m_tuple8->addItem ("ptrk", m_ptot_etof);
211 status = m_tuple8->addItem ("cntr", m_cntr_etof);
212 status = m_tuple8->addItem ("ph", m_ph_etof);
213 status = m_tuple8->addItem ("rhit", m_rhit_etof);
214 status = m_tuple8->addItem ("qual", m_qual_etof);
215 status = m_tuple8->addItem ("te", m_te_etof);
216 status = m_tuple8->addItem ("tmu", m_tmu_etof);
217 status = m_tuple8->addItem ("tpi", m_tpi_etof);
218 status = m_tuple8->addItem ("tk", m_tk_etof);
219 status = m_tuple8->addItem ("tp", m_tp_etof);
220 }
221 else {
222 log << MSG::ERROR << " Cannot book N-tuple:" << long(m_tuple8) << endmsg;
223 return StatusCode::FAILURE;
224 }
225 }
226 } // check Tof:endcap
227
228
229
230 if(m_checkTof == 1) {
231 NTuplePtr nt9(ntupleSvc(), "FILE1/tof1");
232 if ( nt9 ) m_tuple9 = nt9;
233 else {
234 m_tuple9 = ntupleSvc()->book ("FILE1/tof1", CLID_ColumnWiseTuple, "ks N-Tuple example");
235 if ( m_tuple9 ) {
236 status = m_tuple9->addItem ("ptrk", m_ptot_btof1);
237 status = m_tuple9->addItem ("cntr", m_cntr_btof1);
238 status = m_tuple9->addItem ("ph", m_ph_btof1);
239 status = m_tuple9->addItem ("zhit", m_zhit_btof1);
240 status = m_tuple9->addItem ("qual", m_qual_btof1);
241 status = m_tuple9->addItem ("te", m_te_btof1);
242 status = m_tuple9->addItem ("tmu", m_tmu_btof1);
243 status = m_tuple9->addItem ("tpi", m_tpi_btof1);
244 status = m_tuple9->addItem ("tk", m_tk_btof1);
245 status = m_tuple9->addItem ("tp", m_tp_btof1);
246 }
247 else {
248 log << MSG::ERROR << " Cannot book N-tuple:" << long(m_tuple9) << endmsg;
249 return StatusCode::FAILURE;
250 }
251 }
252 } // check Tof:barrel inner Tof
253
254
255 if(m_checkTof == 1) {
256 NTuplePtr nt10(ntupleSvc(), "FILE1/tof2");
257 if ( nt10 ) m_tuple10 = nt10;
258 else {
259 m_tuple10 = ntupleSvc()->book ("FILE1/tof2", CLID_ColumnWiseTuple, "ks N-Tuple example");
260 if ( m_tuple10 ) {
261 status = m_tuple10->addItem ("ptrk", m_ptot_btof2);
262 status = m_tuple10->addItem ("cntr", m_cntr_btof2);
263 status = m_tuple10->addItem ("ph", m_ph_btof2);
264 status = m_tuple10->addItem ("zhit", m_zhit_btof2);
265 status = m_tuple10->addItem ("qual", m_qual_btof2);
266 status = m_tuple10->addItem ("te", m_te_btof2);
267 status = m_tuple10->addItem ("tmu", m_tmu_btof2);
268 status = m_tuple10->addItem ("tpi", m_tpi_btof2);
269 status = m_tuple10->addItem ("tk", m_tk_btof2);
270 status = m_tuple10->addItem ("tp", m_tp_btof2);
271 }
272 else {
273 log << MSG::ERROR << " Cannot book N-tuple:" << long(m_tuple10) << endmsg;
274 return StatusCode::FAILURE;
275 }
276 }
277 } // check Tof:barrel outter Tof
278
279
280 NTuplePtr nt11(ntupleSvc(), "FILE1/pid");
281 if ( nt11 ) m_tuple11 = nt11;
282 else {
283 m_tuple11 = ntupleSvc()->book ("FILE1/pid", CLID_ColumnWiseTuple, "ks N-Tuple example");
284 if ( m_tuple11 ) {
285 status = m_tuple11->addItem ("ptrk", m_ptrk_pid);
286 status = m_tuple11->addItem ("cost", m_cost_pid);
287 status = m_tuple11->addItem ("dedx", m_dedx_pid);
288 status = m_tuple11->addItem ("tof1", m_tof1_pid);
289 status = m_tuple11->addItem ("tof2", m_tof2_pid);
290 status = m_tuple11->addItem ("prob", m_prob_pid);
291 }
292 else {
293 log << MSG::ERROR << " Cannot book N-tuple:" << long(m_tuple11) << endmsg;
294 return StatusCode::FAILURE;
295 }
296 }
297
298
299 //
300 //--------end of book--------
301 //
302
303 log << MSG::INFO << "successfully return from initialize()" <<endmsg;
304 return StatusCode::SUCCESS;
305
306}
307
308// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
309StatusCode Rhopi::execute() {
310
311 //std::cout << "execute()" << std::endl;
312
313 MsgStream log(msgSvc(), name());
314 log << MSG::INFO << "in execute()" << endreq;
315
316 SmartDataPtr<Event::EventHeader> eventHeader(eventSvc(),"/Event/EventHeader");
317 int runNo=eventHeader->runNumber();
318 int event=eventHeader->eventNumber();
319 log << MSG::DEBUG <<"run, evtnum = "
320 << runNo << " , "
321 << event <<endreq;
322 cout<<"event "<<event<<endl;
323 Ncut0++;
324
325 SmartDataPtr<EvtRecEvent> evtRecEvent(eventSvc(), EventModel::EvtRec::EvtRecEvent);
326 // log << MSG::INFO << "get event tag OK" << endreq;
327 log << MSG::DEBUG <<"ncharg, nneu, tottks = "
328 << evtRecEvent->totalCharged() << " , "
329 << evtRecEvent->totalNeutral() << " , "
330 << evtRecEvent->totalTracks() <<endreq;
331
332 SmartDataPtr<EvtRecTrackCol> evtRecTrkCol(eventSvc(), EventModel::EvtRec::EvtRecTrackCol);
333 //
334 // check x0, y0, z0, r0
335 // suggest cut: |z0|<5 && r0<1
336 //
337 Vint iGood, ipip, ipim;
338 iGood.clear();
339 ipip.clear();
340 ipim.clear();
341 Vp4 ppip, ppim;
342 ppip.clear();
343 ppim.clear();
344
345 int nCharge = 0;
346
347 Hep3Vector xorigin(0,0,0);
348
349
350 /*VertexDbSvc* vtxsvc;
351 Gaudi::svcLocator()->service("VertexDbSvc", vtxsvc);
352 if(vtxsvc->isVertexValid()){
353 double* dbv = vtxsvc->PrimaryVertex();
354 double* vv = vtxsvc->SigmaPrimaryVertex();
355 xorigin.setX(dbv[0]);
356 xorigin.setY(dbv[1]);
357 xorigin.setZ(dbv[2]);
358 }*/
359 SmartIF<IVertexDbSvc> m_vtxSvc;
360 m_vtxSvc = serviceLocator()->service("VertexDbSvc");
361 if(m_vtxSvc->isVertexValid()){
362 //if(m_vtxSvc){
363 double* dbv = m_vtxSvc->PrimaryVertex();
364 double* vv = m_vtxSvc->SigmaPrimaryVertex();
365 xorigin.setX(dbv[0]);
366 xorigin.setY(dbv[1]);
367 xorigin.setZ(dbv[2]);
368 }
369
370 for(int i = 0; i < evtRecEvent->totalCharged(); i++){
371 EvtRecTrackIterator itTrk=evtRecTrkCol->begin() + i;
372 if(!(*itTrk)->isMdcTrackValid()) continue;
373 RecMdcTrack *mdcTrk = (*itTrk)->mdcTrack();
374 double pch=mdcTrk->p();
375 double x0=mdcTrk->x();
376 double y0=mdcTrk->y();
377 double z0=mdcTrk->z();
378 double theta0=mdcTrk->theta();
379 double phi0=mdcTrk->helix(1);
380 double xv=xorigin.x();
381 double yv=xorigin.y();
382 double Rxy=(x0-xv)*cos(phi0)+(y0-yv)*sin(phi0);
383 m_vx0 = x0;
384 m_vy0 = y0;
385 m_vz0 = z0;
386 m_vr0 = Rxy;
387
388 HepVector a = mdcTrk->helix();
389 HepSymMatrix Ea = mdcTrk->err();
390 HepPoint3D point0(0.,0.,0.); // the initial point for MDC recosntruction
391 HepPoint3D IP(xorigin[0],xorigin[1],xorigin[2]);
392 VFHelix helixip(point0,a,Ea);
393 helixip.pivot(IP);
394 HepVector vecipa = helixip.a();
395 double Rvxy0=fabs(vecipa[0]); //the nearest distance to IP in xy plane
396 double Rvz0=vecipa[3]; //the nearest distance to IP in z direction
397 double Rvphi0=vecipa[1];
398 m_rvxy0=Rvxy0;
399 m_rvz0=Rvz0;
400 m_rvphi0=Rvphi0;
401
402 m_tuple1->write();
403// if(fabs(z0) >= m_vz0cut) continue;
404// if(fabs(Rxy) >= m_vr0cut) continue;
405
406 if(fabs(Rvz0) >= 10.0) continue;
407 if(fabs(Rvxy0) >= 1.0) continue;
408 if(fabs(cos(theta0)) >= 0.93) continue;
409
410 iGood.push_back(i);
411 nCharge += mdcTrk->charge();
412 }
413
414 //
415 // Finish Good Charged Track Selection
416 //
417 int nGood = iGood.size();
418 log << MSG::DEBUG << "ngood, totcharge = " << nGood << " , " << nCharge << endreq;
419 if((nGood != 2)||(nCharge!=0)){
420 return StatusCode::SUCCESS;
421 }
422 Ncut1++;
423
424 Vint iGam;
425 iGam.clear();
426 for(int i = evtRecEvent->totalCharged(); i< evtRecEvent->totalTracks(); i++) {
427 EvtRecTrackIterator itTrk=evtRecTrkCol->begin() + i;
428 if(!(*itTrk)->isEmcShowerValid()) continue;
429 RecEmcShower *emcTrk = (*itTrk)->emcShower();
430 Hep3Vector emcpos(emcTrk->x(), emcTrk->y(), emcTrk->z());
431 // find the nearest charged track
432 double dthe = 200.;
433 double dphi = 200.;
434 double dang = 200.;
435 for(int j = 0; j < evtRecEvent->totalCharged(); j++) {
436 EvtRecTrackIterator jtTrk = evtRecTrkCol->begin() + j;
437 if(!(*jtTrk)->isExtTrackValid()) continue;
438 RecExtTrack *extTrk = (*jtTrk)->extTrack();
439 if(extTrk->emcVolumeNumber() == -1) continue;
440 Hep3Vector extpos = extTrk->emcPosition();
441 // double ctht = extpos.cosTheta(emcpos);
442 double angd = extpos.angle(emcpos);
443 double thed = extpos.theta() - emcpos.theta();
444 double phid = extpos.deltaPhi(emcpos);
445 thed = fmod(thed+CLHEP::twopi+CLHEP::twopi+pi, CLHEP::twopi) - CLHEP::pi;
446 phid = fmod(phid+CLHEP::twopi+CLHEP::twopi+pi, CLHEP::twopi) - CLHEP::pi;
447 if(angd < dang){
448 dang = angd;
449 dthe = thed;
450 dphi = phid;
451 }
452 }
453 if(dang>=200) continue;
454 double eraw = emcTrk->energy();
455 double theta1 = emcTrk->theta();
456 double time = emcTrk->time();
457 dthe = dthe * 180 / (CLHEP::pi);
458 dphi = dphi * 180 / (CLHEP::pi);
459 dang = dang * 180 / (CLHEP::pi);
460 m_dthe = dthe;
461 m_dphi = dphi;
462 m_dang = dang;
463 m_eraw = eraw;
464 m_tuple2->write();
465
466 if (fabs(cos(theta1))<0.80)
467 {if(eraw <= (m_energyThreshold/2)) continue;}
468 else if (fabs(cos(theta1))>0.86&&fabs(cos(theta1))<0.92)
469 {if(eraw <= m_energyThreshold) continue;}
470 else continue;
471
472
473// if((fabs(dthe) < m_gammaThetaCut) && (fabs(dphi)<m_gammaPhiCut) ) continue;
474 if(fabs(dang) < m_gammaAngleCut) continue;
475
476 if(time < 0 || time > 14) continue;
477
478 //
479 // good photon cut will be set here
480 //
481 iGam.push_back(i);
482 }
483
484 //
485 // Finish Good Photon Selection
486 //
487 int nGam = iGam.size();
488
489 log << MSG::DEBUG << "num Good Photon " << nGam << " , " <<evtRecEvent->totalNeutral()<<endreq;
490 if(nGam<2){
491 return StatusCode::SUCCESS;
492 }
493 Ncut2++;
494
495
496
497 //
498 //
499 // check dedx infomation
500 //
501 //
502
503 if(m_checkDedx == 1) {
504 for(int i = 0; i < nGood; i++) {
505 EvtRecTrackIterator itTrk = evtRecTrkCol->begin() + iGood[i];
506 if(!(*itTrk)->isMdcTrackValid()) continue;
507 if(!(*itTrk)->isMdcDedxValid())continue;
508 RecMdcTrack* mdcTrk = (*itTrk)->mdcTrack();
509 RecMdcDedx* dedxTrk = (*itTrk)->mdcDedx();
510 m_ptrk = mdcTrk->p();
511
512 m_chie = dedxTrk->chiE();
513 m_chimu = dedxTrk->chiMu();
514 m_chipi = dedxTrk->chiPi();
515 m_chik = dedxTrk->chiK();
516 m_chip = dedxTrk->chiP();
517 m_ghit = dedxTrk->numGoodHits();
518 m_thit = dedxTrk->numTotalHits();
519 m_probPH = dedxTrk->probPH();
520 m_normPH = dedxTrk->normPH();
521 m_tuple7->write();
522 }
523 }
524
525 //
526 // check TOF infomation
527 //
528
529
530 if(m_checkTof == 1) {
531 for(int i = 0; i < nGood; i++) {
532 EvtRecTrackIterator itTrk = evtRecTrkCol->begin() + iGood[i];
533 if(!(*itTrk)->isMdcTrackValid()) continue;
534 if(!(*itTrk)->isTofTrackValid()) continue;
535
536 RecMdcTrack * mdcTrk = (*itTrk)->mdcTrack();
537 SmartRefVector<RecTofTrack> tofTrkCol = (*itTrk)->tofTrack();
538
539 double ptrk = mdcTrk->p();
540
541 SmartRefVector<RecTofTrack>::iterator iter_tof = tofTrkCol.begin();
542 for(;iter_tof != tofTrkCol.end(); iter_tof++ ) {
543 TofHitStatus *status = new TofHitStatus;
544 status->setStatus((*iter_tof)->status());
545 if(!(status->is_barrel())){//endcap
546 if( !(status->is_counter()) ) continue; // ?
547 if( status->layer()!=0 ) continue;//layer1
548 double path=(*iter_tof)->path(); // ?
549 double tof = (*iter_tof)->tof();
550 double ph = (*iter_tof)->ph();
551 double rhit = (*iter_tof)->zrhit();
552 double qual = 0.0 + (*iter_tof)->quality();
553 double cntr = 0.0 + (*iter_tof)->tofID();
554 double texp[5];
555 for(int j = 0; j < 5; j++) {
556 double gb = ptrk/xmass[j];
557 double beta = gb/sqrt(1+gb*gb);
558 texp[j] = 10 * path /beta/velc;
559 }
560 m_cntr_etof = cntr;
561 m_ptot_etof = ptrk;
562 m_ph_etof = ph;
563 m_rhit_etof = rhit;
564 m_qual_etof = qual;
565 m_te_etof = tof - texp[0];
566 m_tmu_etof = tof - texp[1];
567 m_tpi_etof = tof - texp[2];
568 m_tk_etof = tof - texp[3];
569 m_tp_etof = tof - texp[4];
570 m_tuple8->write();
571 }
572 else {//barrel
573 if( !(status->is_counter()) ) continue; // ?
574 if(status->layer()==1){ //layer1
575 double path=(*iter_tof)->path(); // ?
576 double tof = (*iter_tof)->tof();
577 double ph = (*iter_tof)->ph();
578 double rhit = (*iter_tof)->zrhit();
579 double qual = 0.0 + (*iter_tof)->quality();
580 double cntr = 0.0 + (*iter_tof)->tofID();
581 double texp[5];
582 for(int j = 0; j < 5; j++) {
583 double gb = ptrk/xmass[j];
584 double beta = gb/sqrt(1+gb*gb);
585 texp[j] = 10 * path /beta/velc;
586 }
587
588 m_cntr_btof1 = cntr;
589 m_ptot_btof1 = ptrk;
590 m_ph_btof1 = ph;
591 m_zhit_btof1 = rhit;
592 m_qual_btof1 = qual;
593 m_te_btof1 = tof - texp[0];
594 m_tmu_btof1 = tof - texp[1];
595 m_tpi_btof1 = tof - texp[2];
596 m_tk_btof1 = tof - texp[3];
597 m_tp_btof1 = tof - texp[4];
598 m_tuple9->write();
599 }
600
601 if(status->layer()==2){//layer2
602 double path=(*iter_tof)->path(); // ?
603 double tof = (*iter_tof)->tof();
604 double ph = (*iter_tof)->ph();
605 double rhit = (*iter_tof)->zrhit();
606 double qual = 0.0 + (*iter_tof)->quality();
607 double cntr = 0.0 + (*iter_tof)->tofID();
608 double texp[5];
609 for(int j = 0; j < 5; j++) {
610 double gb = ptrk/xmass[j];
611 double beta = gb/sqrt(1+gb*gb);
612 texp[j] = 10 * path /beta/velc;
613 }
614
615 m_cntr_btof2 = cntr;
616 m_ptot_btof2 = ptrk;
617 m_ph_btof2 = ph;
618 m_zhit_btof2 = rhit;
619 m_qual_btof2 = qual;
620 m_te_btof2 = tof - texp[0];
621 m_tmu_btof2 = tof - texp[1];
622 m_tpi_btof2 = tof - texp[2];
623 m_tk_btof2 = tof - texp[3];
624 m_tp_btof2 = tof - texp[4];
625 m_tuple10->write();
626 }
627 }
628
629 delete status;
630 }
631 } // loop all charged track
632 } // check tof
633
634
635 //
636 // Assign 4-momentum to each photon
637 //
638
639 Vp4 pGam;
640 pGam.clear();
641 for(int i = 0; i < nGam; i++) {
642 EvtRecTrackIterator itTrk = evtRecTrkCol->begin() + iGam[i];
643 RecEmcShower* emcTrk = (*itTrk)->emcShower();
644 double eraw = emcTrk->energy();
645 double phi = emcTrk->phi();
646 double the = emcTrk->theta();
647 HepLorentzVector ptrk;
648 ptrk.setPx(eraw*sin(the)*cos(phi));
649 ptrk.setPy(eraw*sin(the)*sin(phi));
650 ptrk.setPz(eraw*cos(the));
651 ptrk.setE(eraw);
652
653// ptrk = ptrk.boost(-0.011,0,0);// boost to cms
654
655 pGam.push_back(ptrk);
656 }
657 //cout<<"before pid"<<endl;
658 //
659 // Assign 4-momentum to each charged track
660 //
662 for(int i = 0; i < nGood; i++) {
663 EvtRecTrackIterator itTrk = evtRecTrkCol->begin() + iGood[i];
664 // if(pid) delete pid;
665 pid->init();
666 pid->setMethod(pid->methodProbability());
667// pid->setMethod(pid->methodLikelihood()); //for Likelihood Method
668
669 pid->setChiMinCut(4);
670 pid->setRecTrack(*itTrk);
671// pid->usePidSys(pid->useDedx() | pid->useTof1() | pid->useTof2() | pid->useTofE()); // for BOSS version < 7.0.0
672 pid->usePidSys(pid->useDedx() | pid->useTofCorr()); // use PID sub-system
673 pid->identify(pid->onlyPionKaonProton()); // seperater Pion/Kaon/Proton
674 // pid->identify(pid->onlyPion());
675 // pid->identify(pid->onlyKaon());
676 pid->calculate();
677 if(!(pid->IsPidInfoValid())) continue;
678 RecMdcTrack* mdcTrk = (*itTrk)->mdcTrack();
679 m_ptrk_pid = mdcTrk->p();
680 m_cost_pid = cos(mdcTrk->theta());
681 m_dedx_pid = pid->chiDedx(2);
682 m_tof1_pid = pid->chiTof1(2);
683 m_tof2_pid = pid->chiTof2(2);
684 m_prob_pid = pid->probPion();
685 m_tuple11->write();
686
687 if((pid->probPion() < pid->probProton()) || (pid->probPion() < pid->probKaon())) continue;
688// if(pid->probPion() < 0.001) continue;
689// if(pid->pdf(2)<pid->pdf(3)) continue; // for Likelihood Method(0=electron 1=muon 2=pion 3=kaon 4=proton)
690
691 RecMdcKalTrack* mdcKalTrk = (*itTrk)->mdcKalTrack();//After ParticleID, use RecMdcKalTrack substitute RecMdcTrack
692 RecMdcKalTrack::setPidType (RecMdcKalTrack::pion);//PID can set to electron, muon, pion, kaon and proton;The default setting is pion
693
694 if(mdcKalTrk->charge() >0) {
695 ipip.push_back(iGood[i]);
696 HepLorentzVector ptrk;
697 ptrk.setPx(mdcKalTrk->px());
698 ptrk.setPy(mdcKalTrk->py());
699 ptrk.setPz(mdcKalTrk->pz());
700 double p3 = ptrk.mag();
701 ptrk.setE(sqrt(p3*p3+mpi*mpi));
702
703// ptrk = ptrk.boost(-0.011,0,0);//boost to cms
704
705 ppip.push_back(ptrk);
706 } else {
707 ipim.push_back(iGood[i]);
708 HepLorentzVector ptrk;
709 ptrk.setPx(mdcKalTrk->px());
710 ptrk.setPy(mdcKalTrk->py());
711 ptrk.setPz(mdcKalTrk->pz());
712 double p3 = ptrk.mag();
713 ptrk.setE(sqrt(p3*p3+mpi*mpi));
714
715// ptrk = ptrk.boost(-0.011,0,0);//boost to cms
716
717 ppim.push_back(ptrk);
718 }
719 }
720
721/*
722 for(int i = 0; i < nGood; i++) {//for rhopi without PID
723 EvtRecTrackIterator itTrk = evtRecTrkCol->begin() + iGood[i];
724 RecMdcTrack* mdcTrk = (*itTrk)->mdcTrack();
725 if(mdcTrk->charge() >0) {
726 ipip.push_back(iGood[i]);
727 HepLorentzVector ptrk;
728 ptrk.setPx(mdcTrk->px());
729 ptrk.setPy(mdcTrk->py());
730 ptrk.setPz(mdcTrk->pz());
731 double p3 = ptrk.mag();
732 ptrk.setE(sqrt(p3*p3+mpi*mpi));
733 ppip.push_back(ptrk);
734 } else {
735 ipim.push_back(iGood[i]);
736 HepLorentzVector ptrk;
737 ptrk.setPx(mdcTrk->px());
738 ptrk.setPy(mdcTrk->py());
739 ptrk.setPz(mdcTrk->pz());
740 double p3 = ptrk.mag();
741 ptrk.setE(sqrt(p3*p3+mpi*mpi));
742 ppim.push_back(ptrk);
743 }
744 }// without PID
745*/
746
747 int npip = ipip.size();
748 int npim = ipim.size();
749 if(npip*npim != 1) return SUCCESS;
750
751 Ncut3++;
752
753
754 //
755 // Loop each gamma pair, check ppi0 and pTot
756 //
757
758 HepLorentzVector pTot;
759 for(int i = 0; i < nGam - 1; i++){
760 for(int j = i+1; j < nGam; j++) {
761 HepLorentzVector p2g = pGam[i] + pGam[j];
762 pTot = ppip[0] + ppim[0];
763 pTot += p2g;
764 m_m2gg = p2g.m();
765 m_etot = pTot.e();
766 m_tuple3 -> write();
767
768 }
769 }
770
771
772 RecMdcKalTrack *pipTrk = (*(evtRecTrkCol->begin()+ipip[0]))->mdcKalTrack();
773 RecMdcKalTrack *pimTrk = (*(evtRecTrkCol->begin()+ipim[0]))->mdcKalTrack();
774
775 WTrackParameter wvpipTrk, wvpimTrk;
776 wvpipTrk = WTrackParameter(mpi, pipTrk->getZHelix(), pipTrk->getZError());
777 wvpimTrk = WTrackParameter(mpi, pimTrk->getZHelix(), pimTrk->getZError());
778
779/* Default is pion, for other particles:
780 wvppTrk = WTrackParameter(mp, pipTrk->getZHelixP(), pipTrk->getZErrorP());//proton
781 wvmupTrk = WTrackParameter(mmu, pipTrk->getZHelixMu(), pipTrk->getZErrorMu());//muon
782 wvepTrk = WTrackParameter(me, pipTrk->getZHelixE(), pipTrk->getZErrorE());//electron
783 wvkpTrk = WTrackParameter(mk, pipTrk->getZHelixK(), pipTrk->getZErrorK());//kaon
784*/
785 //
786 // Test vertex fit
787 //
788
789 HepPoint3D vx(0., 0., 0.);
790 HepSymMatrix Evx(3, 0);
791 double bx = 1E+6;
792 double by = 1E+6;
793 double bz = 1E+6;
794 Evx[0][0] = bx*bx;
795 Evx[1][1] = by*by;
796 Evx[2][2] = bz*bz;
797
798 VertexParameter vxpar;
799 vxpar.setVx(vx);
800 vxpar.setEvx(Evx);
801
802 VertexFit* vtxfit = VertexFit::instance();
803 vtxfit->init();
804 vtxfit->AddTrack(0, wvpipTrk);
805 vtxfit->AddTrack(1, wvpimTrk);
806 vtxfit->AddVertex(0, vxpar,0, 1);
807 if(!vtxfit->Fit(0)) return SUCCESS;
808 vtxfit->Swim(0);
809
810 WTrackParameter wpip = vtxfit->wtrk(0);
811 WTrackParameter wpim = vtxfit->wtrk(1);
812
813 //KinematicFit * kmfit = KinematicFit::instance();
815
816 //
817 // Apply Kinematic 4C fit
818 //
819 //cout<<"before 4c"<<endl;
820 if(m_test4C==1) {
821// double ecms = 3.097;
822 HepLorentzVector ecms(0.034,0,0,3.097);
823
824 double chisq = 9999.;
825 int ig1 = -1;
826 int ig2 = -1;
827 for(int i = 0; i < nGam-1; i++) {
828 RecEmcShower *g1Trk = (*(evtRecTrkCol->begin()+iGam[i]))->emcShower();
829 for(int j = i+1; j < nGam; j++) {
830 RecEmcShower *g2Trk = (*(evtRecTrkCol->begin()+iGam[j]))->emcShower();
831 kmfit->init();
832 kmfit->AddTrack(0, wpip);
833 kmfit->AddTrack(1, wpim);
834 kmfit->AddTrack(2, 0.0, g1Trk);
835 kmfit->AddTrack(3, 0.0, g2Trk);
836 kmfit->AddFourMomentum(0, ecms);
837 bool oksq = kmfit->Fit();
838 if(oksq) {
839 double chi2 = kmfit->chisq();
840 if(chi2 < chisq) {
841 chisq = chi2;
842 ig1 = iGam[i];
843 ig2 = iGam[j];
844 }
845 }
846 }
847 }
848
849 if(chisq < 200) {
850
851 RecEmcShower *g1Trk = (*(evtRecTrkCol->begin()+ig1))->emcShower();
852 RecEmcShower *g2Trk = (*(evtRecTrkCol->begin()+ig2))->emcShower();
853 kmfit->init();
854 kmfit->AddTrack(0, wpip);
855 kmfit->AddTrack(1, wpim);
856 kmfit->AddTrack(2, 0.0, g1Trk);
857 kmfit->AddTrack(3, 0.0, g2Trk);
858 kmfit->AddFourMomentum(0, ecms);
859 bool oksq = kmfit->Fit();
860 if(oksq) {
861 HepLorentzVector ppi0 = kmfit->pfit(2) + kmfit->pfit(3);
862 m_mpi0 = ppi0.m();
863 m_chi1 = kmfit->chisq();
864 m_tuple4->write();
865 Ncut4++;
866 }
867 }
868 }
869
870 //
871 // Apply Kinematic 5C Fit
872 //
873
874 // find the best combination over all possible pi+ pi- gamma gamma pair
875 if(m_test5C==1) {
876// double ecms = 3.097;
877 HepLorentzVector ecms(0.034,0,0,3.097);
878 double chisq = 9999.;
879 int ig1 = -1;
880 int ig2 = -1;
881 for(int i = 0; i < nGam-1; i++) {
882 RecEmcShower *g1Trk = (*(evtRecTrkCol->begin()+iGam[i]))->emcShower();
883 for(int j = i+1; j < nGam; j++) {
884 RecEmcShower *g2Trk = (*(evtRecTrkCol->begin()+iGam[j]))->emcShower();
885 kmfit->init();
886 kmfit->AddTrack(0, wpip);
887 kmfit->AddTrack(1, wpim);
888 kmfit->AddTrack(2, 0.0, g1Trk);
889 kmfit->AddTrack(3, 0.0, g2Trk);
890 kmfit->AddResonance(0, 0.135, 2, 3);
891 kmfit->AddFourMomentum(1, ecms);
892 if(!kmfit->Fit(0)) continue;
893 if(!kmfit->Fit(1)) continue;
894 bool oksq = kmfit->Fit();
895 if(oksq) {
896 double chi2 = kmfit->chisq();
897 if(chi2 < chisq) {
898 chisq = chi2;
899 ig1 = iGam[i];
900 ig2 = iGam[j];
901 }
902 }
903 }
904 }
905
906
907 log << MSG::INFO << " chisq = " << chisq <<endreq;
908
909 if(chisq < 200) {
910 RecEmcShower* g1Trk = (*(evtRecTrkCol->begin()+ig1))->emcShower();
911 RecEmcShower* g2Trk = (*(evtRecTrkCol->begin()+ig2))->emcShower();
912
913 kmfit->init();
914 kmfit->AddTrack(0, wpip);
915 kmfit->AddTrack(1, wpim);
916 kmfit->AddTrack(2, 0.0, g1Trk);
917 kmfit->AddTrack(3, 0.0, g2Trk);
918 kmfit->AddResonance(0, 0.135, 2, 3);
919 kmfit->AddFourMomentum(1, ecms);
920 bool oksq = kmfit->Fit();
921 if(oksq){
922 HepLorentzVector ppi0 = kmfit->pfit(2) + kmfit->pfit(3);
923 HepLorentzVector prho0 = kmfit->pfit(0) + kmfit->pfit(1);
924 HepLorentzVector prhop = ppi0 + kmfit->pfit(0);
925 HepLorentzVector prhom = ppi0 + kmfit->pfit(1);
926
927 m_chi2 = kmfit->chisq();
928 m_mrh0 = prho0.m();
929 m_mrhp = prhop.m();
930 m_mrhm = prhom.m();
931 double eg1 = (kmfit->pfit(2)).e();
932 double eg2 = (kmfit->pfit(3)).e();
933 double fcos = abs(eg1-eg2)/ppi0.rho();
934 m_tuple5->write();
935 Ncut5++;
936 //
937 // Measure the photon detection efficiences via
938 // J/psi -> rho0 pi0
939 //
940 if(fabs(prho0.m()-0.770)<0.150) {
941 if(fabs(fcos)<0.99) {
942 m_fcos = (eg1-eg2)/ppi0.rho();
943 m_elow = (eg1 < eg2) ? eg1 : eg2;
944 m_tuple6->write();
945 Ncut6++;
946 }
947 } // rho0 cut
948 } //oksq
949 }
950 }
951 return StatusCode::SUCCESS;
952}
953
954
955// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
956StatusCode Rhopi::finalize() {
957 cout<<"total number: "<<Ncut0<<endl;
958 cout<<"nGood==2, nCharge==0: "<<Ncut1<<endl;
959 cout<<"nGam>=2: "<<Ncut2<<endl;
960 cout<<"Pass Pid: "<<Ncut3<<endl;
961 cout<<"Pass 4C: "<<Ncut4<<endl;
962 cout<<"Pass 5C: "<<Ncut5<<endl;
963 cout<<"J/psi->rho0 pi0: "<<Ncut6<<endl;
964 MsgStream log(msgSvc(), name());
965 log << MSG::INFO << "in finalize()" << endmsg;
966 return StatusCode::SUCCESS;
967}
968
double sin(const BesAngle a)
Definition: BesAngle.h:210
double cos(const BesAngle a)
Definition: BesAngle.h:213
int runNo
Definition: DQA_TO_DB.cxx:12
Double_t time
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
const double mpi
Definition: Gam4pikp.cxx:47
std::vector< int > Vint
Definition: Gam4pikp.cxx:52
int Ncut2
Definition: Ppjrhopi.cxx:53
int Ncut1
Definition: Ppjrhopi.cxx:53
int Ncut0
Definition: Ppjrhopi.cxx:53
int Ncut5
Definition: Ppjrhopi.cxx:53
int Ncut3
Definition: Ppjrhopi.cxx:53
int Ncut6
Definition: Ppjrhopi.cxx:53
int Ncut4
Definition: Ppjrhopi.cxx:53
int Ncut2
Definition: Rhopi.cxx:54
HepGeom::Point3D< double > HepPoint3D
Definition: Rhopi.cxx:34
int Ncut1
Definition: Rhopi.cxx:54
std::vector< HepLorentzVector > Vp4
Definition: Rhopi.cxx:52
int Ncut0
Definition: Rhopi.cxx:54
const double xmass[5]
Definition: Rhopi.cxx:48
int Ncut5
Definition: Rhopi.cxx:54
const double velc
Definition: Rhopi.cxx:50
int Ncut3
Definition: Rhopi.cxx:54
const double mpi
Definition: Rhopi.cxx:47
std::vector< int > Vint
Definition: Rhopi.cxx:51
int Ncut6
Definition: Rhopi.cxx:54
int Ncut4
Definition: Rhopi.cxx:54
INTupleSvc * ntupleSvc()
IMessageSvc * msgSvc()
@ theta1
Definition: TrkKalDeriv.h:24
double theta() const
Definition: DstEmcShower.h:38
double phi() const
Definition: DstEmcShower.h:39
double x() const
Definition: DstEmcShower.h:35
double time() const
Definition: DstEmcShower.h:50
double z() const
Definition: DstEmcShower.h:37
double energy() const
Definition: DstEmcShower.h:45
double y() const
Definition: DstEmcShower.h:36
const Hep3Vector emcPosition() const
Definition: DstExtTrack.h:126
const int emcVolumeNumber() const
Definition: DstExtTrack.h:132
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 px() const
static void setPidType(PidType pidType)
const double pz() const
const double py() const
const int charge() const
const double theta() const
Definition: DstMdcTrack.h:59
const HepSymMatrix err() const
const int charge() const
Definition: DstMdcTrack.h:53
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
void AddFourMomentum(int number, HepLorentzVector p4)
HepLorentzVector pfit(int n) const
void AddResonance(int number, double mres, std::vector< int > tlis)
static KalmanKinematicFit * instance()
int useTofCorr() const
int methodProbability() const
int useDedx() const
int onlyPionKaonProton() 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 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
const HepSymMatrix & getZError() const
Definition: Rhopi.h:10
StatusCode execute()
Definition: Rhopi.cxx:309
StatusCode finalize()
Definition: Rhopi.cxx:956
StatusCode initialize()
Definition: Rhopi.cxx:76
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.
WTrackParameter wtrk(int n) const
Definition: VertexFit.h:79
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
void Swim(int n)
Definition: VertexFit.h:59
bool Fit()
Definition: VertexFit.cxx:301
void setEvx(const HepSymMatrix &eVx)
void setVx(const HepPoint3D &vx)
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
double double * p3
Definition: qcdloop1.h:76
const float pi
Definition: vector3.h:133