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

#include <PndAnaWithTrigger.h>

Inheritance diagram for PndAnaWithTrigger:

Public Member Functions

 PndAnaWithTrigger (double pbarmom, TString outname)
 
 ~PndAnaWithTrigger ()
 
virtual InitStatus Init ()
 
virtual void Exec (Option_t *opt)
 
virtual void Finish ()
 

Private Member Functions

 ClassDef (PndAnaWithTrigger, 1)
 

Private Attributes

int fEvtCount
 
TLorentzVector fIni
 
PndAnalysisfAnalysis
 
TDatabasePDG * fPdg
 
TFile * fFile
 
TString fOutName
 
RhoTuplentp1
 
RhoTuplentp2
 
RhoTuplenmc
 
RhoMassParticleSelectorjpsiMassSel
 
RhoMassParticleSelectorjpsiPreMassSel
 
TClonesArray * fOnlineFilterInfo
 

Detailed Description

Definition at line 24 of file PndAnaWithTrigger.h.

Constructor & Destructor Documentation

PndAnaWithTrigger::PndAnaWithTrigger ( double  pbarmom,
TString  outname 
)

Definition at line 63 of file PndAnaWithTrigger.cxx.

References fIni, fOutName, mp, and sqrt().

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

Definition at line 75 of file PndAnaWithTrigger.cxx.

References fFile.

76 {
77  delete fFile;
78 }

Member Function Documentation

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

Definition at line 132 of file PndAnaWithTrigger.cxx.

References RhoKinFitter::Add4MomConstraint(), all, RhoTuple::Column(), RhoTuple::DumpData(), fAnalysis, fEvtCount, PndAnalysis::FillList(), fIni, RhoKinFitter::Fit(), fOnlineFilterInfo, RhoFitterBase::GetChi2(), PndAnalysis::GetEventInTask(), RhoCandidate::GetFit(), RhoCandidate::GetMcTruth(), PndOnlineFilterInfo::GetNTag(), PndOnlineFilterInfo::GetNTagTotal(), i, jpsiMassSel, jpsiPreMassSel, mclist, PndAnalysis::McTruthMatch(), muminus, muplus, nmc, ntp1, ntp2, RhoCandidate::P4(), piminus, piplus, PndOnlineFilterInfo::Tagged(), and TString.

133 {
134  // *** some variables
135  int i=0,j=0;
136 
137  // *** necessary to read the next event
139 
140  // *** print event counter
141  if (!(++fEvtCount%100)) cout << "evt "<<fEvtCount<<endl;
142 
143  // *******
144  // ******* PUT ANALYSIS CODE HERE
145  // *******
146 
147  // *** fetch trigger info
148  PndOnlineFilterInfo *stInfo = 0;
149  if (fOnlineFilterInfo) stInfo = ( PndOnlineFilterInfo* ) fOnlineFilterInfo->At ( 0 );
150 
151  TString pidalg = "PidChargedProbability";
152 
153  // *** QA tool for simple dumping of analysis results in RhoRuple
154  // *** WIKI: https://panda-wiki.gsi.de/foswiki/bin/view/Computing/PandaRootAnalysisJuly13#PndRhoTupleQA
155  PndRhoTupleQA qa(fAnalysis,fIni.P());
156 
157  // *** RhoCandLists for the analysis
158  RhoCandList muplus, muminus, piplus, piminus, jpsi, psi2s, all, mclist;
159 
160 
161  // *** Select with PID info pidalg and ('All'); type and mass are set
162  fAnalysis->FillList(muplus, "MuonAllPlus", pidalg);
163  fAnalysis->FillList(muminus, "MuonAllMinus", pidalg);
164  fAnalysis->FillList(piplus, "PionAllPlus", pidalg);
165  fAnalysis->FillList(piminus, "PionAllMinus", pidalg);
166 
167 
168  // *** Setup event shape object
169  fAnalysis->FillList(all, "All", pidalg);
170  PndEventShape evsh(all, fIni, 0.05, 0.1);
171 
172 
173  // *** store MC info in ntuple
174  fAnalysis->FillList(mclist, "McTruth");
175 
176  nmc->Column("ev", (Int_t) fEvtCount);
177  qa.qaMcList("", mclist, nmc);
178  nmc->DumpData();
179 
180 
181  // *** combinatorics for J/psi -> mu+ mu-
182  jpsi.Combine(muplus, muminus);
183  jpsi.Select(jpsiPreMassSel);
184  jpsi.SetType(443);
185  int njmct = fAnalysis->McTruthMatch(jpsi); // match the whole list to count #matches (should be only 1)
186 
187  // *** write ntuple for the jpsi reconstruction
188  for (j=0;j<jpsi.GetLength();++j)
189  {
190  // some general info about event (actually written for each candidate!)
191  ntp1->Column("ev", (Float_t) i);
192  ntp1->Column("cand", (Float_t) j);
193  ntp1->Column("ncand", (Float_t) jpsi.GetLength());
194  ntp1->Column("nmct", (Float_t) njmct);
195 
196  // *** store info from trigger
197  if (stInfo)
198  {
199  ntp1->Column("sttrig", (Int_t) stInfo->Tagged()); // event triggered
200  ntp1->Column("stntot", (Int_t) stInfo->GetNTagTotal()); // total number of triggered candidates from all active lines
201  ntp1->Column("stn200", (Int_t) stInfo->GetNTag(200)); // number of triggered candidates from J/psi->e+ e- line
202  ntp1->Column("stn201", (Int_t) stInfo->GetNTag(201)); // number of triggered candidates from J/psi->mu+ mu- line
203  }
204 
205  // store info about initial 4-vector
206  qa.qaP4("beam", fIni, ntp1);
207 
208  // dump information about composite candidate tree recursively (see analysis/AnalysisTools/PndRhoTupleQA)
209  qa.qaComp("j", jpsi[j], ntp1);
210 
211  // dump info about event shapes
212  qa.qaEventShapeShort("es",&evsh, ntp1);
213 
214  // *** store the 4-vector of the truth matched candidate (or a dummy, if not matched to keep ntuple consistent)
215  RhoCandidate *truth = jpsi[j]->GetMcTruth();
216  TLorentzVector lv;
217  if (truth) lv = truth->P4();
218  qa.qaP4("trj", lv, ntp1);
219 
220  ntp1->DumpData();
221  }
222 
223 
224  // *** some rough mass selection on J/psi
225  jpsi.Select(jpsiMassSel);
226 
227  // *** combinatorics for psi(2S) -> J/psi pi+ pi-
228  psi2s.Combine(jpsi, piplus, piminus);
229  psi2s.SetType(88880);
230  int npsimct = fAnalysis->McTruthMatch(psi2s); // match the whole list to count #matches (should be only 1)
231 
232  // *** write ntuple for the psi(2S) reconstruction
233  for (j=0;j<psi2s.GetLength();++j)
234  {
235  // some general info about event (actually written for each candidate!)
236  ntp2->Column("ev", (Float_t) i);
237  ntp2->Column("cand", (Float_t) j);
238  ntp2->Column("ncand", (Float_t) psi2s.GetLength());
239  ntp2->Column("nmct", (Float_t) npsimct);
240 
241  // *** store info from trigger
242  if (stInfo)
243  {
244  ntp2->Column("sttrig", (Int_t) stInfo->Tagged()); // event triggered
245  ntp2->Column("stntot", (Int_t) stInfo->GetNTagTotal()); // total number of triggered candidates from all active lines
246  ntp2->Column("stn200", (Int_t) stInfo->GetNTag(200)); // number of triggered candidates from J/psi->e+ e- line
247  ntp2->Column("stn201", (Int_t) stInfo->GetNTag(201)); // number of triggered candidates from J/psi->mu+ mu- line
248  }
249 
250  RhoKinFitter kinfit(psi2s[j]);
251  kinfit.Add4MomConstraint(fIni);
252  kinfit.Fit();
253 
254  RhoCandidate *psifit=psi2s[j]->GetFit();
255  // store info about initial 4-vector
256  qa.qaP4("beam", fIni, ntp2);
257 
258  // dump information about composite candidate tree recursively (see analysis/AnalysisTools/PndRhoTupleQA)
259  qa.qaComp("psi", psi2s[j], ntp2);
260  qa.qaComp("fpsi",psifit, ntp2);
261  ntp2->Column("fchi2", (Float_t) kinfit.GetChi2());
262 
263  // dump info about event shapes
264  qa.qaEventShapeShort("es",&evsh, ntp2);
265 
266  // *** store the 4-vector of the truth matched candidate (or a dummy, if not matched to keep ntuple consistent)
267  RhoCandidate *truth = psi2s[j]->GetMcTruth();
268  TLorentzVector lv;
269  if (truth) lv = truth->P4();
270  qa.qaP4("trpsi", lv, ntp2);
271 
272  ntp2->DumpData();
273  }
274 
275 }
CandList mclist
RhoMassParticleSelector * jpsiPreMassSel
Int_t i
Definition: run_full.C:25
Bool_t FillList(RhoCandList &l, TString listkey="All", TString pidTcaNames="", int trackHypothesis=-1)
CandList piplus
RhoMassParticleSelector * jpsiMassSel
CandList muplus
void GetEventInTask()
TLorentzVector fIni
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
RhoCandidate * GetMcTruth() const
Definition: RhoCandidate.h:437
void DumpData()
Definition: RhoTuple.cxx:391
CandList piminus
RhoCandidate * GetFit() const
Definition: RhoCandidate.h:293
Bool_t McTruthMatch(RhoCandidate *cand, Int_t level=2, bool verbose=false)
CandList muminus
TClonesArray * fOnlineFilterInfo
PndAnalysis * fAnalysis
void PndAnaWithTrigger::Finish ( )
virtual

Definition at line 278 of file PndAnaWithTrigger.cxx.

References fFile.

