FairRoot/PandaRoot
softtrigger_kin5.C
Go to the documentation of this file.
1 //#include "SimpleCand.C"
2 #include <iostream>
3 #include <vector>
4 #include <map>
5 
6 #include "TFile.h"
7 #include "TTree.h"
8 #include "TRandom3.h"
9 #include "TDatabasePDG.h"
10 #include "TParticlePDG.h"
11 #include "TCanvas.h"
12 #include "TH1F.h"
13 #include "TH2F.h"
14 #include "TLegend.h"
15 #include "TParticle.h"
16 #include "TClonesArray.h"
17 #include "TF1.h"
18 #include "TStyle.h"
19 #include "TNtuple.h"
20 #include "TLatex.h"
21 #include "TVector.h"
22 #include "TMatrixD.h"
23 #include "TSystem.h"
24 #include "TMatrixDEigen.h"
25 #include "TMath.h"
26 
27 #include "EventShape.h"
28 
29 
30 using std::cout;
31 using std::endl;
32 
33 TRandom3 fRand;
34 TLorentzVector fIni;
35 
36 // *** some program parameters
37 const int fMAX = 100; // *** max tracks per event
38 int nbins = 300; // *** bins of mass histograms
39 double low = 0.0; // *** lower edge of histo
40 double high = 5.0; // *** higher edge of histo (will be modified acc. to resonance mass)
41 
42 // *** some contants for neutrals
43 const double fNeutEff = 0.9; // *** efficiency for neutrals
44 const double fdE = 0.03; // *** neutral energy resolution
45 const double fNdtht = 0.003; // *** neutrals tht resolution
46 const double fNdphi = 0.003; // *** neutrals phi resolution
47 const double fNeutThrsh = 0.02; // *** energy threshold of neutrals
48 
49 // *** constants for particles
50 int fPi0Pdg = 111; // *** PDG code for pi0
51 double fPi0Mass = 0.13497; // *** center of mass cut window
52 double fPi0win = 0.020;
53 
54 int fKsPdg = 310; // *** PDG code for pi0
55 double fKsMass = 0.497614; // *** center of mass cut window
56 double fKswin = 0.025;
57 
58 int fDPdg = 411; // *** PDG code to match D
59 double fDMass = 1.8693; // *** center of mass cut window
60 double fDsig = 0.88; // *** slope of linear width function
61 double fDwin =0.25;
62 
63 int fD0Pdg = 421; // *** PDG code to match D0
64 double fD0Mass = 1.8645; // *** center of mass cut window
65 double fD0sig = 1.04; // *** slope of linear width function
66 double fD0win =0.25;
67 
68 int fJpsiPdg = 443; // *** PDG code to match J/psi
69 double fJMass = 3.09687; // *** center of mass cut window
70 double fJsig = 2.24; // *** slope of linear width function
71 double fJsig3pi = 1.5; // *** slope of linear width function for 3pi mode
72 double fJ3pioff = 0.027; // *** offset of linear function for 3pi mode
73 double fJwin =0.7;
74 
75 int fLamPdg = 3122; // *** PDG code to match Lam
76 double fLamMass = 1.115684; // *** center of mass cut window
77 double fLamsig = 0.08; // *** slope of linear width function
78 double fLamwin =0.2;
79 
80 int fDsPdg = 431; // *** PDG code to match Ds
81 double fDsMass = 1.9686; // *** center of mass cut window
82 double fDssig = 0.75; // *** slope of linear width function
83 double fDswin =0.2;
84 
85 int fPhiPdg = 333; // *** PDG code to match Phi
86 double fPhiMass = 1.0195; // *** center of mass cut window
87 double fPhisig = 0.15; // *** slope of linear width function
88 double fPhiwin =0.2;
89 
90 int fEtacPdg = 441; // *** PDG code to match etac
91 double fEtacMass = 2.9798; // *** center of mass cut window
92 double fEtacsig = 0.15; // *** slope of linear width function
93 double fEtacwin =0.15;
94 
95 int fLamcPdg = 4122; // *** PDG code to match Lambda_c
96 double fLamcMass = 2.2849; // *** center of mass cut window
97 double fLamcsig = 0.15; // *** slope of linear width function
98 double fLamcwin =0.2;
99 
100 
101 // *** for event shape variables (global)
102 double _sph, _pla, _apl, _thr;
103 TVector3 _sphAx, _secAx, _thrAx;
104 
105 // *** max order for Fox-Wolfram Moments
106 const int _nfwmom = 5;
107 // *** result container for computation (global)
109 
110 // *** pid probability table; will be modified according to purity parameter value
111 // *** true pid: e mu pi K p classified as [%]
112 double pidtab [25] = { 80, 5, 5, 5, 5, // e
113  5, 80, 5, 5, 5, // mu
114  5, 5, 80, 5, 5, // pi
115  5, 5, 5, 80, 5, // K
116  5, 5, 5, 5, 80
117  }; // p
118 
119 // *** table determined by Donghee 11/2012 for P>0.1
120 double pidreal [25] = { 95, 6, 10, 7, 10, // e
121  5, 86, 27, 6, 7, // mu
122  9, 29, 85, 18, 16, // pi
123  9, 21, 40, 67, 24, // K
124  12, 12, 20, 16, 90
125  }; // p
126 
127 // *** global stable particle masses used in pid selection; set in init()
128 double mE, mMu, mPi, mK, mP, mPi0;
129 
130 TF1 *f_res=0;
131 
132 // *** this map contains mapping from PDG code to indices 0...4 (for pidtab)
133 // *** 11 -> 0, 13 -> 1, 211 -> 2, 321 -> 3, 2212 -> 4
134 std::map<int,int> pdgidx;
135 
136 TDatabasePDG *pdg_db;
137 
138 typedef std::vector<SimpleCand> CandList;
139 
140 // *** all global particle lists we need to do combinatorics
142 
143 // *** the global pid lists
145 
146 // *** boost full list with -boost and store in boostlist
147 void boostList(CandList &list, CandList &boostlist, TVector3 boost)
148 {
149  int N=list.size();
150  boostlist.clear();
151 
152  for (unsigned int i=0;i<N;i++)
153  {
154  SimpleCand c=list[i];
155  TLorentzVector lv(c.P4());
156  lv.Boost(-boost);
157  c.SetP4(lv);
158  boostlist.push_back(c);
159  }
160 }
161 
162 // *** aux for thrust; eps(x)=sign(x)
163 double eps(TVector3 v1, TVector3 v2)
164 {
165  double prod = v1.Dot(v2);
166  return prod/fabs(prod);
167 }
168 
169 // *** routine to compute thrust vector and thrust
170 double compThrust(CandList& list, TVector3 &axis)
171 {
172  TVector3 n0;
173  int N=list.size();
174 
175  int i,j;
176  double pmax=0;
177 
178  //find starting vector as maximum momentum vector
179  for (i=0;i<N;++i)
180  {
181  if (list[i].P4().Vect().Mag()>pmax)
182  {
183  n0=list[i].P4().Vect();
184  pmax=list[i].P4().Vect().Mag();
185  }
186  }
187  //n0=list[0].P4().Vect();
188 
189  TVector3 nnew(0,0,0);
190  // find thrust axis (10 iterations)
191  for (i=0;i<10;++i)
192  {
193  for (j=0;j<N;++j)
194  {
195  //nnew+=eps(n0,list[j].P4().Vect())*list[j].P4().Vect();
196  nnew+=n0*list[j].P4().Vect()>0 ? list[j].P4().Vect() : -list[j].P4().Vect();
197  }
198  n0=nnew.Unit();
199  }
200 
201  axis = n0;
202  double thr=0, sum=0;
203  for (i=0;i<N;++i)
204  {
205  thr+=fabs(axis.Dot(list[i].P4().Vect()));
206  sum+=list[i].P4().Vect().Mag();
207  }
208  return thr/sum;
209 }
210 
211 // *** Legendre Function; auxiliary function for Fox-Wolfram-Moments
212 double Legendre( int l, double x )
213 {
214  assert(fabs(x) <= 1.);
215 
216  if (l==0) return 1.;
217 
218  double pmm = 1.;
219  double pmmp1 = x;
220 
221  if (l>1)
222  {
223  for(int ll=2; ll<=l; ll++)
224  {
225  double pll = (x * (2 * ll - 1) * pmmp1 - (ll - 1) * pmm) / (ll);
226  pmm = pmmp1;
227  pmmp1 = pll;
228  }
229  }
230  return pmmp1;
231 }
232 
233 // *** compute relative Fox-Wolfram-Moments R up to 4th order (see _nfwmom contant above)
235 {
236  int n=list.size();
237  if( n==0 ) return;
238 
239  double s = 0.;
240  double _sumarray[_nfwmom]={0,0,0,0};
241  int i,j,l;
242 
243  for (i=0; i<n-1; ++i)
244  {
245  // this candidate's 3-momentum
246  TVector3 p1(list[i].P4().Vect());
247  double pmag1 = p1.Mag();
248 
249  // loop over other candidates, starting at the next one in the list
250 
251  for (j=i+1; j<n; ++j)
252  {
253  // this candidate's 3-momentum
254  TVector3 p2(list[j].P4().Vect());
255  double pmag2 = p2.Mag();
256 
257  // the cosine of the angle between the two candidates
258  double cosPhi = cos ( p1.Angle(p2) );
259 
260  // the contribution of this pair of track
261  // (note the factor 2 : the pair enters the sum twice)
262  for( l=0; l<_nfwmom; l++ )
263  _sumarray[l] += 2 * pmag1 * pmag2 * Legendre( l, cosPhi );
264 
265  }
266  // contribution of this track
267  for( l=0; l<_nfwmom; l++ )
268  _sumarray[l] += pmag1 * pmag1 * Legendre( l, 1. );
269 
270  // total energy
271  s += list[i].P4().Energy();
272  }
273 
274  // well ...
275  if( s<=0. ) return;
276  double ssq=s*s;
277  // normalize Fox Wolfram Moments
278  for( i=0; i<_nfwmom; i++)
279  _foxwolfR[i] = _sumarray[i]/ssq ;
280 }
281 
282 // *** routine to compute sphericity, planarity and aplanarity
284 {
285  int N=list.size();
286 
287  double stot=0, sxx=0, sxy=0, sxz=0, syy=0, syz=0, szz=0;
288  int i;
289 
290  for (i=0;i<N;++i)
291  {
292  TVector3 v(list[i].P4().Vect());
293 
294  sxx += v.X()*v.X(); sxy += v.X()*v.Y(); sxz += v.X()*v.Z();
295  syy += v.Y()*v.Y(); syz += v.Y()*v.Z();
296  szz += v.Z()*v.Z();
297 
298  stot += v.Mag2();
299  }
300 
301  TMatrixD sm(3,3);
302  sm(0,0) = sxx/stot; sm(0,1) = sxy/stot; sm(0,2) = sxz/stot;
303  sm(1,0) = sxy/stot; sm(1,1) = syy/stot; sm(1,2) = syz/stot;
304  sm(2,0) = sxz/stot; sm(2,1) = syz/stot; sm(2,2) = szz/stot;
305 
306  TMatrixDEigen ei(sm);
307  TMatrixD eiv=ei.GetEigenValues();
308 
309  _sph = 1.5 * (eiv(1,1) + eiv(2,2));
310  _apl = 1.5 * eiv(2,2);
311  _pla = eiv(1,1) - eiv(2,2);
312 
313  return _sph;
314 }
315 
316 // *** check PDG code is final state particle
317 bool isFinalState(int n)
318 {
319  n=abs(n);
320  return (n==11 || n==13 || n==211 || n==321 || n==2212 || n==22);
321 }
322 
323 // *** online momentum resoltion based on talk of Yutie Lang, PANDA CM 12/2011, TRK oder DAQ session
324 // *** values (pt[GeV/c], dp/p[%]): (1.0, 3.8) (2.0, 5.5) (5.0, 12.2)
325 // *** see init routine for function parameters
326 double dponline(TLorentzVector &v)
327 {
328  return f_res->Eval(v.Pt())/100.;
329 }
330 
331 // *** computes the pid table for given efficiency and purity
332 void makePidTable(double pideff, double misid)
333 {
334  cout<< "Using PID table:"<<endl;
335 
336  if (pideff<0 && misid<0)
337  {
338  for (int i=0;i<5;++i)
339  {
340  for (int j=0;j<5;++j)
341  {
342  //pidtab[j+5*i] = pidreal[j+5*i];
343  printf("%4.0f",pidtab[j+5*i]);
344  //cout << pidtab[i+5*j]<<" ";
345  }
346  cout<< endl;
347  }
348  return;
349  }
350 
351  if (pideff<0) pideff=0;
352  if (pideff>100) pideff=100;
353  if (misid>100) misid=100;
354 
355  // *** neff is the efficiency for particles of wrong type
356  //double neff = pideff*(100-pidpur)/(4*pidpur);
357  //if (neff>100) neff=100;
358 
359  for (int i=0;i<5;++i)
360  {
361  for (int j=0;j<5;++j)
362  {
363  if (i==j) pidtab[i+5*j] = pideff;
364  else pidtab[i+5*j] = misid;
365  cout << pidtab[i+5*j]<<" ";
366  }
367  cout<< endl;
368  }
369 }
370 
371 // *** initialization of some globals
372 void init(double pideff, double misid)
373 {
374  // *** setting the mapping of PDG code to tab index
375  pdgidx.clear();
376  pdgidx[11]=0;
377  pdgidx[13]=1;
378  pdgidx[211]=2;
379  pdgidx[321]=3;
380  pdgidx[2212]=4;
381 
382  pdg_db=TDatabasePDG::Instance();
383 
384  // *** create pid table
385  makePidTable(pideff, misid);
386 
387  // *** gets particle masses from pdg database and set them globally
388  mE = pdg_db->GetParticle(11)->Mass();
389  mMu = pdg_db->GetParticle(13)->Mass();
390  mPi = pdg_db->GetParticle(211)->Mass();
391  mK = pdg_db->GetParticle(321)->Mass();
392  mP = pdg_db->GetParticle(2212)->Mass();
393 
394  mPi0= pdg_db->GetParticle(111)->Mass();
395 
396  // *** resolution fcn based on Yutie Lang results; see routine dponline above
397  f_res=new TF1("f_res","pol2(0)");
398  f_res->SetParameters(2.367,1.3,0.133);
399  f_res->SetRange(0,10);
400 }
401 
402 // *** prints single candidate of type SimpleCand
404 {
405  cout <<c.Id()<<" pid:"<<c.Pid()<<" ("<<c.P4().X()<<","<<c.P4().Y()<<","<<c.P4().Z()<<";"<<c.P4().T()<<") ";
406  cout <<c.Charge()<<" "<<c.P4().M();
407  cout <<" moth:"<<c.MotherIdx()<<" nsib:"<<c.NSiblings()<<" nfs:"<<c.NFS();
408  cout <<" mcpid:"<<c.McPid()<<" mct:"<<c.Mct()<<endl;
409 }
410 
411 // *** prints list of candidates
412 void printList(CandList &l, TString tit="")
413 {
414  if (tit!="") cout << tit<<endl;
415  for (unsigned int j=0;j<l.size();++j) printCand(l[j]);
416  cout <<endl;
417 }
418 
419 // *** prints all global lists
421 {
422  printList(mclist, "MC List");
423  printList(recolist, "Reco List");
424 
425  printList(eplus, "e+ List");
426  printList(eminus, "e- List");
427 
428  printList(muplus, "mu+ List");
429  printList(muminus, "mu- List");
430 
431  printList(piplus, "pi+ List");
432  printList(piminus, "pi- List");
433 
434  printList(kplus, "k+ List");
435  printList(kminus, "k- List");
436 
437  printList(pplus, "p+ List");
438  printList(pminus, "p- List");
439 
440  printList(gam, "gamma List");
441 }
442 
443 // *** counts how many combinations (=composite particles in l) fullfill the mass criterion
444 int massCrit(CandList &l, double m, double w)
445 {
446  int cnt = 0;
447  for (unsigned int i=0;i<l.size();++i)
448  {
449  double mass=l[i].P4().M();
450  if (mass>m-w && mass<m+w) cnt++;
451  }
452 
453  return cnt;
454 }
455 
456 // *** fill the mass histos
457 int fillHisto(CandList &l, TH1 *hall, TH1* hbg=0, double m =0, double w=0, TH1* hallsel=0, TH1* hbgsel=0, TH1* hsig=0, TH1* hmm=0, TH2* hmvsmm=0)
458 {
459  int cnt =0;
460  for (unsigned int i=0;i<l.size();++i)
461  {
462  double mass = l[i].P4().M();
463  double miss = (fIni-l[i].P4()).M();
464  hall->Fill(mass);
465 
466  //missing mass
467  if (hmm) hmm->Fill(miss);
468  if (hmvsmm) hmvsmm->Fill(mass, miss);
469 
470  bool mct = l[i].Mct();
471  if (hbg && !mct) hbg->Fill(mass);
472  if (hsig && mct) hsig->Fill(mass);
473  if (mass>m-w && mass<m+w)
474  {
475  cnt++;
476  if (hallsel) hallsel->Fill(mass);
477  if (hbgsel && !mct) hbgsel->Fill(mass);
478  }
479  }
480  return cnt;
481 }
482 
483 int combine(CandList &l1, CandList &l2, CandList &out, int matchPdg=0)
484 {
485  out.clear();
486 
487  bool sameList = (&l1 == &l2);
488 
489  for (unsigned int i=0;i<l1.size();++i)
490  {
491  for (unsigned int j=0;j<l2.size();++j)
492  {
493  if (sameList && j<=i) continue;
494  if (l1[i].Marker() & l2[j].Marker()) continue; // *** overlap
495  //if (l1[i].EvtId()!=l2[j].EvtId()) continue; // *** tracks from different events
496  SimpleCand c;
497  c.SetP4(l1[i].P4()+l2[j].P4());
498  c.SetCharge(l1[i].Charge()+l2[j].Charge());
499  c.SetNFS(l1[i].NFS()+l2[j].NFS());
500  c.SetMarker(l1[i].Marker() | l2[j].Marker());
501 
502  int m1 = l1[i].MotherIdx();
503  int m2 = l2[j].MotherIdx();
504  int ev1 = l1[i].EvtId();
505  int ev2 = l2[j].EvtId();
506 
507  if ( ev1==ev2 && // *** particles are from same event
508  m1>=0 && m1 == m2 && // *** there is a common mother of all particles
509  (l1[i].McPid()||l1[i].Mct()) && // *** all particles match (FS: PID, Comp:mct)
510  (l2[j].McPid()||l2[j].Mct()) && // *** "
511  2==mclist[m1].NDau() && // *** mother has as many daughters as we combine particles
512  matchPdg == mclist[m1].Pid() ) // *** mother has correct PDG
513  {
514  c.SetId(m1);
515  c.SetMcPid();
516  c.SetMct();
517  c.SetMotherIdx(mclist[m1].MotherIdx());
518  }
519 
520  out.push_back(c);
521  }
522  }
523  return out.size();
524 }
525 
526 int combine(CandList &l1, CandList &l2, CandList &l3, CandList &out, int matchPdg=0)
527 {
528  out.clear();
529 
530  bool same12 = (&l1 == &l2);
531  bool same13 = (&l1 == &l3);
532  bool same23 = (&l2 == &l3);
533 
534  for (unsigned int i1=0;i1<l1.size();++i1)
535  {
536  for (unsigned int i2=0;i2<l2.size();++i2)
537  {
538  if (same12 && i2<i1) continue;
539  if (l1[i1].Marker() & l2[i2].Marker()) continue;
540 
541  for (unsigned int i3=0;i3<l3.size();++i3)
542  {
543  if ((same13 && i3<i1) || (same23 && i3<i2)) continue;
544 
545  if ( l1[i1].Marker() & l3[i3].Marker() ) continue;
546  if ( l2[i2].Marker() & l3[i3].Marker() ) continue;
547 
548  //if (l1[i1].EvtId()!=l2[i2].EvtId() ||
549  // l1[i1].EvtId()!=l3[i3].EvtId() ||
550  // l2[i2].EvtId()!=l3[i3].EvtId() ) continue; // *** tracks from different events
551 
552  SimpleCand c;
553  c.SetP4( l1[i1].P4() + l2[i2].P4() + l3[i3].P4());
554  c.SetCharge( l1[i1].Charge() + l2[i2].Charge() + l3[i3].Charge() );
555  c.SetNFS( l1[i1].NFS() + l2[i2].NFS() + l3[i3].NFS());
556  c.SetMarker( l1[i1].Marker() | l2[i2].Marker() | l3[i3].Marker() );
557 
558  int m1 = l1[i1].MotherIdx();
559  int m2 = l2[i2].MotherIdx();
560  int m3 = l3[i3].MotherIdx();
561 
562  int ev1 = l1[i1].EvtId();
563  int ev2 = l2[i2].EvtId();
564  int ev3 = l3[i3].EvtId();
565 
566  if ( ev1==ev2 && ev2==ev3 && // *** particles are from same event
567  m1>=0 && m1 == m2 && m2 == m3 && // *** there is a common mother of all particles
568  (l1[i1].McPid()||l1[i1].Mct()) && // *** all particles match (FS: PID, Comp:mct)
569  (l2[i2].McPid()||l2[i2].Mct()) && // *** "
570  (l3[i3].McPid()||l3[i3].Mct()) && // *** "
571  3==mclist[m1].NDau() && // *** mother has as many daughters as we combine particles
572  matchPdg == mclist[m1].Pid() ) // *** mother has correct PDG
573  {
574  c.SetId(m1);
575  c.SetMcPid();
576  c.SetMct();
577  c.SetMotherIdx(mclist[m1].MotherIdx());
578  }
579 
580  out.push_back(c);
581  }
582  }
583  }
584 
585  return out.size();
586 }
587 
588 void smearMom(SimpleCand &c, double dpr=0.05, double dtht=0.000, double dphi=0.000)
589 {
590  TLorentzVector p4=c.P4();
591 
592  // *** use resolution function
593  if (dpr<0) dpr = dponline(p4);
594 
595  double p = p4.P();
596  double m = p4.M();
597 
598  double dp = p*dpr;
599  double newp = -1;
600 
601  if (dp==0) newp=p;
602  while (newp<0) newp = p4.P() + fRand.Gaus(0.0,dp);
603 
604  p4.SetVectM( p4.Vect().Unit()*newp, m );
605 
606  double newtht = p4.Theta();
607  if (dtht!=0) newtht += fRand.Gaus(0.,dtht);
608  p4.SetTheta(newtht);
609 
610  double newphi = p4.Phi();
611  if (dphi!=0) newphi += fRand.Gaus(0.,dphi);
612  p4.SetPhi(newphi);
613 
614  c.SetP4(p4);
615 
616 }
617 
618 bool isDetected(double eff=0.90)
619 {
620  // very simple efficiency cut; might be dependent on candidate itself later
621  if (fRand.Rndm()<eff) return true;
622  return false;
623 }
624 
625 void makeRecoCands(CandList &mc, CandList &reco, double eff=0.90, double dp=0.05, double dtht=0.001, double dphi=0.001)
626 {
627  reco.clear();
628 
629  for (unsigned int i=0;i<mc.size();i++)
630  {
631  SimpleCand c=mc[i];
632  int pid = abs(c.Pid());
633  if (!isFinalState(pid)) continue;
634 
635  if (pid!=22 && isDetected(eff)) // *** charged track
636  {
637  smearMom(c,dp,dtht,dphi);
638  // set pion mass as default
639  c.SetMass(0.13957);
640  reco.push_back(c);
641  }
642  else if (pid==22 && isDetected(fNeutEff) && c.P4().E()>fNeutThrsh) // *** neutral
643  {
644  smearMom(c,fdE, fNdtht, fNdphi);
645  reco.push_back(c);
646  }
647  }
648 }
649 
650 // **** A simple mass selector
651 void selectMass(CandList &in, CandList &out, double mmean=0, double mw=0)
652 {
653  out.clear();
654 
655  for (unsigned int i=0;i<in.size();++i)
656  if (fabs(in[i].P4().M()-mmean)<mw) out.push_back(in[i]);
657 
658 }
659 
660 // *** perform pid selection according to the pidtable
661 void select(CandList &in, CandList &out, int chrg=0, int pdg=0, double mass=0.139)
662 {
663  out.clear();
664 
665  TParticlePDG *part = pdg_db->GetParticle(pdg);
666 
667  //pdg *= pdg<20 ? -chrg : chrg;
668 
669  if (part) mass = part->Mass();
670  if (mass<0) mass = 0.13957;
671 
672  for (unsigned int i=0;i<in.size();++i)
673  {
674  SimpleCand c=in[i];
675 
676  // select a specific charge in any case
677  if (c.Charge()==chrg)
678  {
679  c.SetMass(mass);
680  if (c.Pid()==pdg) c.SetMcPid();
681  if (pdg==0 || pdg==22) out.push_back(c); // do we want select a pdg code? 0=only charge
682  else
683  {
684  // invalid pdg code?
685  if (pdgidx.find(abs(pdg))==pdgidx.end()) continue;
686 
687  int selidx = pdgidx[abs(pdg)];
688  int partidx = pdgidx[abs(c.Pid())];
689  double prob = pidtab[5*selidx + partidx]/100.;
690  //cout <<pdg <<" "<<c.Pid()<<" : "<<selidx<<" "<<partidx<<" P:"<<prob<<endl;
691 
692  if (fRand.Rndm()<prob) out.push_back(c);
693  }
694  }
695  }
696 }
697 
699 {
700  // *** prepare globally defined PID lists
701  select(l, eplus, 1, -11, mE);
702  select(l, eminus, -1, 11, mE);
703 
704  select(l, muplus, 1, -13, mMu);
705  select(l, muminus, -1, 13, mMu);
706 
707  select(l, piplus, 1, 211, mPi);
708  select(l, piminus, -1, -211, mPi);
709 
710  select(l, kplus, 1, 321, mK);
711  select(l, kminus, -1, -321, mK);
712 
713  select(l, pplus, 1, 2212, mP);
714  select(l, pminus, -1, -2212, mP);
715 
716  select(l, gam, 0, 22, 0);
717 }
718 
719 void makeIni4Vector(TLorentzVector &l, double s)
720 {
721  double X = (s*s-2*mP*mP)/(2*mP);
722  double p = sqrt(X*X - mP*mP);
723  l.SetXYZM(0,0,p,s);
724 }
725 
726 
727 int softtrigger_kin5(TString fsig, int mode=0, double sqrts=3.77, int nev=1000, bool evcut=false, double dp=0.03,
728  double trkeff=100., double pideff=95., double misid =5., double P_mix=0.00 )
729 {
730  fRand.SetSeed();
731 
732  init(pideff, misid);
733 
734  unsigned int i,j,k;
735 
736  makeIni4Vector(fIni, sqrts);
737 
738  //gStyle->SetOptStat(1111);
739 
740  // *** open file and get tree
741  TFile *f=new TFile(fsig,"READ");
742  if (f->IsZombie())
743  {
744  cout <<"File not found:"<<fsig<<endl;
745  return;
746  }
747 
748  bool dpm=fsig.Contains("DPM"); // *** DPM events
749  bool bkg=dpm || fsig.Contains("BKG"); // *** Background events from EvtGen OR DPM events
750  TString treename = dpm ? "data" : "ntp"; // *** set treename
751  TTree *t=(TTree*)f->Get(treename);
752 
753  if (t==0)
754  {
755  cout <<"Tree '"<<treename<<"' not found in '"<<fsig<<"'."<<endl;
756  return;
757  }
758 
759  // *** how many events ?
760  int Nevt = t->GetEntriesFast();
761  if (nev==0) nev = Nevt;
762 
763  double mmean=0, mwsig=0, mwoff=0;
764  int matchPDG=0;
765 
766  // *** choose histogram scale according to resonance
767  high = mmean+2;
768 
769  double min=0, max=5.5;
770 
771  TString name = "ntp_";
772  if (!bkg) name+="sig_mode_";
773  else name+="bkg_mode_";
774  name+=mode; name+=".root";
775 
776  TFile *ggg=new TFile(name, "RECREATE");
777 
778  //branches
779  // 0 1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33
780  TString branches = "sig:mode:npart:nneut:nchrg:ne:nmu:npi:nk:npr:pmax:ptmax:np05:np10:np15:np20:np25:np30:np40:np50:npt05:npt10:npt15:npt20:npt25:npt30:npt40:npt50:sph:apl:pla:thr:sne:scp:";
781 
782  // 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56
783  branches += "fw1:fw2:fw3:pmaxl:ptmaxl:np05l:np10l:np15l:np20l:np25l:np30l:np40l:np50l:ne10:ne20:nmu10:nmu20:npi10:npi20:nk10:nk20:npr10:npr20:";
784 
785  // 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79
786  branches += "sumpt:sumptl:sumptc:sumptcl:sumetn:sumetnl:sumpc:sumpcl:sumen:sumenl:cir:fw4:fw5:pmin:ptmin:pminl:prapmax:sumpt05:sumpt10:sumpc05:sumpc10:sumpc05l:sumpc10l:";
787 
788  // 80 81 82 83 84 85 86 87 88 89
789  branches += "sumen05:sumen10:sumen05l:sumen10l:npi0:nks0:npi005:nks005:npi010:nks010";
790 
791  TNtuple *ntpev = new TNtuple("ntpev","ntpev", branches);
792 
793  Float_t ntpval[100];
794 
795  // *** prepare the tree to loop through
796 
797  Int_t nTrk,Id[fMAX],DF[fMAX],DL[fMAX];
798  Double_t px[fMAX],py[fMAX],pz[fMAX],E[fMAX],m[fMAX];
799  TClonesArray* DpmPart=0;
800 
801  if (dpm)
802  {
803  DpmPart=new TClonesArray("TParticle",100);
804  t->SetBranchAddress("Npart",&nTrk);
805  t->SetBranchAddress("Particles", &DpmPart);
806  }
807  else
808  {
809  t->SetBranchAddress("nTrk",&nTrk);
810  t->SetBranchAddress("Id",Id);
811  t->SetBranchAddress("DF",DF);
812  t->SetBranchAddress("DL",DL);
813  t->SetBranchAddress("px",px);
814  t->SetBranchAddress("py",py);
815  t->SetBranchAddress("pz",pz);
816  t->SetBranchAddress("E",E);
817  t->SetBranchAddress("m",m);
818  }
819 
820 
821  // *** the event loop
822 
823  int mixevts=1;
824  int currmix=0;
825 
826  // how many events are tagged by individual algo?
827  bool channelselect[6];
828  double channelcount[6];
829  TString channelname[6]={"J/psi","D+-","Ds","Phi","D0","Lamc"};
830 
831  for (i=0;i<6;++i) channelcount[i]=0.;
832 
833  if (fRand.Rndm()<P_mix) mixevts=2;
834 
835  for (i=0;i<nev;++i)
836  {
837  t->GetEntry(i);
838 
839  if (i%5000==0) cout <<"ev "<<i<<endl;
840 
841  // *** read in from file the mc candidate
842  if (0==currmix)
843  {
844  mclist.clear(); // *** contains all MC cands (also non-final-states)
845  mcFSlist.clear(); // *** contains only final state MC cands
846  }
847 
848  bool good = true;
849 
850  if (dpm)
851  {
852  for (j=0;j<nTrk;++j)
853  {
854  TParticle *pt=(TParticle*)DpmPart->At(j);
855  Id[j] = pt->GetPdgCode();
856  int pdg = abs(Id[j]);
857  px[j] = pt->Px();
858  py[j] = pt->Py();
859  pz[j] = pt->Pz();
860  E[j] = pt->Energy();
861  m[j] = pt->GetMass();
862  DF[j] = 0;
863  DL[j] = -1;
864  }
865  }
866 
867  for (j=0;j<nTrk;++j)
868  {
869  int chrg = 0;
870  int pdg = abs(Id[j]);
871  TParticlePDG *p=pdg_db->GetParticle(Id[j]);
872  if (p) chrg = p->Charge();
873  if (abs(chrg)>1) chrg/=3;
874 
875  SimpleCand c(px[j], py[j], pz[j], E[j], Id[j], j, chrg);
876  c.SetDau(DF[j], DL[j]);
877 
878  for (k=0;k<j;++k) if (j>=DF[k] && j<=DL[k])
879  {
880  c.SetMotherIdx(k);
881  c.SetNSiblings(DL[k]-DF[k]+1);
882  }
883 
884  c.SetEvtId(i);
885 
886  mclist.push_back(c);
887 
888  // final state particle?
889  if (isFinalState(pdg)) mcFSlist.push_back(c);
890  }
891 
892  if (++currmix<mixevts) continue;
893 
894  currmix=0;
895 
896  // *** the boost of the system
897  TVector3 boost = fIni.BoostVector();
898 
899  // *** prepare the toy reco: smear momentum and apply track efficiency (both flat at the moment)
900  makeRecoCands(mcFSlist, recolist, trkeff/100., dp );
901 
902 /* int nneut = 0;
903  int nchrg = 0;
904 
905  for (j=0;j<recolist.size();++j)
906  {
907  int chrg = recolist[j].Charge();
908  if (abs(chrg)==0) nneut++; else nchrg++;
909  }
910 */
911  boostList(recolist, recocmslist, boost);
912 
913  // *** prepare all pid list according to the pidtable probabilities
915 
916  EventShape evs(recolist, fIni);
917 
918  // *** the lists we need for combinatorics
919  CandList comb1;
920  CandList LPi0, LKs; // pi0 -> gg (mass cut), KS -> pi+ pi- (mass cut)
921 
922  // *** do combinatorics (pi0 -> gg)
923  combine(gam, gam, comb1, fPi0Pdg);
924  //fillHisto(comb1, h_pi0, 0, fPi0Mass, fPi0win, h_pi0sel, 0, h_pi0sig );
925  selectMass(comb1, LPi0, fPi0Mass, fPi0win);
926 
927  // *** do combinatorics (pi0 -> gg)
928  combine(piplus, piminus, comb1, fKsPdg);
929  //fillHisto(comb1, h_ks, 0, fKsMass, fKswin, h_kssel, 0, h_kssig );
930  selectMass(comb1, LKs, fKsMass, fKswin);
931 
932  // *** Make Event Selection
933 // if (evcut)
934 // {
935 // EventShape evshape(recolist, fIni);
936 // if (evshape.ChrgPtSumLab()<1.5) continue;
937 // }
938 
939  // 0 1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33
940  // "sig:mode:npart:nneut:nchrg:ne:nmu:npi:nk:npr:pmax:ptmax:np05:np10:np15:np20:np25:np30:np40:np50:npt05:npt10:npt15:npt20:npt25:npt30:npt40:npt50:sph:apl:pla:thr:sne:scp";
941 
942  // 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56
943  // "fw1:fw2:fw3:pmaxl:ptmaxl:np05l:np10l:np15l:np20l:np25l:np30l:np40l:np50l:ne10:ne20:nmu10:nmu20:npi10:npi20:nk10:nk20:npr10:npr20");
944 
945  // 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79
946  // "sumpt:sumptl:sumptc:sumptcl:sumetn:sumetnl:sumpc:sumpcl:sumen:sumenl:cir:fw4:fw5:pmin:ptmin:pminl:prapmax:sumpt05:sumpt10:sumpc05:sumpc10:sumpc05l:sumpc10l:";
947 
948  // 80 81 82 83 84 85 86 87 88 89
949  // "sumen05:sumen10:sumen05l:sumen10l:npi0:nks0:npi005:nks005:npi010:nks010";
950 
951  ntpval[0] = !bkg; // signal or background events
952  ntpval[1] = mode; // decay mode (parameter of script - set by hand)
953 
954  ntpval[2] = evs.NParticles(); // #particles (multiplicity)
955  ntpval[3] = evs.NNeutral(); // #neutrals cands
956  ntpval[4] = evs.NCharged(); // #charged cands
957 
958  ntpval[5] = eplus.size()+eminus.size(); // #electron candidates
959  ntpval[6] = muplus.size()+muminus.size(); // #muon candidates
960  ntpval[7] = piplus.size()+piminus.size(); // #pion candidates
961  ntpval[8] = kplus.size()+kminus.size(); // #kaon candidates
962  ntpval[9] = pplus.size()+pminus.size(); // #proton candidates
963  ntpval[84] = LPi0.size(); // #pi0 candidates
964  ntpval[85] = LKs.size(); // #KS0 candidates
965 
966  ntpval[37] = evs.PmaxLab(); // maximum momentum (lab)
967  ntpval[38] = evs.Ptmax(); // maximum pt (lab)
968  ntpval[72] = evs.PminLab(); // minimum momentum (lab)
969  ntpval[71] = evs.Ptmin(); // minimum pt (lab)
970  ntpval[73] = evs.PRapmax(); // maximum pseudo rapidity (lab)
971 
972  ntpval[39] = evs.MultPminLab(0.5); // #particles with p>0.5 GeV (lab)
973  ntpval[40] = evs.MultPminLab(1.0); // #particles with p>1.0 GeV (lab)
974  ntpval[41] = evs.MultPminLab(1.5); // #particles with p>1.5 GeV (lab)
975  ntpval[42] = evs.MultPminLab(2.0); // #particles with p>2.0 GeV (lab)
976  ntpval[43] = evs.MultPminLab(2.5); // #particles with p>2.5 GeV (lab)
977  ntpval[44] = evs.MultPminLab(3.0); // #particles with p>3.0 GeV (lab)
978  ntpval[45] = evs.MultPminLab(4.0); // #particles with p>4.0 GeV (lab)
979  ntpval[46] = evs.MultPminLab(5.0); // #particles with p>5.0 GeV (lab)
980 
981  ntpval[58] = evs.PtSumLab(); // sum pt of all particles (lab)
982  ntpval[60] = evs.ChrgPtSumLab(); // sum pt of chrg particles (lab)
983  ntpval[62] = evs.NeutEtSumLab(); // sum Et of neut particles (lab)
984  ntpval[64] = evs.ChrgPSumLab(); // sum p of chrg particles (lab)
985  ntpval[66] = evs.NeutESumLab(); // sum E of neut particles (lab)
986 
987  ntpval[78] = evs.SumChrgPminLab(0.5); // sum p of charged particles with p>0.5 (lab)
988  ntpval[79] = evs.SumChrgPminLab(1.0); // sum p of charged particles with p>1.0 (lab)
989 
990  ntpval[82] = evs.SumNeutEminLab(0.5); // sum E of neutral particles with E>0.5 (lab)
991  ntpval[83] = evs.SumNeutEminLab(1.0); // sum E of neutral particles with E>1.0 (lab)
992 
993  ntpval[58] = evs.PtSumLab(); // sum pt of all particles (lab)
994  ntpval[74] = evs.SumPtminLab(0.5); // sum pt of all particles with pt>0.5 GeV
995  ntpval[75] = evs.SumPtminLab(1.0); // sum pt of all particles with pt>1.0 GeV
996 
997  // all following values are in center-of-momentum system
998  ntpval[10] = evs.PmaxCms(); // maximum momentum
999  ntpval[11] = evs.Ptmax(); // maximum pt
1000  ntpval[70] = evs.PminCms(); // minimum momentum
1001 
1002  ntpval[12] = evs.MultPminCms(0.5); // #particles with p>0.5 GeV
1003  ntpval[13] = evs.MultPminCms(1.0); // #particles with p>1.0 GeV
1004  ntpval[14] = evs.MultPminCms(1.5); // #particles with p>1.5 GeV
1005  ntpval[15] = evs.MultPminCms(2.0); // #particles with p>2.0 GeV
1006  ntpval[16] = evs.MultPminCms(2.5); // #particles with p>2.5 GeV
1007  ntpval[17] = evs.MultPminCms(3.0); // #particles with p>3.0 GeV
1008  ntpval[18] = evs.MultPminCms(4.0); // #particles with p>4.0 GeV
1009  ntpval[19] = evs.MultPminCms(5.0); // #particles with p>5.0 GeV
1010 
1011  ntpval[57] = evs.PtSumCms(); // sum pt of all particles
1012  ntpval[59] = evs.ChrgPtSumCms(); // sum pt of chrg particles
1013  ntpval[61] = evs.NeutEtSumCms(); // sum Et of neut particles
1014  ntpval[63] = evs.ChrgPSumCms(); // sum p of chrg particles
1015  ntpval[65] = evs.NeutESumCms(); // sum E of neut particles
1016 
1017  ntpval[76] = evs.SumChrgPminCms(0.5); // sum p of charged particles with p>0.5 (cms)
1018  ntpval[77] = evs.SumChrgPminCms(1.0); // sum p of charged particles with p>1.0 (cms)
1019 
1020  ntpval[80] = evs.SumNeutEminCms(0.5); // sum E of neutral particles with E>0.5 (cms)
1021  ntpval[81] = evs.SumNeutEminCms(1.0); // sum E of neutral particles with E>1.0 (cms)
1022 
1023  ntpval[20] = 0.0; // #particles with pt>0.5 GeV
1024  ntpval[21] = 0.0; // #particles with pt>1.0 GeV
1025  ntpval[22] = 0.0; // #particles with pt>1.5 GeV
1026  ntpval[23] = 0.0; // #particles with pt>2.0 GeV
1027  ntpval[24] = 0.0; // #particles with pt>2.5 GeV
1028  ntpval[25] = 0.0; // #particles with pt>3.0 GeV
1029  ntpval[26] = 0.0; // #particles with pt>4.0 GeV
1030  ntpval[27] = 0.0; // #particles with pt>5.0 GeV
1031 
1032  ntpval[47] = 0.0; // #electrons with p>1.0
1033  ntpval[48] = 0.0; // #electrons with p>2.0
1034 
1035  ntpval[49] = 0.0; // #muons with p>1.0
1036  ntpval[50] = 0.0; // #muons with p>2.0
1037 
1038  ntpval[51] = 0.0; // #pions with p>1.0
1039  ntpval[52] = 0.0; // #pions with p>2.0
1040 
1041  ntpval[53] = 0.0; // #kaons with p>1.0
1042  ntpval[54] = 0.0; // #kaons with p>2.0
1043 
1044  ntpval[55] = 0.0; // #protons with p>1.0
1045  ntpval[56] = 0.0; // #protons with p>2.0
1046 
1047  ntpval[86] = 0.0; // #pi0s with p>0.5
1048  ntpval[88] = 0.0; // #pi0s with p>1.0
1049 
1050  ntpval[87] = 0.0; // #KSs with p>0.5
1051  ntpval[89] = 0.0; // #KSs with p>1.0
1052 
1053  //compSphericity(recocmslist);
1054 
1055  //TVector3 thrAxis;
1056  //_thr = compThrust(recocmslist, thrAxis);
1057 
1058  // quantities based on momentum tensor
1059  ntpval[28] = evs.Sphericity(); // sphericity = 3/2 * (lam_2+lam_3)
1060  ntpval[29] = evs.Aplanarity(); // aplanarity = 3/2 * lam_3
1061  ntpval[30] = evs.Planarity(); // planarity = lam_2 - lam_3
1062  ntpval[67] = evs.Circularity(); // circularity = 2 * min(lam_1, lam_2)/(lam_1 + lam_2)
1063  ntpval[31] = evs.Thrust(); // thrust
1064 
1065  //compFoxWolfMoms(recocmslist);
1066 
1067  ntpval[34] = evs.FoxWolfMomR(1); // Fox Wolfram Moment R1 = H1/H0
1068  ntpval[35] = evs.FoxWolfMomR(2); // Fox Wolfram Moment R2 = H2/H0
1069  ntpval[36] = evs.FoxWolfMomR(3); // Fox Wolfram Moment R3 = H3/H0
1070  ntpval[68] = evs.FoxWolfMomR(4); // Fox Wolfram Moment R4 = H4/H0
1071  ntpval[69] = evs.FoxWolfMomR(5); // Fox Wolfram Moment R5 = H5/H0
1072 
1073  // count electrons
1074  for (j=0; j<eplus.size(); ++j)
1075  {
1076  if (eplus[j].P4().P()>1.0) ntpval[47]+=1.0;
1077  if (eplus[j].P4().P()>2.0) ntpval[48]+=1.0;
1078  }
1079  for (j=0; j<eminus.size(); ++j)
1080  {
1081  if (eminus[j].P4().P()>1.0) ntpval[47]+=1.0;
1082  if (eminus[j].P4().P()>2.0) ntpval[48]+=1.0;
1083  }
1084 
1085  // count muons
1086  for (j=0; j<muplus.size(); ++j)
1087  {
1088  if (muplus[j].P4().P()>1.0) ntpval[49]+=1.0;
1089  if (muplus[j].P4().P()>2.0) ntpval[50]+=1.0;
1090  }
1091  for (j=0; j<muminus.size(); ++j)
1092  {
1093  if (muminus[j].P4().P()>1.0) ntpval[49]+=1.0;
1094  if (muminus[j].P4().P()>2.0) ntpval[50]+=1.0;
1095  }
1096 
1097  // count pions
1098  for (j=0; j<piplus.size(); ++j)
1099  {
1100  if (piplus[j].P4().P()>1.0) ntpval[51]+=1.0;
1101  if (piplus[j].P4().P()>2.0) ntpval[52]+=1.0;
1102  }
1103  for (j=0; j<piminus.size(); ++j)
1104  {
1105  if (piminus[j].P4().P()>1.0) ntpval[51]+=1.0;
1106  if (piminus[j].P4().P()>2.0) ntpval[52]+=1.0;
1107  }
1108 
1109  // count kaons
1110  for (j=0; j<kplus.size(); ++j)
1111  {
1112  if (kplus[j].P4().P()>1.0) ntpval[53]+=1.0;
1113  if (kplus[j].P4().P()>2.0) ntpval[54]+=1.0;
1114  }
1115  for (j=0; j<kminus.size(); ++j)
1116  {
1117  if (kminus[j].P4().P()>1.0) ntpval[53]+=1.0;
1118  if (kminus[j].P4().P()>2.0) ntpval[54]+=1.0;
1119  }
1120 
1121  // count protons
1122  for (j=0; j<pplus.size(); ++j)
1123  {
1124  if (pplus[j].P4().P()>1.0) ntpval[55]+=1.0;
1125  if (pplus[j].P4().P()>2.0) ntpval[56]+=1.0;
1126  }
1127  for (j=0; j<pminus.size(); ++j)
1128  {
1129  if (pminus[j].P4().P()>1.0) ntpval[55]+=1.0;
1130  if (pminus[j].P4().P()>2.0) ntpval[56]+=1.0;
1131  }
1132 
1133  // count pi0
1134  for (j=0; j<LPi0.size(); ++j)
1135  {
1136  if (LPi0[j].P4().P()>0.5) ntpval[86]+=1.0;
1137  if (LPi0[j].P4().P()>1.0) ntpval[88]+=1.0;
1138  }
1139  // count KS
1140  for (j=0; j<LKs.size(); ++j)
1141  {
1142  if (LKs[j].P4().P()>0.5) ntpval[87]+=1.0;
1143  if (LKs[j].P4().P()>1.0) ntpval[89]+=1.0;
1144  }
1145 
1146 
1147 
1148  for (j=0;j<recolist.size();++j)
1149  {
1150  // *** candidate in lab system
1151  SimpleCand cnd = recolist[j];
1152  TLorentzVector lv = cnd.P4();
1153 
1154  // *** candidate in cms system
1155  SimpleCand cndcms = recocmslist[j];
1156  TLorentzVector lvcm = cndcms.P4();
1157 
1158  // compute in LAB
1159 
1160 
1161  if (lvcm.Pt()>0.5) ntpval[20]+=1.0;
1162  if (lvcm.Pt()>1.0) ntpval[21]+=1.0;
1163  if (lvcm.Pt()>1.5) ntpval[22]+=1.0;
1164  if (lvcm.Pt()>2.0) ntpval[23]+=1.0;
1165  if (lvcm.Pt()>2.5) ntpval[24]+=1.0;
1166  if (lvcm.Pt()>3.0) ntpval[25]+=1.0;
1167  if (lvcm.Pt()>4.0) ntpval[26]+=1.0;
1168  if (lvcm.Pt()>5.0) ntpval[27]+=1.0;
1169 
1170  }
1171 
1172  ntpev->Fill(ntpval);
1173 
1174  mixevts=1;
1175  if (fRand.Rndm()<P_mix) mixevts=2;
1176  }
1177 
1178 
1179  ggg->cd();
1180  ntpev->Write();
1181  ggg->Close();
1182 
1183  return 0;
1184 }
double ChrgPtSumCms() const
Definition: EventShape.h:39
double fDwin
int fLamPdg
const double fNdtht
CandList pminus
void select(CandList &in, CandList &out, int chrg=0, int pdg=0, double mass=0.139)
int fPhiPdg
CandList mclist
Double_t p
Definition: anasim.C:58
printf("RealTime=%f seconds, CpuTime=%f seconds\n", rtime, ctime)
void printList(CandList &l, TString tit="")
double fKswin
friend F32vec4 cos(const F32vec4 &a)
Definition: P4_F32vec4.h:112
double SumNeutEminLab(double emin)
Definition: EventShape.h:746
int NNeutral() const
Definition: EventShape.h:17
double mK
void SetP4(TLorentzVector lv)
Definition: SimpleCand.C:20
void SetEvtId(int ev)
Definition: SimpleCand.C:35
double _thr
int fJpsiPdg
double NeutESumLab() const
Definition: EventShape.h:31
const double fNeutEff
double mP
void makePidSelection(CandList &l)
double NeutEtSumCms() const
Definition: EventShape.h:37
Int_t i
Definition: run_full.C:25
__m128 m
Definition: P4_F32vec4.h:28
int low
Definition: anaMvdDigi.C:51
int NCharged() const
Definition: EventShape.h:16
double fPi0win
TRandom3 fRand
Definition: full_core_ntp.C:28
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
double fDsMass
double _sph
TDatabasePDG * pdg_db
int NParticles() const
Definition: EventShape.h:15
TLorentzVector s
Definition: Pnd2DStar.C:50
int n
int pid()
const int fMAX
Definition: full_core_ntp.C:32
CandList kminus
int MultPminLab(double pmin)
Definition: EventShape.h:471
double Aplanarity()
Definition: EventShape.h:299
double mE
CandList recocmslist
double fJ3pioff
CandList piplus
double fJsig
double fLamcsig
int NFS() const
Definition: SimpleCand.C:47
double Ptmin() const
Definition: EventShape.h:25
friend F32vec4 max(const F32vec4 &a, const F32vec4 &b)
Definition: P4_F32vec4.h:26
void makePidTable(double pideff, double misid)
bool Mct() const
Definition: SimpleCand.C:48
double fLamwin
int fEtacPdg
void makeRecoCands(CandList &mc, CandList &reco, double eff=0.90, double dp=0.05, double dtht=0.001, double dphi=0.001)
void SetMcPid(bool mct=true)
Definition: SimpleCand.C:30
int combine(CandList &l1, CandList &l2, CandList &out, int matchPdg=0)
CandList eplus
int fillHisto(CandList &l, TH1 *hall, TH1 *hbg=0, double m=0, double w=0, TH1 *hallsel=0, TH1 *hbgsel=0, TH1 *hsig=0, TH1 *hmm=0, TH2 *hmvsmm=0)
int Pid() const
Definition: SimpleCand.C:40
std::vector< SimpleCand > CandList
Definition: SimpleCand.C:76
double SumNeutEminCms(double emin)
Definition: EventShape.h:768
double fDswin
__m128 v
Definition: P4_F32vec4.h:4
std::map< int, int > pdgidx
int fDPdg
double PmaxLab() const
Definition: EventShape.h:20
double dponline(TLorentzVector &v)
double fDMass
double compSphericity(CandList &list)
void SetMarker(unsigned int i)
Definition: SimpleCand.C:33
int fLamcPdg
double PtSumLab() const
Definition: EventShape.h:29
double NeutEtSumLab() const
Definition: EventShape.h:30
static int init
Definition: ranlxd.cxx:374
void SetCharge(int ch)
Definition: SimpleCand.C:27
double fKsMass
void SetMct(bool mct=true)
Definition: SimpleCand.C:29
TString pt(TString pts, TString exts="px py pz")
Definition: invexp.C:133
int NSiblings() const
Definition: SimpleCand.C:46
CandList muplus
bool isFinalState(int n)
Int_t mode
Definition: autocutx.C:47
double fJsig3pi
const double fdE
int massCrit(CandList &l, double m, double w)
TString m2(TString pts, TString exts="e px py pz")
Definition: invexp.C:117
double SumChrgPminLab(double pmin)
Definition: EventShape.h:792
double _foxwolfR[_nfwmom]
double _apl
double pidreal[25]
double compThrust(CandList &list, TVector3 &axis)
TVector3 _secAx
double PRapmax() const
Definition: EventShape.h:26
void SetMass(double mass)
Definition: SimpleCand.C:22
Double_t
double eps(TVector3 v1, TVector3 v2)
double fLamsig
double pidtab[25]
void SetNSiblings(int n)
Definition: SimpleCand.C:31
int Id() const
Definition: SimpleCand.C:42
double NeutESumCms() const
Definition: EventShape.h:38
void SetId(int id)
Definition: SimpleCand.C:26
double fLamcMass
double fPhiwin
double fLamcwin
double fPhisig
double PminCms() const
Definition: EventShape.h:23
double fEtacMass
TFile * f
Definition: bump_analys.C:12
TVector3 _sphAx
double fPhiMass
double mPi0
CandList pplus
int nbins
Definition: full_core_ntp.C:33
double Sphericity()
Definition: EventShape.h:281
double fDsig
friend F32vec4 min(const F32vec4 &a, const F32vec4 &b)
Definition: P4_F32vec4.h:25
double Circularity()
Definition: EventShape.h:308
axis
Definition: PndRadMapPlane.h:9
int MultPminCms(double pmin)
Definition: EventShape.h:493
double Planarity()
Definition: EventShape.h:290
TLorentzVector P4() const
Definition: SimpleCand.C:39
CandList recolist
double fEtacsig
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
TFile * out
Definition: reco_muo.C:20
const double fNeutThrsh
void smearMom(SimpleCand &c, double dpr=0.05, double dtht=0.000, double dphi=0.000)
CandList gam
int fD0Pdg
TString name
double fEtacwin
void boostList(CandList &list, CandList &boostlist, TVector3 boost)
int softtrigger_kin5(TString fsig, int mode=0, double sqrts=3.77, int nev=1000, bool evcut=false, double dp=0.03, double trkeff=100., double pideff=95., double misid=5., double P_mix=0.00)
double _pla
TPad * p2
Definition: hist-t7.C:117
double ChrgPSumCms() const
Definition: EventShape.h:40
double X
Definition: anaLmdDigi.C:68
TVector3 _thrAx
double Thrust()
Definition: EventShape.h:324
int fPi0Pdg
Double_t x
double FoxWolfMomR(int order)
Definition: EventShape.h:462
int reco()
double PmaxCms() const
Definition: EventShape.h:21
double fPi0Mass
CandList piminus
TVector3 v1
Definition: bump_analys.C:40
Int_t cnt
Definition: hist-t7.C:106
CandList mcFSlist
double high
Definition: full_core_ntp.C:35
TVector3 v2
Definition: bump_analys.C:40
double fJMass
TPad * p1
Definition: hist-t7.C:116
double fLamMass
void printCand(TLorentzVector l, TVector3 p)
double PtSumCms() const
Definition: EventShape.h:36
TF1 * f_res
TLorentzVector fIni
Definition: full_core_ntp.C:29
TTree * t
Definition: bump_analys.C:13
void selectMass(CandList &in, CandList &out, double mmean=0, double mw=0)
double PminLab() const
Definition: EventShape.h:22
int fDsPdg
double mPi
double SumChrgPminCms(double pmin)
Definition: EventShape.h:814
void printAllLists()
CandList muminus
int fKsPdg
double ChrgPSumLab() const
Definition: EventShape.h:33
double fD0sig
const double fNdphi
double mMu
void SetMotherIdx(int idx)
Definition: SimpleCand.C:28
bool isDetected(double eff=0.90)
double fDssig
void SetNFS(int n)
Definition: SimpleCand.C:32
double ChrgPtSumLab() const
Definition: EventShape.h:32
double fD0Mass
void compFoxWolfMoms(CandList &list)
PndMCTrack * mct
Definition: hist-t7.C:147
void SetDau(int n, int m=0)
Definition: SimpleCand.C:34
CandList kplus
double fJwin
double fD0win
int MotherIdx() const
Definition: SimpleCand.C:45
const int _nfwmom
TMatrixT< double > TMatrixD
Definition: PndLmdDim.h:52
int Charge() const
Definition: SimpleCand.C:44
double pz[39]
Definition: pipisigmas.h:14
bool McPid() const
Definition: SimpleCand.C:49
void makeIni4Vector(TLorentzVector &l, double s)
double Legendre(int l, double x)
double Ptmax() const
Definition: EventShape.h:24
CandList eminus
double SumPtminLab(double ptmin)
Definition: EventShape.h:700