BOSS 7.1.2
BESIII Offline Software System
Loading...
Searching...
No Matches
ParticleID.cxx
Go to the documentation of this file.
1#include <iostream>
2// #include <cmath>
3#include <cstdlib>
4
6
7#ifndef BEAN
9#include "GaudiKernel/Bootstrap.h"
10#include "GaudiKernel/ISvcLocator.h"
11#include "GaudiKernel/ISvcLocator.h"
12#include "GaudiKernel/IDataProviderSvc.h"
13#include "GaudiKernel/SmartDataPtr.h"
15#endif
16
17//
18// Author: K.L. He & L. L. Wang & Gang.Qin, 01/07/2007 created
19
20ParticleID * ParticleID::m_pointer = 0;
21
23 if(!m_pointer) m_pointer = new ParticleID();
24 return m_pointer;
25}
26
28
29 if(IsDedxInfoUsed()) {
30 if(!m_dedxpid) m_dedxpid = DedxPID::instance();
31 m_dedxpid->init();
32 }
33
34 if(IsTofInfoUsed()|IsTof1InfoUsed()|IsTof2InfoUsed()) {
35 if(!m_tofpid) m_tofpid = TofPID::instance();
36 m_tofpid->init();
37 }
38
39 if(IsTofCorrInfoUsed()) {
40 if(!m_tofcorrpid) m_tofcorrpid = TofCorrPID::instance();
41 // m_tofcorrpid->init();
42 }
43
44 if(IsEmcInfoUsed()) {
45 if(!m_emcpid) m_emcpid = EmcPID::instance();
46 m_emcpid->init();
47 }
48
49 if(IsMucInfoUsed()) {
50 if(!m_mucpid) m_mucpid = MucPID::instance();
51 m_mucpid->init();
52 }
53
54 // global info.
55 m_pidsys = 0;
56 m_pidcase = 0;
57 m_method = 0;
58 m_TotalLikelihood =0;
59 m_discard =1;
60 // probability
61 m_ndof = 0;
62 m_nhitcut=5;
63 // chi cut
64 setChiMinCut( 4.0 );
65 setChiMaxCut( 6.0 );
66 for(int i = 0; i < 5; i++) {
67 m_chisq[i] = 9999.;
68 m_prob[i] = -1.0;
69 }
70}
71
72
73ParticleID::ParticleID() : ParticleIDBase() {
74 m_dedxpid = 0;
75 m_tofpid = 0;
76 m_tofepid = 0;
77 m_tofqpid = 0;
78 m_tofcpid = 0;
79 m_tofcorrpid = 0;
80 m_emcpid = 0;
81 m_mucpid = 0;
82
83
84
85}
86
87
89 // if(m_dedxpid) delete m_dedxpid;
90 // if(m_tof1pid) delete m_tof1pid;
91 // if(m_tof2pid) delete m_tof2pid;
92 // if(m_tofepid) delete m_tofepid;
93 // if(m_tofqpid) delete m_tofqpid;
94 // if(m_emcpid) delete m_emcpid;
95}
96
98#ifdef BEAN
99 cout << " please use ParticleID::calculate(run) ! " << endl;
100 exit(1);
101}
102
103
104void ParticleID::calculate(int run)
105{
106#endif
107 int nhitcutpid=getNhitCut();
108
109#ifndef BEAN
110 IDataProviderSvc* m_eventSvc;
111 Gaudi::svcLocator()->service("EventDataSvc", m_eventSvc, true);
112
113 SmartDataPtr<Event::EventHeader> eventHeaderpid(m_eventSvc,"/Event/EventHeader");
114 int runpid=eventHeaderpid->runNumber();
115 int eventpid=eventHeaderpid->eventNumber();
116 // cout<<"runpid="<<runpid<<endl;
117 // cout<<"eventpid="<<eventpid<<endl;
118#else
119 int runpid=run;
120#endif
121
122 EvtRecTrack* recTrk = PidTrk();
123 // int runnum=getRunNo();
124 // cout<<"runnum="<<runnum<<endl;
125 // if user did not specify sub sys, sepcify the default value
126 if(m_pidsys == 0) {
127 m_pidsys = useDedx() | useTof() | useTofE() | useEmc() | useMuc() | useTofQ() | useTofC() | useTofCorr();
128 }
129 // if user did not set the seperate case, set the default value
130
131 if(m_pidcase == 0 ) {
132 m_pidcase = all();
133 }
134 //dedx sys
135 if(IsDedxInfoUsed()) {
136 if(!m_dedxpid) m_dedxpid = DedxPID::instance();
137 m_dedxpid->init();
138 m_dedxpid->setRunNo(runpid);
139 m_dedxpid->setNhitCutDx(nhitcutpid);
140 m_dedxpid->setRecTrack(recTrk);
141 m_dedxpid->setChiMinCut(chiMinCut());
142 m_dedxpid->setPdfMinSigmaCut(pdfMinSigmaCut());
143 m_dedxpid->calculate();
144 }
145
146 // tof1 and tof2 sys
147 if(IsTofInfoUsed()|IsTof1InfoUsed()|IsTof2InfoUsed()|IsTofCInfoUsed())
148 {
149 if(IsTofCInfoUsed())
150 {
151 if(!m_tofcpid) m_tofcpid = TofCPID::instance();
152 m_tofcpid->init();
153 m_tofcpid->setRunNo(runpid);
154 m_tofcpid->setRecTrack(recTrk);
155 m_tofcpid->setChiMinCut(chiMinCut());
156 m_tofcpid->setPdfMinSigmaCut(pdfMinSigmaCut());
157 m_tofcpid->calculate();
158 }
159 else
160 {
161 if(!m_tofpid) m_tofpid = TofPID::instance();
162 m_tofpid->init();
163 m_tofpid->setRecTrack(recTrk);
164 m_tofpid->setChiMinCut(chiMinCut());
166 m_tofpid->calculate();
167 }
168
169 }
170
171
172 // tof secondary correction sys
173 if(IsTofCorrInfoUsed()) {
174 if(!m_tofcorrpid) m_tofcorrpid = TofCorrPID::instance();
175 m_tofcorrpid->setRunNo(runpid);
176 m_tofcorrpid->init();
177 m_tofcorrpid->setRecTrack(recTrk);
178 m_tofcorrpid->setChiMinCut(chiMinCut());
179 m_tofcorrpid->setChiMaxCut(chiMaxCut());
180 m_tofcorrpid->setPdfMinSigmaCut(pdfMinSigmaCut());
181 m_tofcorrpid->calculate();
182 }
183
184 /*
185 // tof1 sys
186 if(IsTof1InfoUsed()){
187 if(!m_tof1pid) m_tof1pid = Tof1PID::instance();
188 m_tof1pid->init();
189 m_tof1pid->setRecTrack(recTrk);
190 m_tof1pid->setChiMinCut(4);
191 m_tof1pid->setPdfMinSigmaCut(4);
192 m_tof1pid->calculate();
193 }
194
195 // tof2 sys
196 if(IsTof2InfoUsed()){
197 if(!m_tof2pid) m_tof2pid = Tof2PID::instance();
198 m_tof2pid->init();
199 m_tof2pid->setRecTrack(recTrk);
200 m_tof2pid->setChiMinCut(4);
201 m_tof2pid->setPdfMinSigmaCut(4);
202 m_tof2pid->calculate();
203 }
204
205 */
206 // tofq sys
207 if(IsTofQInfoUsed()) {
208 if(!m_tofqpid) m_tofqpid = TofQPID::instance();
209 m_tofqpid->init();
210 m_tofqpid->setRecTrack(recTrk);
211 m_tofqpid->setChiMinCut(chiMinCut());
212 m_tofqpid->calculate();
213 }
214
215 // endcap tof sys
216 if(IsTofEInfoUsed()&&(!IsTofCorrInfoUsed())) {
217 if(!m_tofepid) m_tofepid = TofEPID::instance();
218 m_tofepid->init();
219 m_tofepid->setRecTrack(recTrk);
220 m_tofepid->setChiMinCut(chiMinCut());
221 m_tofepid->setPdfMinSigmaCut(pdfMinSigmaCut());
222 m_tofepid->calculate();
223 }
224 // emc sys
225 if(IsEmcInfoUsed()) {
226 if(!m_emcpid) m_emcpid = EmcPID::instance();
227 m_emcpid->init();
228 m_emcpid->setRecTrack(recTrk);
229 m_emcpid->setChiMinCut(chiMinCut());
230 m_emcpid->calculate();
231 }
232
233 // muc sys
234 if(IsMucInfoUsed()) {
235 if(!m_mucpid) m_mucpid = MucPID::instance();
236 m_mucpid->init();
237 m_mucpid->setRecTrack(recTrk);
238 m_mucpid->setChiMinCut(chiMinCut());
239 m_mucpid->calculate();
240 }
241 // probability method
242 if((m_method & methodProbability()) == methodProbability())
243 if(particleIDCalculation() < 0) m_ndof = 0;
244 // std::cout<<"m_ndof="<<m_ndof<<std::endl;
245
246 //likelihood method
247 if((m_method & methodLikelihood()) == methodLikelihood())
248 if(LikelihoodCalculation() < 0) m_discard =0;
249 // neuron network
250 if((m_method & methodNeuronNetwork()) == methodNeuronNetwork())
251 // m_neuronPid = neuronPIDCalculation();
252 if(LikelihoodCalculation() < 0) m_discard =0;
253
254}
255
256int ParticleID ::particleIDCalculation() {
257 int irc = -1;
258 bool valid = IsDedxInfoValid() || IsTofInfoValid()||IsTofEInfoValid()
259 || IsTofQInfoValid() || IsEmcInfoValid() || IsMucInfoValid()
260 || IsTofCInfoValid() || IsTofCorrInfoValid();
261
262 if(!valid) return irc;
263
264 double chisq[5];
265 bool pidcase[5];
266 for(int i = 0; i < 5; i++) {
267 chisq[i] = 0;
268 pidcase[i] = false;
269 }
270
271 if((m_pidcase & onlyElectron()) == onlyElectron()) pidcase[0] = true;
272 if((m_pidcase & onlyMuon()) == onlyMuon()) pidcase[1] = true;
273 if((m_pidcase & onlyPion()) == onlyPion()) pidcase[2] = true;
274 if((m_pidcase & onlyKaon()) == onlyKaon()) pidcase[3] = true;
275 if((m_pidcase & onlyProton()) == onlyProton()) pidcase[4] = true;
276
277 //
278 // dEdx PID
279 //
280 if(IsDedxInfoUsed()) {
281 if(IsDedxInfoValid()) {
282 bool okpid = false;
283 for(int i = 0; i < 5; i++) {
284 if(pidcase[i] && (fabs(chiDedx(i)) < m_dedxpid->chiMinCut()))
285 if(!okpid) okpid = true;
286 }
287 if(okpid) {
288 m_ndof++;
289 for(int i = 0; i < 5; i++) chisq[i] += chiDedx(i)*chiDedx(i);
290
291
292 }
293 } // dE/dx
294 }
295 //
296 // Barrel TOF
297 //
298
299 if(IsTofInfoUsed()|IsTof1InfoUsed()|IsTof2InfoUsed() | IsTofCInfoUsed())
300 { if(IsTofCInfoUsed())
301 {
302 if(IsTofCInfoValid()) {
303 bool okpid = false;
304 for(int i = 0; i < 5; i++) {
305 if(pidcase[i] && (fabs(chiTofC(i)) < m_tofcpid->chiMinCut()))
306 if(!okpid) okpid = true;
307 }
308 if(okpid) {
309 m_ndof++;
310 for(int i = 0; i < 5; i++) chisq[i] += chiTofC(i)*chiTofC(i);
311 }
312 } // TOF1
313 }
314 else {
315 if(IsTofInfoValid()) {
316 bool okpid = false;
317 for(int i = 0; i < 5; i++) {
318 if(pidcase[i] && (fabs(chiTof(i)) < m_tofpid->chiMinCut()))
319 if(!okpid) okpid = true;
320 }
321 if(okpid) {
322 m_ndof++;
323 for(int i = 0; i < 5; i++) chisq[i] += chiTof(i)*chiTof(i);
324 }
325 } // TOF1
326
327
328 //
329 // EndCap Tof
330 //
331
332 if(IsTofEInfoUsed()) {
333 if(IsTofEInfoValid()) {
334 bool okpid = false;
335 for(int i = 0; i < 5; i++) {
336 if(pidcase[i] && (fabs(chiTofE(i)) < m_tofepid->chiMinCut()))
337 if(!okpid) okpid = true;
338 }
339 if(okpid) {
340 m_ndof++;
341 for(int i = 0; i < 5; i++) chisq[i] += chiTofE(i)*chiTofE(i);
342 }
343 } // EndCap TOF
344 }
345
346 }
347 }
348
349 // Secondary TOF correction
350 if(IsTofCorrInfoUsed()) {
351 if(IsTofCorrInfoValid()) {
352 bool okpid = false;
353 for(int i = 0; i < 5; i++) {
354 if(pidcase[i] && ( chiTofCorr(i) < m_tofcorrpid->chiMaxCut() ) && ( chiTofCorr(i) > ( 0.0 - m_tofcorrpid->chiMinCut() ) ) )
355 // if(pidcase[i] && ( chiTofCorr(i) < 6.0 && chiTofCorr(i) > -4.0 ) )
356 if(!okpid) okpid = true;
357 }
358 if(okpid) {
359 m_ndof++;
360 for(int i = 0; i < 5; i++) chisq[i] += chiTofCorr(i)*chiTofCorr(i);
361 }
362 }
363 }
364
365 //
366 // Barrel TOF Q
367 //
368
369 if(IsTofQInfoUsed()) {
370 if(IsTofQInfoValid()) {
371 bool okpid = false;
372 for(int i = 0; i < 5; i++) {
373 if(pidcase[i] && (fabs(chiTofQ(i)) < m_tofqpid->chiMinCut()))
374 if(!okpid) okpid = true;
375 }
376 if(okpid) {
377 m_ndof++;
378 for(int i = 0; i < 5; i++) chisq[i] += chiTofQ(i)*chiTofQ(i);
379 }
380 } // TofQ
381 }
382
383 // Muc Pid
384 if(IsMucInfoUsed()) {
385 if(IsMucInfoValid()) {
386 bool okpid = false;
387 for(int i = 0; i < 5; i++) {
388 if(pidcase[i] && (fabs(chiMuc(i)) < m_mucpid->chiMinCut()))
389 if(!okpid) okpid = true;
390 }
391 if(okpid) {
392 m_ndof++;
393 for(int i = 0; i < 5; i++) chisq[i] += chiMuc(i)*chiMuc(i);
394 }
395 } // Muc Pid
396 }
397
398
399 // Emc PID
400 if(IsEmcInfoUsed()) {
401 if(IsEmcInfoValid()) {
402 bool okpid = false;
403 for(int i = 0; i < 5; i++) {
404 if(pidcase[i] && (fabs(chiEmc(i)) < m_emcpid->chiMinCut()))
405 if(!okpid) okpid = true;
406 }
407 if(okpid) {
408 m_ndof++;
409 for(int i = 0; i < 5; i++) chisq[i] += chiEmc(i)*chiEmc(i);
410 }
411 } // Emc Pid
412 }
413
414
415 if(m_ndof <= 0) return irc;
416
417
418 for(int i = 0; i < 5; i++) {
419 m_chisq[i] = chisq[i];
420 m_prob[i] = probCalculate(chisq[i], m_ndof);
421 }
422
423
424 irc = 0;
425 return irc;
426}
427
428
429
430
431
433 int irc = -1;
434
436 if(!valid) return irc;
437 double pdf[5];
438 bool pidcase[5];
439 for(int i = 0; i < 5; i++) {
440 pdf[i] = 1;
441 pidcase[i] = false;
442 }
443
444 if((m_pidcase & onlyElectron()) == onlyElectron()) pidcase[0] = true;
445 if((m_pidcase & onlyMuon()) == onlyMuon()) pidcase[1] = true;
446 if((m_pidcase & onlyPion()) == onlyPion()) pidcase[2] = true;
447 if((m_pidcase & onlyKaon()) == onlyKaon()) pidcase[3] = true;
448 if((m_pidcase & onlyProton()) == onlyProton()) pidcase[4] = true;
449
450 for(int i = 0; i < 5; i++) {
451 if(pidcase[i]==0)
452 pdf[i]=0;
453 }
454
455 //
456 // dEdx PID
457 //
458 if(IsDedxInfoUsed()) {
459 if(IsDedxInfoValid()) {
460 bool okpid = false;
461 for(int i = 0; i < 5; i++) {
462 if(pidcase[i] && pdfCalculate(chiDedx(i),1) > pdfCalculate(m_dedxpid->pdfMinSigmaCut(),1.5))
463 if(!okpid) okpid = true;
464 }
465 if(okpid) {
466 m_ndof++;
467 for(int i = 0; i < 5; i++) {
468 pdf[i] *= pdfDedx(i);
469 }
470 }
471 } // dE/dx
472 }
473
474
475 //
476 // Barrel TOF
477 //
478 if(IsTofInfoUsed()|IsTof1InfoUsed()|IsTof2InfoUsed()|IsTofCInfoUsed())
479 { if(IsTofCInfoUsed())
480 {
481
482 if(IsTofCInfoValid()) {
483 bool okpid = false;
484 for(int i = 0; i < 5; i++) {
485 if(pidcase[i] && pdfCalculate(chiTof(i),1) > pdfCalculate(m_tofcpid->pdfMinSigmaCut(),1.5))
486 if(!okpid) okpid = true;
487 }
488 if(okpid) {
489 m_ndof++;
490 for(int i = 0; i < 5; i++) {
491 pdf[i] *= pdfTofC(i);
492 }
493 }
494 } // TOF
495 }
496
497 else {
498 if(IsTofInfoValid()) {
499 bool okpid = false;
500 for(int i = 0; i < 5; i++) {
501 if(pidcase[i] && pdfCalculate(chiTof(i),1) > pdfCalculate(m_tofpid->pdfMinSigmaCut(),1.5))
502 if(!okpid) okpid = true;
503 }
504 if(okpid) {
505 m_ndof++;
506 for(int i = 0; i < 5; i++) {
507 pdf[i] *= pdfTof(i);
508 }
509 }
510 } // TOF
511
512
513
514 //
515 // EndCap Tof
516 //
517
518 if(IsTofEInfoUsed()) {
519 if(IsTofEInfoValid()) {
520 bool okpid = false;
521 for(int i = 0; i < 5; i++) {
522 if(pidcase[i]&& pdfCalculate(chiTofE(i),1) > pdfCalculate(m_tofepid->pdfMinSigmaCut(),1.5))
523 if(!okpid) okpid = true;
524 }
525 if(okpid) {
526 // m_ndof++;
527 // for(int i = 0; i < 5; i++) pdf[i] *= pdfTofE(i);
528 }
529 } // EndCap TOF
530 }
531 }
532
533 }
534
535 // Secondary TOF correction
536 if(IsTofCorrInfoUsed()) {
537 if(IsTofCorrInfoValid()) {
538 bool okpid = false;
539 for(int i = 0; i < 5; i++) {
540 if(pidcase[i] && pdfCalculate(chiTofCorr(i),1) > pdfCalculate(m_tofcorrpid->pdfMinSigmaCut(),1.5))
541 if(!okpid) okpid = true;
542 }
543 if(okpid) {
544 m_ndof++;
545 for(int i = 0; i < 5; i++) {
546 pdf[i] *= pdfTofCorr(i);
547 }
548 }
549 }
550 }
551
552
553 //
554 // Barrel TOF Q
555 //
556
557 if(IsTofQInfoValid()) {
558 bool okpid = false;
559 for(int i = 0; i < 5; i++) {
560 if(pidcase[i])
561 if(!okpid) okpid = true;
562 }
563 if(okpid) {
564 // m_ndof++;
565 for(int i = 0; i < 5; i++) pdf[i] *= pdfTofQ(i);
566 }
567 } // TofQ
568
569 //
570 // Emc PID
571 //
572 if(IsEmcInfoUsed()) {
573 if(IsEmcInfoValid()) {
574 bool okpid = false;
575 for(int i = 0; i < 5; i++) {
576 if(pidcase[i]&&pdfEmc(i)>0)
577 if(!okpid) okpid = true;
578 }
579 if(okpid) {
580 m_ndof++;
581 for(int i = 0; i < 5; i++) {
582 pdf[i] *= pdfEmc(i);
583 }
584 } // Emc Pid
585 }
586 }
587 if(IsMucInfoUsed()) {
588 if(IsMucInfoValid()) {
589 bool okpid = false;
590 for(int i = 0; i < 5; i++) {
591 if(pidcase[i]&&pdfMuc(i)>0)
592 if(!okpid) okpid = true;
593 }
594 if(okpid) {
595 m_ndof++;
596 for(int i = 0; i < 5; i++) {
597 pdf[i] *= pdfMuc(i);
598 }
599 }
600 } // Emc Pid
601 }
602
603
604
605 if(m_ndof <= 0) return irc;
606 for(int i = 0; i < 5; i++) {
607 m_pdf[i] = pdf[i];
608 m_TotalLikelihood += pdf[i];
609 }
610 for(int i = 0; i < 5; i++) {
611 m_likelihoodfraction[i] = pdf[i]/m_TotalLikelihood;
612 }
613
614
615 irc = 0;
616 return irc;
617}
618
void init()
Definition DedxPID.cxx:26
void setNhitCutDx(const int nhitcuthdx=5)
Definition DedxPID.h:35
void calculate()
Definition DedxPID.cxx:42
static DedxPID * instance()
Definition DedxPID.cxx:17
void calculate()
Definition EmcPID.cxx:153
void init()
Definition EmcPID.cxx:95
static EmcPID * instance()
Definition EmcPID.cxx:23
void init()
Definition MucPID.cxx:67
static MucPID * instance()
Definition MucPID.cxx:22
void calculate()
Definition MucPID.cxx:101
double chiMinCut() const
int useTofCorr() const
void setPdfMinSigmaCut(const double pdf=4)
EvtRecTrack * PidTrk() const
int useTofE() const
int onlyProton() const
double chiMaxCut() const
int useTof() const
int useTofQ() const
int methodProbability() const
int useDedx() const
int useTofC() const
double pdfCalculate(double offset, double sigma)
void setRunNo(const double runh=8093)
int onlyMuon() const
double pdfMinSigmaCut() const
int onlyKaon() const
int onlyElectron() const
int onlyPion() const
int methodNeuronNetwork() const
void setChiMaxCut(const double chi=6)
int useEmc() const
void setChiMinCut(const double chi=4)
int all() const
int useMuc() const
int methodLikelihood() const
void setRecTrack(EvtRecTrack *trk)
double pdfEmc(int n)
double pdfTofC(int n)
int getNhitCut() const
Definition ParticleID.h:88
bool IsTofEInfoValid() const
bool IsEmcInfoValid() const
double pdfDedx(int n)
int LikelihoodCalculation()
double chiTofE(int n) const
bool IsTofCorrInfoValid() const
bool IsTofCInfoValid() const
bool IsTofInfoValid() const
double pdfTof(int n)
bool IsTofQInfoValid() const
int particleIDCalculation()
static ParticleID * instance()
bool IsDedxInfoValid() const
double pdfMuc(int n)
double pdfTofCorr(int n)
void calculate()
void init()
bool IsMucInfoValid() const
double pdf(int n) const
Definition ParticleID.h:118
double chiTof(int n) const
double pdfTofQ(int n)
double chiTofCorr(int n) const
double chiDedx(int n) const
void calculate()
Definition TofCPID.cxx:55
void init()
Definition TofCPID.cxx:43
static TofCPID * instance()
Definition TofCPID.cxx:19
void init()
void calculate()
static TofCorrPID * instance()
void init()
Definition TofEPID.cxx:25
void calculate()
Definition TofEPID.cxx:39
static TofEPID * instance()
Definition TofEPID.cxx:16
void init()
Definition TofPID.cxx:24
void calculate()
Definition TofPID.cxx:48
static TofPID * instance()
Definition TofPID.cxx:15
void init()
Definition TofQPID.cxx:22
void calculate()
Definition TofQPID.cxx:34
static TofQPID * instance()
Definition TofQPID.cxx:13