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