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

#include <PndKFParticleFinderQA.h>

Inheritance diagram for PndKFParticleFinderQA:

Public Member Functions

 PndKFParticleFinderQA (const char *name="PndKFParticleFinderQA", Int_t iVerbose=0, KFParticleTopoReconstructor *tr=0, TString outFileName="PndKFParticleFinderQA.root")
 
 ~PndKFParticleFinderQA ()
 
void SetEffFileName (const TString &name)
 
void SetMCTrackBranchName (const TString &name)
 
void SetChargedTrackBranchName (const TString &name)
 
void SetNeutralTrackBranchName (const TString &name)
 
virtual InitStatus Init ()
 
virtual void Exec (Option_t *opt)
 
virtual void Finish ()
 
void SetPrintEffFrequency (Int_t n)
 
void SaveParticles (Bool_t b=1)
 
void SaveMCParticles (Bool_t b=1)
 

Private Member Functions

const PndKFParticleFinderQAoperator= (const PndKFParticleFinderQA &)
 
 PndKFParticleFinderQA (const PndKFParticleFinderQA &)
 
void FindClosestMCTrackToBump (const int trackId, int &closestTrack, float &drMin, const float *rReco, const std::vector< std::vector< int > > &mcDaughters)
 
void FindEmcClusterMother (const int iDaughter, int &iMother)
 
void WriteHistosCurFile (TObject *obj)
 
 ClassDef (PndKFParticleFinderQA, 1)
 

Private Attributes

TString fMCTracksBranchName
 
TString fChargedTrackBranchName
 Name of the input TCA with MC tracks. More...
 
TString fNeutralTrackBranchName
 Name of the input TCA with charged tracks. More...
 
TClonesArray * fMCTrackArray
 Name of the input TCA with neutral tracks. More...
 
TClonesArray * fChargedTrackArray
 
TClonesArray * fNeutralTrackArray
 
TClonesArray * fEmcBumps
 
TClonesArray * fRecParticles
 
TClonesArray * fMCParticles
 
TClonesArray * fMatchParticles
 
Bool_t fSaveParticles
 
Bool_t fSaveMCParticles
 
TString fOutFileName
 
TFile * fOutFile
 
TString fEfffileName
 
KFTopoPerformance * fTopoPerformance
 
Int_t fPrintFrequency
 
Int_t fNEvents
 
Double_t fTime [5]
 

Detailed Description

Definition at line 19 of file PndKFParticleFinderQA.h.

Constructor & Destructor Documentation

PndKFParticleFinderQA::PndKFParticleFinderQA ( const char *  name = "PndKFParticleFinderQA",
Int_t  iVerbose = 0,
KFParticleTopoReconstructor *  tr = 0,
TString  outFileName = "PndKFParticleFinderQA.root" 
)

Definition at line 33 of file PndKFParticleFinderQA.cxx.

References Double_t, fOutFile, fOutFileName, fTime, fTopoPerformance, i, KFPartEfficiencies::partPDG, and KFPartEfficiencies::partTitle.