279 {
280 
281  // *******
282  // ******* STORE YOUR HISTOS AND TUPLES
283  // *******
284 
285  fFile->Write();
286  fFile->Close();
287 
288 }
InitStatus PndAnaWithTrigger::Init ( )
virtual

Definition at line 83 of file PndAnaWithTrigger.cxx.

References fAnalysis, fEvtCount, fFile, fOnlineFilterInfo, fOutName, fPdg, RhoTuple::GetInternalTree(), jpsiMassSel, jpsiPreMassSel, PndAnalysis::McMatchAllowPhotos(), nmc, ntp1, and ntp2.

84 {
85  // *** initialize PndAnalysis object
86  fAnalysis = new PndAnalysis();
88 
89  // *** reset the event counter
90  fEvtCount = 0;
91 
92  // *******
93  // ******* PREPARE/CREATE THE STUFF YOU NEED
94  // *******
95 
96  fPdg = TDatabasePDG::Instance();
97 
98  // *** Mass selector for the jpsi cands
99  double m0_jpsi = fPdg->GetParticle("J/psi")->Mass(); // Get nominal PDG mass of the J/psi
100  jpsiMassSel = new RhoMassParticleSelector("jpsi" ,m0_jpsi,0.5);
101  jpsiPreMassSel = new RhoMassParticleSelector("jpsipre",m0_jpsi,2.0);
102 
103  // ***
104  // *** Prepare RhoTuple output
105  // ***
106  TDirectory *dir = gDirectory;
107  fFile=new TFile(fOutName,"RECREATE");
108  fFile->cd();
109 
110  // *** create some ntuples
111  ntp1 = new RhoTuple("ntp1", "jpsi analysis");
112  ntp2 = new RhoTuple("ntp2", "psi(2S) analysis");
113  nmc = new RhoTuple("nmc", "mctruth info");
114 
115  // assign RhoTuples to outfile
116  if (ntp1) ntp1->GetInternalTree()->SetDirectory(gDirectory);
117  if (ntp2) ntp2->GetInternalTree()->SetDirectory(gDirectory);
118  if (nmc) nmc->GetInternalTree()->SetDirectory(gDirectory);
119 
120  // *** restore original gDirectory
121  dir->cd();
122 
123  // *** Connect to the Online Filter Info
124  fOnlineFilterInfo = ( TClonesArray* ) FairRootManager::Instance()->GetObject ( "OnlineFilterInfo" );
125 
126  return kSUCCESS;
127 }
TDatabasePDG * fPdg
RhoMassParticleSelector * jpsiPreMassSel
RhoMassParticleSelector * jpsiMassSel
void McMatchAllowPhotos(int maxn=1, double thresh=0.05)
Definition: PndAnalysis.h:62
TTree * GetInternalTree()
Definition: RhoTuple.h:207
TClonesArray * fOnlineFilterInfo
PndAnalysis * fAnalysis

Member Data Documentation

PndAnalysis* PndAnaWithTrigger::fAnalysis
private

Definition at line 54 of file PndAnaWithTrigger.h.

Referenced by Exec(), and Init().

int PndAnaWithTrigger::fEvtCount
private

Definition at line 48 of file PndAnaWithTrigger.h.

Referenced by Exec(), and Init().

TFile* PndAnaWithTrigger::fFile
private

Definition at line 56 of file PndAnaWithTrigger.h.

Referenced by Finish(), Init(), and ~PndAnaWithTrigger().

TLorentzVector PndAnaWithTrigger::fIni
private

Definition at line 51 of file PndAnaWithTrigger.h.

Referenced by Exec(), and PndAnaWithTrigger().

TClonesArray* PndAnaWithTrigger::fOnlineFilterInfo
private

Definition at line 72 of file PndAnaWithTrigger.h.

Referenced by Exec(), and Init().

TString PndAnaWithTrigger::fOutName
private

Definition at line 57 of file PndAnaWithTrigger.h.

Referenced by Init(), and PndAnaWithTrigger().

TDatabasePDG* PndAnaWithTrigger::fPdg
private

Definition at line 55 of file PndAnaWithTrigger.h.

Referenced by Init().

RhoMassParticleSelector* PndAnaWithTrigger::jpsiMassSel
private

Definition at line 68 of file PndAnaWithTrigger.h.

Referenced by Exec(), and Init().

RhoMassParticleSelector* PndAnaWithTrigger::jpsiPreMassSel
private

Definition at line 69 of file PndAnaWithTrigger.h.

Referenced by Exec(), and Init().

RhoTuple* PndAnaWithTrigger::nmc
private

Definition at line 66 of file PndAnaWithTrigger.h.

Referenced by Exec(), and Init().

RhoTuple* PndAnaWithTrigger::ntp1
private

Definition at line 64 of file PndAnaWithTrigger.h.

Referenced by Exec(), and Init().

RhoTuple* PndAnaWithTrigger::ntp2
private

Definition at line 65 of file PndAnaWithTrigger.h.

Referenced by Exec(), and Init().


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