FairRoot/PandaRoot
Public Member Functions | Private Member Functions | Private Attributes | List of all members
PndTutAnaTaskD0 Class Reference

#include <PndTutAnaTaskD0.h>

Inheritance diagram for PndTutAnaTaskD0:

Public Member Functions

 PndTutAnaTaskD0 (double pbarmom)
 
 ~PndTutAnaTaskD0 ()
 
virtual InitStatus Init ()
 
virtual void Exec (Option_t *opt)
 
virtual void Finish ()
 

Private Member Functions

virtual void SetParContainers ()
 
 ClassDef (PndTutAnaTaskD0, 1)
 

Private Attributes

int fEvtCount
 
TLorentzVector fIni
 
PndAnalysisfAnalysis
 
RhoTuplend0
 
RhoMassParticleSelectorfD0Sel
 
RhoVtxPocafVtxPoca
 

Detailed Description

Definition at line 22 of file PndTutAnaTaskD0.h.

Constructor & Destructor Documentation

PndTutAnaTaskD0::PndTutAnaTaskD0 ( double  pbarmom)

Definition at line 63 of file PndTutAnaTaskD0.cxx.

References fIni, mp, and sqrt().

63  :
64  FairTask("Panda Tutorial Analysis Task")
65 {
66  double mp=0.938272;
67  fIni.SetXYZT(0,0,pbarmom, sqrt(pbarmom*pbarmom+mp*mp)+mp);
68 }
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
TLorentzVector fIni
static const double mp
Definition: mzparameters.h:11
PndTutAnaTaskD0::~PndTutAnaTaskD0 ( )

Definition at line 73 of file PndTutAnaTaskD0.cxx.

73 { }

Member Function Documentation

PndTutAnaTaskD0::ClassDef ( PndTutAnaTaskD0  ,
 
)
private
void PndTutAnaTaskD0::Exec ( Option_t *  opt)
virtual

Definition at line 109 of file PndTutAnaTaskD0.cxx.

References RhoCandList::Append(), RhoCandidate::Charge(), RhoTuple::Column(), RhoCandList::Combine(), d0, RhoTuple::DumpData(), fAnalysis, fD0Sel, fEvtCount, PndAnalysis::FillList(), fIni, fVtxPoca, PndAnalysis::GetEvent(), RhoCandList::GetLength(), RhoCandidate::GetPidInfo(), RhoVtxPoca::GetPocaVtx(), RhoCandidate::GetRecoCandidate(), i, mct, PndAnalysis::McTruthMatch(), nd0, P, RhoCandidate::P(), RhoCandidate::P3(), RhoCandidate::P4(), RhoCandList::Select(), RhoCandList::SetType(), and TString.

110 {
111  // *** some variables
112  int i=0,j=0, k=0, l=0;
113 
114  // *** necessary to read the next event
115  fAnalysis->GetEvent();
116 
117  // *** print event counter
118  if (!(++fEvtCount%100)) cout << "evt "<<fEvtCount<<endl;
119 
120  // *******
121  // ******* PUT ANALYSIS CODE HERE
122  // *******
123  RhoCandList kp, km, pip, pim, d0, d0b;
124  TString pidsel = "PidAlgoEmcBayes;PidAlgoDrc;PidAlgoDisc;PidAlgoStt;PidAlgoMdtHardCuts";
125 
126  fAnalysis->FillList( kp, "KaonAllPlus", pidsel);
127  fAnalysis->FillList( km, "KaonAllMinus", pidsel);
128  fAnalysis->FillList( pip, "PionAllPlus", pidsel);
129  fAnalysis->FillList( pim, "PionAllMinus", pidsel);
130 
131  d0.Combine(km, pip);
132  d0.SetType(421);
133  d0b.Combine(kp, pim);
134  d0b.SetType(-421);
135  d0.Append(d0b);
136 
137  d0.Select(fD0Sel);
138 
139  for (j=0;j<d0.GetLength();++j)
140  {
141  RhoCandidate *dka = d0[j]->Daughter(0);
142  RhoCandidate *dpi = d0[j]->Daughter(1);
143 
146 
147  TLorentzVector d0cm = d0[j]->P4();
148  TLorentzVector kacm = dka->P4();
149  TLorentzVector picm = dpi->P4();
150  d0cm.Boost(-fIni.BoostVector());
151  kacm.Boost(-fIni.BoostVector());
152  picm.Boost(-fIni.BoostVector());
153 
154  bool mct = fAnalysis->McTruthMatch(d0[j]);
155 
156  TVector3 vtx;
157  double qavtx = fVtxPoca->GetPocaVtx(vtx, d0[j]);
158 
159  nd0->Column("ev", (Float_t) i);
160  nd0->Column("cand", (Float_t) j);
161  nd0->Column("pdg", (Float_t) d0[j]->PdgCode());
162  nd0->Column("mct", (Float_t) mct);
163 
164  nd0->Column("d0m", (Float_t) d0[j]->M());
165  nd0->Column("d0p", (Float_t) d0[j]->P());
166  nd0->Column("d0tht", (Float_t) d0[j]->P3().Theta());
167  nd0->Column("d0phi", (Float_t) d0[j]->P3().Phi());
168  nd0->Column("d0pt", (Float_t) d0[j]->P3().Pt());
169  nd0->Column("d0pcm", (Float_t) d0cm.P());
170  nd0->Column("d0thtcm", (Float_t) d0cm.Theta());
171  nd0->Column("d0miss", (Float_t) (fIni-(d0[j]->P4())).M());
172 
173  nd0->Column("d0vx", (Float_t) vtx.X());
174  nd0->Column("d0vy", (Float_t) vtx.Y());
175  nd0->Column("d0vz", (Float_t) vtx.Z());
176  nd0->Column("d0vqa", (Float_t) qavtx);
177 
178  nd0->Column("kch", (Float_t) dka->Charge());
179  nd0->Column("kp", (Float_t) dka->P());
180  nd0->Column("ktht", (Float_t) dka->P3().Theta());
181  nd0->Column("kphi", (Float_t) dka->P3().Phi());
182  nd0->Column("kpt", (Float_t) dka->P3().Pt());
183  nd0->Column("kpid", (Float_t) dka->GetPidInfo(3));
184  nd0->Column("kapcm", (Float_t) kacm.P());
185  nd0->Column("kathtcm", (Float_t) kacm.Theta());
186 
187  nd0->Column("pich", (Float_t) dpi->Charge());
188  nd0->Column("pip", (Float_t) dpi->P());
189  nd0->Column("pitht", (Float_t) dpi->P3().Theta());
190  nd0->Column("piphi", (Float_t) dpi->P3().Phi());
191  nd0->Column("pit", (Float_t) dpi->P3().Pt());
192  nd0->Column("pipid", (Float_t) dpi->GetPidInfo(2));
193  nd0->Column("pipcm", (Float_t) picm.P());
194  nd0->Column("pithtcm", (Float_t) picm.Theta());
195 
196  nd0->DumpData();
197  }
198 
199 }
void Append(const RhoCandidate *c)
Definition: RhoCandList.h:52
RhoVtxPoca * fVtxPoca
Int_t i
Definition: run_full.C:25
TLorentzVector fIni
RhoMassParticleSelector * fD0Sel
Int_t GetLength() const
Definition: RhoCandList.h:46
Bool_t FillList(RhoCandList &l, TString listkey="All", TString pidTcaNames="", int trackHypothesis=-1)
PndPidCandidate * GetRecoCandidate() const
Definition: RhoCandidate.h:376
Double_t d0
Definition: checkhelixhit.C:59
void Combine(RhoCandList &l1, RhoCandList &l2)
Double_t GetPocaVtx(TVector3 &vertex, RhoCandidate *composite)
Definition: RhoVtxPoca.cxx:27
void Select(RhoParticleSelectorBase *pidmgr)
TLorentzVector P4() const
Definition: RhoCandidate.h:195
void Column(const char *label, Bool_t value, Bool_t defval=0, const char *block=0)
Definition: RhoTuple.cxx:56
void SetType(const TParticlePDG *pdt, int start=0)
Double_t P() const
void DumpData()
Definition: RhoTuple.cxx:391
GeV c P
Double_t Charge() const
Definition: RhoCandidate.h:184
Bool_t McTruthMatch(RhoCandidate *cand, Int_t level=2, bool verbose=false)
Int_t GetEvent(Int_t n=-1)
PndMCTrack * mct
Definition: hist-t7.C:147
double GetPidInfo(int hypo)
TVector3 P3() const
Definition: RhoCandidate.h:199
PndAnalysis * fAnalysis
void PndTutAnaTaskD0::Finish ( )
virtual

