FairRoot/PandaRoot
PndFtfGenerator.cxx
Go to the documentation of this file.
1 // -------------------------------------------------------------------------
2 // ----- PndFtfGenerator source file -----
3 // -------------------------------------------------------------------------
4 
5 
6 #include <iostream>
7 #include "TClonesArray.h"
8 #include "TFile.h"
9 #include "TLorentzVector.h"
10 #include "TTree.h"
11 #include "TVector3.h"
12 #include "TParticle.h"
13 #include "PndFtfGenerator.h"
14 #include "FairPrimaryGenerator.h"
15 
16 
17 
18 
19 using namespace std;
20 
21 // ----- Default constructor ------------------------------------------
23  iEvent = 0;
24  fInputFile = NULL;
25  fInputTree = NULL;
26 }
27 // ------------------------------------------------------------------------
28 
29 
30 
31 // ----- Standard constructor -----------------------------------------
32 PndFtfGenerator::PndFtfGenerator(const Char_t* fileName) {
33  iEvent = 0;
34  fFileName = fileName;
35  fInputFile = new TFile(fFileName);
36  fInputTree = (TTree*) fInputFile->Get("data");
37  fParticles = new TClonesArray("TParticle",100);
38  fInputTree->SetBranchAddress("Particles", &fParticles);
39 }
40 // ------------------------------------------------------------------------
41 
42 
43 
44 // ----- Destructor ---------------------------------------------------
46  CloseInput();
47 }
48 // ------------------------------------------------------------------------
49 
50 
51 
52 // ----- Public method ReadEvent --------------------------------------
53 Bool_t PndFtfGenerator::ReadEvent(FairPrimaryGenerator* primGen) {
54 
55  // Check for input file
56  if ( ! fInputFile ) {
57  cout << "-E PndFtfGenerator: Input file nor open!" << endl;
58  return kFALSE;
59  }
60 
61  // Check for number of events in input file
62  if ( iEvent > fInputTree->GetEntries() ) {
63  cout << "-E PndFtfGenerator: No more events in input file!" << endl;
64  CloseInput();
65  return kFALSE;
66  }
67  TFile *g=gFile;
68  fInputFile->cd();
69  fInputTree->GetEntry(iEvent++);
70  g->cd();
71 
72  // Get number of particle in TClonesrray
73  Int_t nParts = fParticles->GetEntriesFast();
74 
75  // Loop over particles in TClonesArray
76  for (Int_t iPart=0; iPart < nParts; iPart++) {
77  TParticle* part = (TParticle*) fParticles->At(iPart);
78  Int_t pdgType = part->GetPdgCode();
79 
80  // Check if particle type is known to database
81  if ( ! pdgType ) {
82  cout << "-W PndFtfGenerator: Unknown type " << part->GetPdgCode()
83  << ", skipping particle." << endl;
84  continue;
85  }
86 
87  Double_t px = part->Px();
88  Double_t py = part->Py();
89  Double_t pz = part->Pz();
90 
91  Double_t vx = part->Vx();
92  Double_t vy = part->Vy();
93  Double_t vz = part->Vz();
94 
95  // Give track to PrimaryGenerator
96  primGen->AddTrack(pdgType, px, py, pz, vx, vy, vz);
97 
98  } // Loop over particle in event
99 
100  return kTRUE;
101 
102 }
103 // ------------------------------------------------------------------------
104 
105 
106 
107 
108 // ----- Private method CloseInput ------------------------------------
110  if ( fInputFile ) {
111  cout << "-I PndFtfGenerator: Closing input file " << fFileName
112  << endl;
113  fInputFile->Close();
114  delete fInputFile;
115  }
116  fInputFile = NULL;
117 }
118 // ------------------------------------------------------------------------
119 
120 
121 
virtual ~PndFtfGenerator()
TFile * g
FairPrimaryGenerator * primGen
Definition: sim_emc_apd.C:81
Int_t fParticles
Definition: run_full.C:23
Double_t
void CloseInput()
Particle array from PLUTO.
ClassImp(PndAnaContFact)
virtual Bool_t ReadEvent(FairPrimaryGenerator *primGen)
double pz[39]
Definition: pipisigmas.h:14