33  :
37  fOutFileName(outFileName), fOutFile(0), fEfffileName("Efficiency.txt"), fTopoPerformance(0), fPrintFrequency(100), fNEvents(0)
38 {
39  for(Int_t i=0; i<5; i++)
40  fTime[i] = 0;
41 
42  fTopoPerformance = new KFTopoPerformance;
43  fTopoPerformance->SetTopoReconstructor(tr);
44 
45  TFile* curFile = gFile;
46  TDirectory* curDirectory = gDirectory;
47 
48  if(!(fOutFileName == ""))
49  fOutFile = new TFile(fOutFileName.Data(),"RECREATE");
50  else
51  fOutFile = gFile;
52  fTopoPerformance->CreateHistos("KFTopoReconstructor",fOutFile);
53 
54  gFile = curFile;
55  gDirectory = curDirectory;
56 
57  //register MC particles from the KFPartEfficiency in the TDatabasePDG
59  int index[4] = {48, 49, 50, 63};
60  {
61  for(int iPall=0; iPall<4; iPall++)
62  {
63  int iParticle = index[iPall];
64  Double_t lifetime = eff.partLifeTime[iParticle]; // lifetime
65  Double_t mass = eff.partMass[iParticle];
66  Int_t PDG = eff.partPDG[iParticle];
67  Double_t charge = eff.partCharge[iParticle];
68  Double_t width = 6.58211928e-22/lifetime;
69 
70  TDatabasePDG::Instance()->AddParticle(eff.partTitle[iParticle], eff.partTitle[iParticle],
71  mass, kFALSE, width, charge, "hadron", PDG);
72  }
73  }
74 }
Int_t i
Definition: run_full.C:25
TString fChargedTrackBranchName
Name of the input TCA with MC tracks.
Double_t
TString partTitle[nParticles]
TString name
KFTopoPerformance * fTopoPerformance
TClonesArray * fNeutralTrackArray
Int_t iVerbose
TString fNeutralTrackBranchName
Name of the input TCA with charged tracks.
TClonesArray * fMCTrackArray
Name of the input TCA with neutral tracks.
TClonesArray * fChargedTrackArray
int partPDG[nParticles]
PndKFParticleFinderQA::~PndKFParticleFinderQA ( )

Definition at line 76 of file PndKFParticleFinderQA.cxx.

References fMatchParticles, fMCParticles, fRecParticles, fSaveMCParticles, fSaveParticles, and fTopoPerformance.

77 {
79 
80  if(fSaveParticles)
81  fRecParticles->Delete();
83  {
84  fMCParticles->Delete();
85  fMatchParticles->Delete();
86  }
87 }
KFTopoPerformance * fTopoPerformance
PndKFParticleFinderQA::PndKFParticleFinderQA ( const PndKFParticleFinderQA )
private

Member Function Documentation

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

Definition at line 144 of file PndKFParticleFinderQA.cxx.

References CAMath::Abs(), Double_t, fChargedTrackArray, fEmcBumps, FindClosestMCTrackToBump(), FindEmcClusterMother(), fMatchParticles, fMCParticles, fMCTrackArray, fNeutralTrackArray, fNEvents, fPrintFrequency, fRecParticles, fSaveMCParticles, fSaveParticles, fTime, fTopoPerformance, PndPidCandidate::GetEmcIndex(), GetEntriesFast(), PndPidCandidate::GetMcIndex(), PndMCTrack::GetMomentum(), PndMCTrack::GetMotherID(), PndMCTrack::GetPdgCode(), PndMCTrack::GetStartVertex(), nEvents, p, R, KFParticleMatch::SetMatch(), KFParticleMatch::SetMatchType(), CAMath::Sqrt(), PndEmcCluster::x(), PndEmcCluster::y(), and PndEmcCluster::z().

145 {
146  if(fSaveParticles)
147  fRecParticles->Delete();
148  if(fSaveMCParticles)
149  {
150  fMCParticles->Delete();
151  fMatchParticles->Delete();
152  }
153 
154  static int nEvents = 0;
155  nEvents++;
156 
157  Int_t nMCTracks = fMCTrackArray->GetEntriesFast();
158  vector<KFMCTrack> mcTracks(nMCTracks);
159  vector< vector<int> > mcDaughters(nMCTracks);
160  for(Int_t iMC=0; iMC<nMCTracks; iMC++)
161  {
162  PndMCTrack *pndMCTrack = (PndMCTrack*)fMCTrackArray->At(iMC);
163 
164  mcTracks[iMC].SetX ( pndMCTrack->GetStartVertex().X() );
165  mcTracks[iMC].SetY ( pndMCTrack->GetStartVertex().Y() );
166  mcTracks[iMC].SetZ ( pndMCTrack->GetStartVertex().Z() );
167  mcTracks[iMC].SetPx( pndMCTrack->GetMomentum().X() );
168  mcTracks[iMC].SetPy( pndMCTrack->GetMomentum().Y() );
169  mcTracks[iMC].SetPz( pndMCTrack->GetMomentum().Z() );
170 
171  Int_t pdg = pndMCTrack->GetPdgCode();
172  Double_t q=1;
173  if ( pdg < 20000 )
174  q = TDatabasePDG::Instance()->GetParticle(pdg)->Charge()/3.0;
175  else q = 0;
176  Double_t p = pndMCTrack->GetMomentum().Mag();
177 
178  mcTracks[iMC].SetMotherId(pndMCTrack->GetMotherID());
179  mcTracks[iMC].SetQP(q/p);
180  mcTracks[iMC].SetPDG(pdg);
181  mcTracks[iMC].SetNMCPoints(0);
182 
183  if( (mcTracks[iMC].R() > 10) && (pdg == 22) )
184  mcTracks[iMC].SetOutOfDetector();
185 
186  if(pndMCTrack->GetMotherID() > -1)
187  mcDaughters[pndMCTrack->GetMotherID()].push_back(iMC);
188  }
189 
190 // for(Int_t iMC=0; iMC<nMCTracks; iMC++)
191 // {
192 // PndMCTrack *pndMCTrack = (PndMCTrack*)fMCTrackArray->At(iMC);
193 // Int_t pdg = pndMCTrack->GetPdgCode();
194 //
195 // if( (mcTracks[iMC].R() < 10) && (pdg == 22) )
196 // {
197 // PndMCTrack *pndMCTrack1 = (PndMCTrack*)fMCTrackArray->At(mcDaughters[iMC][0]);
198 // std::cout << "imc " << iMC << " r " << pndMCTrack1->GetStartVertex().X() << " " << pndMCTrack1->GetStartVertex().Y() << " " << pndMCTrack1->GetStartVertex().Z() << " pdg " << pdg <<std::endl;
199 // }
200 // }
201 
202  Int_t ntrackMatches=fChargedTrackArray->GetEntriesFast();
204  ntrackMatches += fNeutralTrackArray->GetEntriesFast();
205 
206  vector<int> trackMatch(ntrackMatches, -1);
207 static int Nghosts = 0;
208  for(int iTr=0; iTr<ntrackMatches; iTr++)
209  {
210  PndPidCandidate *inTrack;
211  if(iTr<fChargedTrackArray->GetEntriesFast())
212  inTrack = (PndPidCandidate*)fChargedTrackArray->At(iTr);
213  else
214  inTrack = (PndPidCandidate*)fNeutralTrackArray->At(iTr - fChargedTrackArray->GetEntriesFast());
215  if(!inTrack) continue;
216  Int_t trId = inTrack->GetMcIndex();
217 
218  if(trId >= 0)
219  {
220  PndMCTrack *pndMCTrack = (PndMCTrack*)fMCTrackArray->At(trId);
221  if(iTr<fChargedTrackArray->GetEntriesFast() && TMath::Abs(pndMCTrack->GetPdgCode()) > 2500)
222  continue;
223  //Correct linking between EMC clusters and MC tracks
224  float rMC[3] = { static_cast<float>(pndMCTrack->GetStartVertex().X()),
225  static_cast<float>(pndMCTrack->GetStartVertex().Y()),
226  static_cast<float>(pndMCTrack->GetStartVertex().Z()) };
227  float R = TMath::Sqrt(rMC[0]*rMC[0] + rMC[1]*rMC[1] + rMC[2]*rMC[2]);
228 
229  //if cluster points to some mother or grandmother of the gamma
230  if(TMath::Abs(pndMCTrack->GetPdgCode()) != 22 && iTr>=fChargedTrackArray->GetEntriesFast())
231  {
232  float drMin = 1.e9f;
233  int closestTrack = -1;
234  int iMother = -1;
235 
236  int iBump = inTrack->GetEmcIndex();
237  PndEmcBump *inBump = (PndEmcBump*) fEmcBumps->At(iBump);
238  float rReco[3] = { static_cast<float>(inBump->x()),
239  static_cast<float>(inBump->y()),
240  static_cast<float>(inBump->z()) };
241  FindClosestMCTrackToBump(trId, closestTrack, drMin, rReco, mcDaughters);
242  FindEmcClusterMother(closestTrack, iMother);
243 
244  if(iMother < 0)
245  trId = -1;
246  else
247  {
248  PndMCTrack* mother = (PndMCTrack*)fMCTrackArray->At( iMother );
249  if(mother->GetPdgCode() == 22)
250  trId = iMother;
251  else
252  trId = -1;
253  }
254  }
255  }
256 
257  if(trId < 0) continue;
258  mcTracks[trId].SetReconstructed();
259  trackMatch[iTr] = trId;
260  }
261 
262  fTopoPerformance->SetMCTracks(mcTracks);
263  fTopoPerformance->SetTrackMatch(trackMatch);
264 
265  fTopoPerformance->CheckMCTracks();
266  fTopoPerformance->MatchTracks();
267  fTopoPerformance->FillHistos();
268 
269  fNEvents++;
270  fTime[4] += fTopoPerformance->GetTopoReconstructor()->Time();
271  for(int iT=0; iT<4; iT++)
272  fTime[iT] += fTopoPerformance->GetTopoReconstructor()->StatTime(iT);
273  if(fNEvents%fPrintFrequency == 0)
274  {
275  std::cout << "Topo reconstruction time"
276  << " Real = " << std::setw( 10 ) << fTime[4]/fNEvents * 1.e3 << " ms" << std::endl;
277  std::cout << " Init " << fTime[0]/fNEvents * 1.e3 << " ms" << std::endl;
278  std::cout << " PV Finder " << fTime[1]/fNEvents * 1.e3 << " ms" << std::endl;
279  std::cout << " Sort Tracks " << fTime[2]/fNEvents * 1.e3 << " ms" << std::endl;
280  std::cout << " KF Particle Finder " << fTime[3]/fNEvents * 1.e3 << " ms" << std::endl;
281  }
282 
283  // save particles to a ROOT file
284  if(fSaveParticles)
285  {
286  for(unsigned int iP=0; iP < fTopoPerformance->GetTopoReconstructor()->GetParticles().size(); iP++)
287  {
288  new((*fRecParticles)[iP]) KFParticle(fTopoPerformance->GetTopoReconstructor()->GetParticles()[iP]);
289  }
290  }
291 
292  if(fSaveMCParticles)
293  {
294  for(unsigned int iP=0; iP < fTopoPerformance->GetTopoReconstructor()->GetParticles().size(); iP++)
295  {
296  new((*fMatchParticles)[iP]) KFParticleMatch();
298 
299  Short_t matchType = 0;
300  int iMCPart = -1;
301  if(!(fTopoPerformance->ParticlesMatch()[iP].IsMatchedWithPdg())) //background
302  {
303  if(fTopoPerformance->ParticlesMatch()[iP].IsMatched())
304  {
305  iMCPart = fTopoPerformance->ParticlesMatch()[iP].GetBestMatchWithPdg();
306  matchType = 1;
307  }
308  }
309  else
310  {
311  iMCPart = fTopoPerformance->ParticlesMatch()[iP].GetBestMatchWithPdg();
312  matchType = 2;
313  }
314 
315  p->SetMatch(iMCPart);
316  p->SetMatchType(matchType);
317  }
318 
319  for(unsigned int iP=0; iP < fTopoPerformance->MCParticles().size(); iP++)
320  {
321  new((*fMCParticles)[iP]) KFMCParticle(fTopoPerformance->MCParticles()[iP]);
322  }
323  }
324 }
void SetMatchType(Short_t i)
static T Sqrt(const T &x)
Definition: PndCAMath.h:37
Int_t GetPdgCode() const
Definition: PndMCTrack.h:73
TVector3 GetMomentum() const
Definition: PndMCTrack.h:78
void SetMatch(Int_t i)
cout<< "POINTs for new FwEndCap == "<< tsim-> GetEntriesFast()
Double_t p
Definition: anasim.C:58
Int_t GetMcIndex() const
Int_t GetEmcIndex() const
static T Abs(const T &x)
Definition: PndCAMath.h:39
Double_t z() const
Double_t
Int_t nEvents
Definition: hit_dirc.C:11
Double_t x() const
KFTopoPerformance * fTopoPerformance
TClonesArray * fNeutralTrackArray
void FindClosestMCTrackToBump(const int trackId, int &closestTrack, float &drMin, const float *rReco, const std::vector< std::vector< int > > &mcDaughters)
TVector3 GetStartVertex() const
Definition: PndMCTrack.h:76
Int_t GetMotherID() const
Definition: PndMCTrack.h:74
TClonesArray * fMCTrackArray
Name of the input TCA with neutral tracks.
represents a reconstructed (splitted) emc cluster
Definition: PndEmcBump.h:34
Double_t R
Definition: checkhelixhit.C:61
Double_t y() const
TClonesArray * fChargedTrackArray
void FindEmcClusterMother(const int iDaughter, int &iMother)
void PndKFParticleFinderQA::FindClosestMCTrackToBump ( const int  trackId,
int &  closestTrack,
float &  drMin,
const float *  rReco,
const std::vector< std::vector< int > > &  mcDaughters 
)
private

Definition at line 349 of file PndKFParticleFinderQA.cxx.

References dx, dy, dz, fMCTrackArray, PndMCTrack::GetStartVertex(), and sqrt().

Referenced by Exec().

351 {
352  if( (trackId<0) || (trackId>=mcDaughters.size()) ) return;
353 
354  for(int iDaughter=0; iDaughter<mcDaughters[trackId].size(); iDaughter++)
355  {
356  int daughterIndex = mcDaughters[trackId][iDaughter];
357  PndMCTrack* daughter = (PndMCTrack*)fMCTrackArray->At( daughterIndex );
358  float rMC[3] = { static_cast<float>(daughter->GetStartVertex().X()),
359  static_cast<float>(daughter->GetStartVertex().Y()),
360  static_cast<float>(daughter->GetStartVertex().Z()) };
361 
362  float dx = rReco[0] - rMC[0];
363  float dy = rReco[1] - rMC[1];
364  float dz = rReco[2] - rMC[2];
365  float dr = sqrt(dx*dx + dy*dy + dz*dz);
366  if(dr < drMin)
367  {
368  drMin = dr;
369  closestTrack = daughterIndex;
370  }
371 
372  if(mcDaughters[daughterIndex].size() > 0)
373  FindClosestMCTrackToBump( daughterIndex, closestTrack, drMin, rReco, mcDaughters);
374  }
375 }
double dy
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
double dx
void FindClosestMCTrackToBump(const int trackId, int &closestTrack, float &drMin, const float *rReco, const std::vector< std::vector< int > > &mcDaughters)
TVector3 GetStartVertex() const
Definition: PndMCTrack.h:76
TClonesArray * fMCTrackArray
Name of the input TCA with neutral tracks.
void PndKFParticleFinderQA::FindEmcClusterMother ( const int  iDaughter,
int &  iMother 
)
private

Definition at line 377 of file PndKFParticleFinderQA.cxx.

References fMCTrackArray, PndMCTrack::GetMotherID(), PndMCTrack::GetStartVertex(), R, and CAMath::Sqrt().

Referenced by Exec().

378 {
379  PndMCTrack* daughter = (PndMCTrack*)fMCTrackArray->At( iDaughter );
380  float rMC[3] = { static_cast<float>(daughter->GetStartVertex().X()),
381  static_cast<float>(daughter->GetStartVertex().Y()),
382  static_cast<float>(daughter->GetStartVertex().Z()) };
383  float R = TMath::Sqrt(rMC[0]*rMC[0] + rMC[1]*rMC[1] + rMC[2]*rMC[2]);
384 
385 // while( !(R<80. && daughter->GetPdgCode() == 22) )
386  while( R > 10. )
387  {
388  iMother = daughter->GetMotherID();
389  if(iMother < 0) break;
390 
391  daughter = (PndMCTrack*)fMCTrackArray->At( iMother );
392  rMC[0] = daughter->GetStartVertex().X();
393  rMC[1] = daughter->GetStartVertex().Y();
394  rMC[2] = daughter->GetStartVertex().Z();
395 
396  R = TMath::Sqrt(rMC[0]*rMC[0] + rMC[1]*rMC[1] + rMC[2]*rMC[2]);
397  }
398 }
static T Sqrt(const T &x)
Definition: PndCAMath.h:37
TVector3 GetStartVertex() const
Definition: PndMCTrack.h:76
Int_t GetMotherID() const
Definition: PndMCTrack.h:74
TClonesArray * fMCTrackArray
Name of the input TCA with neutral tracks.
Double_t R
Definition: checkhelixhit.C:61
void PndKFParticleFinderQA::Finish ( )
virtual