Definition at line 202 of file PndTutAnaTaskD0.cxx.

References RhoTuple::GetInternalTree(), and nd0.

203 {
204 
205  // *******
206  // ******* STORE YOUR HISTOS AND TUPLES
207  // *******
208  nd0->GetInternalTree()->Write();
209 }
TTree * GetInternalTree()
Definition: RhoTuple.h:207
InitStatus PndTutAnaTaskD0::Init ( )
virtual

Definition at line 78 of file PndTutAnaTaskD0.cxx.

References fAnalysis, fD0Sel, fEvtCount, fVtxPoca, and nd0.

79 {
80  // *** initialize PndAnalysis object
81  fAnalysis = new PndAnalysis();
82 
83  // *** reset the event counter
84  fEvtCount = 0;
85 
86  // *******
87  // ******* PREPARE/CREATE THE STUFF YOU NEED
88  // *******
89  fD0Sel = new RhoMassParticleSelector("D0Sel", 1.864, 1.0);
90  nd0 = new RhoTuple("nd0","my D0 tuple");
91  fVtxPoca = new RhoVtxPoca();
92 
93  return kSUCCESS;
94 }
RhoVtxPoca * fVtxPoca
RhoMassParticleSelector * fD0Sel
PndAnalysis * fAnalysis
void PndTutAnaTaskD0::SetParContainers ( )
privatevirtual

Definition at line 98 of file PndTutAnaTaskD0.cxx.

References run.

99 {
100  // *** Get run and runtime database
101  FairRun* run = FairRun::Instance();
102  if ( ! run ) Fatal("SetParContainers", "No analysis run");
103 }
Int_t run
Definition: autocutx.C:47

Member Data Documentation

PndAnalysis* PndTutAnaTaskD0::fAnalysis
private

Definition at line 52 of file PndTutAnaTaskD0.h.

Referenced by Exec(), and Init().

RhoMassParticleSelector* PndTutAnaTaskD0::fD0Sel
private

Definition at line 59 of file PndTutAnaTaskD0.h.

Referenced by Exec(), and Init().

int PndTutAnaTaskD0::fEvtCount
private

Definition at line 46 of file PndTutAnaTaskD0.h.

Referenced by Exec(), and Init().

TLorentzVector PndTutAnaTaskD0::fIni
private

Definition at line 49 of file PndTutAnaTaskD0.h.

Referenced by Exec(), and PndTutAnaTaskD0().

RhoVtxPoca* PndTutAnaTaskD0::fVtxPoca
private

Definition at line 60 of file PndTutAnaTaskD0.h.

Referenced by Exec(), and Init().

RhoTuple* PndTutAnaTaskD0::nd0
private

Definition at line 58 of file PndTutAnaTaskD0.h.

Referenced by Exec(), Finish(), and Init().


The documentation for this class was generated from the following files: