FairRoot/PandaRoot
pgenerators/FileReaders/PndHypBupGenerator.cxx
Go to the documentation of this file.
1 // -------------------------------------------------------------------------
2 // ----- PndHypBupGenerator source file -----
3 // ----- Created by A. Sanchez -----
4 // -------------------------------------------------------------------------
5 
6 
7 #include <iostream>
8 #include "TClonesArray.h"
9 
10 #include "TFile.h"
11 #include "TF1.h"
12 #include "TRandom.h"
13 #include "TLorentzVector.h"
14 #include "TTree.h"
15 #include "TVector3.h"
16 #include "THParticle.h" //detectors/hyp
17 #include "PndHypBupGenerator.h"
18 #include "FairPrimaryGenerator.h"
19 
20 using namespace std;
21 
22 
23 
24 // ----- Default constructor ------------------------------------------
26  iEvent = 0;
27  fInputFile = NULL;
28  fInputTree = NULL;
29 }
30 // ------------------------------------------------------------------------
31 
32 
33 
34 // ----- Standard constructor -----------------------------------------
35 PndHypBupGenerator::PndHypBupGenerator(const char* fileName) {
36 
37 
38  fFileName = fileName;
39  iEvent = 0;
40  fInputFile = new TFile(fFileName);
41  fInputTree = (TTree*) fInputFile->Get("data");
42  fParticles = new TClonesArray("THParticle",100);
43  fInputTree->SetBranchAddress("Particles", &fParticles);
44 
45 
46 
47 }
48 // ------------------------------------------------------------------------
49 
50 
51 
52 // ----- Destructor ---------------------------------------------------
54  CloseInput();
55 }
56 // ------------------------------------------------------------------------
57 
58 
59 
60 // ----- Public method ReadEvent --------------------------------------
61 Bool_t PndHypBupGenerator::ReadEvent(FairPrimaryGenerator* primGen) {
62 
63  // Check for input file
64  if ( ! fInputFile ) {
65  cout << "-E- PndHypBupGenerator: Input file not open!" << endl;
66  return kFALSE;
67  }
68  // Check for number of events in input file
69  if ( iEvent > fInputTree->GetEntries() ) {
70  cout << "-E PndDpmGenerator: No more events in input file!" << endl;
71  CloseInput();
72  return kFALSE;
73  }
74  TFile *g=gFile;
75  fInputFile->cd();
76  fInputTree->GetEntry(iEvent++);
77  g->cd();
78 
79  // Get number of particle in TClonesrray
80  Int_t nParts = fParticles->GetEntriesFast();
81 
82  // Loop over particles in TClonesArray
83 
84  for (Int_t iPart=0; iPart < nParts; iPart++) {
85  THParticle* part = (THParticle*) fParticles->At(iPart);
86  Int_t pdgType = part->GetPdgCode();
87 
88  // Check if particle type is known to database
89  if ( ! pdgType ) {
90  cout << "-W PndDpmGenerator: Unknown type " << part->GetPdgCode()
91  << ", skipping particle." << endl;
92  continue;
93  }
94 
95  Double_t px = part->Px();
96  Double_t py = part->Py();
97  Double_t pz = part->Pz();
98 
99  // Double_t vx = part->Vx();
100  // Double_t vy = part->Vy();
101  // Double_t vz = part->Vz();//default vertex 0,0,-76.5
102 
103  Double_t vx ;
104  Double_t vy ;
105  Double_t vz ;
106 
107  vx = part->Vx();
108  vy = part->Vy();
109  vz = part->Vz();
110 
111  /* cout<<" track "<<iPart<<" "
112  <<" pdgCode "<<pdgType<<" "
113  <<" energy "<<part->Energy()<<endl;*/
114 
115  //parametrisation of secondary vertex in second.target.
116  //reference 0.5 mm of silicon old(geo)
117 
118  /* TF1* Zdist = new TF1("Zdist","gaus(0)",-77,-73.5);
119  Zdist->SetParameters(10.1355,-75.51,0.76);
120  vz = Zdist->GetRandom();
121 
122  TF1* totaly = new TF1("totaly","gaus(0)+gaus(3)+gaus(6)",-2.,2);
123  totaly->SetParameters(8.9131,-0.635,0.6219,1.201,87.17,194,10.7,0.934,0.300);
124  TF1* totalx = new TF1("totalx","gaus(0)+gaus(3)+gaus(6)",-2.,2);
125  totalx->SetParameters(15.111,-1,-0.27,1,23,45,13,0.5,0.7);
126 
127  vx = totalx->GetRandom();
128  vy = totaly->GetRandom();
129  */
130  //if(pdgType==1020030090){
131 
132  cout<<" Secondary vertex for stopping Xi - "<<pdgType<<" "<<vx<<" "
133  <<vy<<" "<<vz<<endl;
134 
135  // Give track to PrimaryGenerator
136  //primGen->AddTrack(pdgType, px, py, pz, 0., 0., -76.5);
137  primGen->AddTrack(pdgType, px, py, pz, vx, vy,vz);
138  //}
139 
140 
141  }
142 
143 
144  return kTRUE;
145 }
146 // ------------------------------------------------------------------------
147 
148 
149 
150 // ----- Private method CloseInput ------------------------------------
152 
153  if ( fInputFile) {
154  cout << "-I- CbmHypBupGenerator: Closing input file "
155  << fFileName << endl;
156  fInputFile->Close();
157 
158  delete fInputFile;
159  }
160 
161  fInputFile = NULL;
162  }
163 
164 // ------------------------------------------------------------------------
165 
166 
167 
168 // ------------------------------------------------------------------------
170 
Double_t Vy() const
Definition: THParticle.h:70
Double_t Vz() const
Definition: THParticle.h:71
TFile * g
FairPrimaryGenerator * primGen
Definition: sim_emc_apd.C:81
virtual Bool_t ReadEvent(FairPrimaryGenerator *primGen)
Int_t fParticles
Definition: run_full.C:23
Double_t
Double_t Vx() const
Definition: THParticle.h:69
ClassImp(PndAnaContFact)
Int_t GetPdgCode() const
Definition: THParticle.h:60
double pz[39]
Definition: pipisigmas.h:14