Definition at line 326 of file PndKFParticleFinderQA.cxx.

References fEfffileName, fOutFile, fOutFileName, fTopoPerformance, out, and WriteHistosCurFile().

327 {
328  TDirectory *curr = gDirectory;
329  TFile *currentFile = gFile;
330  // Open output file and write histograms
331 
332  fOutFile->cd();
333  WriteHistosCurFile(fTopoPerformance->GetHistosDirectory());
334 // if (histodirmod!=NULL) WriteHistos(histodirmod);
335  //WriteHistos(gDirectory);
336  if(!(fOutFileName == ""))
337  {
338  fOutFile->Close();
339  fOutFile->Delete();
340  }
341  gFile = currentFile;
342  gDirectory = curr;
343 
344  std::fstream eff(fEfffileName.Data(),fstream::out);
345  eff << fTopoPerformance->fParteff;
346  eff.close();
347 }
TFile * out
Definition: reco_muo.C:20
KFTopoPerformance * fTopoPerformance
void WriteHistosCurFile(TObject *obj)
InitStatus PndKFParticleFinderQA::Init ( )
virtual

Definition at line 89 of file PndKFParticleFinderQA.cxx.

References fChargedTrackArray, fChargedTrackBranchName, fEmcBumps, fMatchParticles, fMCParticles, fMCTrackArray, fMCTracksBranchName, fNeutralTrackArray, fNeutralTrackBranchName, fRecParticles, fSaveMCParticles, and fSaveParticles.

90 {
91  //Get ROOT Manager
92  FairRootManager* ioman= FairRootManager::Instance();
93 
94  if(ioman==0)
95  {
96  Error("PndKFParticleFinderQA::Init","RootManager not instantiated!");
97  return kERROR;
98  }
99 
100  //MC Tracks
101  fMCTrackArray=(TClonesArray*) ioman->GetObject(fMCTracksBranchName);
102  if(fMCTrackArray==0)
103  {
104  Error("StandaloneHitGenerator::Init","mc track array not found!");
105  return kERROR;
106  }
107 
108  //Charged Tracks
109  fChargedTrackArray=(TClonesArray*) ioman->GetObject(fChargedTrackBranchName);
110  if(fChargedTrackArray==0)
111  {
112  Error("StandaloneHitGenerator::Init","Charged tracks array not found!");
113  return kERROR;
114  }
115 
116  //Neutral Tracks
117  fNeutralTrackArray=(TClonesArray*) ioman->GetObject(fNeutralTrackBranchName);
118  if(fNeutralTrackArray==0)
119  {
120  Error("StandaloneHitGenerator::Init","Neutral tracks match array not found! Continue without gammas.");
121  }
122  fEmcBumps = (TClonesArray*) ioman->GetObject("EmcBump");
123 
124  if(fSaveParticles)
125  {
126  // create and register TClonesArray with output reco particles
127  fRecParticles = new TClonesArray("KFParticle",100);
128  ioman->Register("RecoParticles", "KFParticle", fRecParticles, kTRUE);
129  }
130 
131  if(fSaveMCParticles)
132  {
133  // create and register TClonesArray with output MC particles
134  fMCParticles = new TClonesArray("KFMCParticle",100);
135  ioman->Register("KFMCParticles", "KFParticle", fMCParticles, kTRUE);
136 
137  // create and register TClonesArray with matching between reco and MC particles
138  fMatchParticles = new TClonesArray("KFParticleMatch",100);
139  ioman->Register("KFParticleMatch", "KFParticle", fMatchParticles, kTRUE);
140  }
141  return kSUCCESS;
142 }
TString fChargedTrackBranchName
Name of the input TCA with MC tracks.
TClonesArray * fNeutralTrackArray
TString fNeutralTrackBranchName
Name of the input TCA with charged tracks.
TClonesArray * fMCTrackArray
Name of the input TCA with neutral tracks.
TClonesArray * fChargedTrackArray
const PndKFParticleFinderQA& PndKFParticleFinderQA::operator= ( const PndKFParticleFinderQA )
private
void PndKFParticleFinderQA::SaveMCParticles ( Bool_t  b = 1)
inline

Definition at line 39 of file PndKFParticleFinderQA.h.

References b, and fSaveMCParticles.

39 { fSaveMCParticles = b; }
TTree * b
void PndKFParticleFinderQA::SaveParticles ( Bool_t  b = 1)
inline

Definition at line 38 of file PndKFParticleFinderQA.h.

References b, and fSaveParticles.

38 { fSaveParticles = b; }
TTree * b
void PndKFParticleFinderQA::SetChargedTrackBranchName ( const TString name)
inline

Definition at line 29 of file PndKFParticleFinderQA.h.

References fChargedTrackBranchName, and name.

Referenced by kfparticle().

TString fChargedTrackBranchName
Name of the input TCA with MC tracks.
TString name
void PndKFParticleFinderQA::SetEffFileName ( const TString name)
inline

Definition at line 27 of file PndKFParticleFinderQA.h.

References fEfffileName, and name.

27 { fEfffileName = name; }
TString name
void PndKFParticleFinderQA::SetMCTrackBranchName ( const TString name)
inline

Definition at line 28 of file PndKFParticleFinderQA.h.

References fMCTracksBranchName, and name.

Referenced by kfparticle().

TString name
void PndKFParticleFinderQA::SetNeutralTrackBranchName ( const TString name)
inline

Definition at line 30 of file PndKFParticleFinderQA.h.

References fNeutralTrackBranchName, and name.

Referenced by kfparticle().

TString name
TString fNeutralTrackBranchName
Name of the input TCA with charged tracks.
void PndKFParticleFinderQA::SetPrintEffFrequency ( Int_t  n)

Definition at line 417 of file PndKFParticleFinderQA.cxx.

References fPrintFrequency, fTopoPerformance, and n.

Referenced by kfparticle().

418 {
419  fTopoPerformance->SetPrintEffFrequency(n);
420  fPrintFrequency = n;
421 }
int n
KFTopoPerformance * fTopoPerformance
void PndKFParticleFinderQA::WriteHistosCurFile ( TObject *  obj)
private

Definition at line 400 of file PndKFParticleFinderQA.cxx.

References obj.

Referenced by Finish().

400  {
401  if( !obj->IsFolder() ) obj->Write();
402  else{
403  TDirectory *cur = gDirectory;
404  TFile *currentFile = gFile;
405 
406  TDirectory *sub = cur->GetDirectory(obj->GetName());
407  sub->cd();
408  TList *listSub = (static_cast<TDirectory*>(obj))->GetList();
409  TIter it(listSub);
410  while( TObject *obj1=it() ) WriteHistosCurFile(obj1);
411  cur->cd();
412  gFile = currentFile;
413  gDirectory = cur;
414  }
415 }
void WriteHistosCurFile(TObject *obj)
PndAnaPidSelector *& obj

Member Data Documentation

TClonesArray* PndKFParticleFinderQA::fChargedTrackArray
private

Definition at line 59 of file PndKFParticleFinderQA.h.

Referenced by Exec(), and Init().

TString PndKFParticleFinderQA::fChargedTrackBranchName
private

Name of the input TCA with MC tracks.

Definition at line 54 of file PndKFParticleFinderQA.h.

Referenced by Init(), and SetChargedTrackBranchName().

TString PndKFParticleFinderQA::fEfffileName
private

Definition at line 73 of file PndKFParticleFinderQA.h.

Referenced by Finish(), and SetEffFileName().

TClonesArray* PndKFParticleFinderQA::fEmcBumps
private

Definition at line 61 of file PndKFParticleFinderQA.h.

Referenced by Exec(), and Init().

TClonesArray* PndKFParticleFinderQA::fMatchParticles
private

Definition at line 65 of file PndKFParticleFinderQA.h.

Referenced by Exec(), Init(), and ~PndKFParticleFinderQA().

TClonesArray* PndKFParticleFinderQA::fMCParticles
private

Definition at line 64 of file PndKFParticleFinderQA.h.

Referenced by Exec(), Init(), and ~PndKFParticleFinderQA().

TClonesArray* PndKFParticleFinderQA::fMCTrackArray
private

Name of the input TCA with neutral tracks.

Definition at line 58 of file PndKFParticleFinderQA.h.

Referenced by Exec(), FindClosestMCTrackToBump(), FindEmcClusterMother(), and Init().

TString PndKFParticleFinderQA::fMCTracksBranchName
private

Definition at line 53 of file PndKFParticleFinderQA.h.

Referenced by Init(), and SetMCTrackBranchName().

TClonesArray* PndKFParticleFinderQA::fNeutralTrackArray
private

Definition at line 60 of file PndKFParticleFinderQA.h.

Referenced by Exec(), and Init().

TString PndKFParticleFinderQA::fNeutralTrackBranchName
private

Name of the input TCA with charged tracks.

Definition at line 55 of file PndKFParticleFinderQA.h.

Referenced by Init(), and SetNeutralTrackBranchName().

Int_t PndKFParticleFinderQA::fNEvents
private

Definition at line 78 of file PndKFParticleFinderQA.h.

Referenced by Exec().

TFile* PndKFParticleFinderQA::fOutFile
private

Definition at line 72 of file PndKFParticleFinderQA.h.

Referenced by Finish(), and PndKFParticleFinderQA().

TString PndKFParticleFinderQA::fOutFileName
private

Definition at line 71 of file PndKFParticleFinderQA.h.

Referenced by Finish(), and PndKFParticleFinderQA().

Int_t PndKFParticleFinderQA::fPrintFrequency
private

Definition at line 77 of file PndKFParticleFinderQA.h.

Referenced by Exec(), and SetPrintEffFrequency().

TClonesArray* PndKFParticleFinderQA::fRecParticles
private

Definition at line 63 of file PndKFParticleFinderQA.h.

Referenced by Exec(), Init(), and ~PndKFParticleFinderQA().

Bool_t PndKFParticleFinderQA::fSaveMCParticles
private

Definition at line 68 of file PndKFParticleFinderQA.h.

Referenced by Exec(), Init(), SaveMCParticles(), and ~PndKFParticleFinderQA().

Bool_t PndKFParticleFinderQA::fSaveParticles
private

Definition at line 67 of file PndKFParticleFinderQA.h.

Referenced by Exec(), Init(), SaveParticles(), and ~PndKFParticleFinderQA().

Double_t PndKFParticleFinderQA::fTime[5]
private

Definition at line 79 of file PndKFParticleFinderQA.h.

Referenced by Exec(), and PndKFParticleFinderQA().

KFTopoPerformance* PndKFParticleFinderQA::fTopoPerformance
private

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