FairRoot
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
FairPrimaryGenerator.h
Go to the documentation of this file.
1 /********************************************************************************
2  * Copyright (C) 2014 GSI Helmholtzzentrum fuer Schwerionenforschung GmbH *
3  * *
4  * This software is distributed under the terms of the *
5  * GNU Lesser General Public Licence (LGPL) version 3, *
6  * copied verbatim in the file "LICENSE" *
7  ********************************************************************************/
8 
23 #ifndef FAIRPRIMARYGENERATOR_H
24 #define FAIRPRIMARYGENERATOR_H
25 
26 #include "FairGenerator.h" // for FairGenerator
27 
28 #include <Rtypes.h> // for Double_t, Bool_t, Int_t, etc
29 #include <TMCProcess.h>
30 #include <TNamed.h> // for TNamed
31 #include <TObjArray.h> // for TObjArray
32 #include <TVector3.h> // for TVector3
33 #include <iostream> // for cout
34 
35 class FairGenericStack;
36 class FairMCEventHeader;
37 class TIterator;
38 
39 class FairPrimaryGenerator : public TNamed
40 {
41 
42  public:
45 
47  FairPrimaryGenerator(const char *name, const char *title = "FAIR Generator");
48 
50  virtual ~FairPrimaryGenerator();
51 
53  virtual Bool_t Init();
54 
56  void AddGenerator(FairGenerator *generator)
57  {
58  if (!fGenList) {
59  std::cout << "Empty fGenList pointer ! " << std::endl;
60  return;
61  }
62  fGenList->Add(generator);
63  }
64 
65  void SetEventNr(Int_t evtNr) { fEventNr = evtNr; }
66 
74  virtual Bool_t GenerateEvent(FairGenericStack *pStack);
75 
83  virtual void AddTrack(Int_t pdgid,
84  Double_t px,
85  Double_t py,
86  Double_t pz,
87  Double_t vx,
88  Double_t vy,
89  Double_t vz,
90  Int_t parent = -1,
91  Bool_t wanttracking = true,
92  Double_t e = -9e9,
93  Double_t tof = 0.,
94  Double_t weight = 0.,
95  TMCProcess proc = kPPrimary);
96 
99 
106  void SetBeam(Double_t beamX0, Double_t beamY0, Double_t beamSigmaX, Double_t beamSigmaY);
107 
114  void SetBeamAngle(Double_t beamAngleX0, Double_t beamAngleY0, Double_t beamAngleSigmaX, Double_t beamAngleSigmaY);
115 
122  void SetEventPlane(Double_t phiMin, Double_t phiMax);
123 
128  void SetTarget(Double_t targetZ, Double_t targetDz);
129 
136  void SetMultTarget(Int_t nroftargets, Double_t *targetZ, Double_t targetDz);
137 
139  void SmearVertexZ(Bool_t flag);
140  void SmearGausVertexZ(Bool_t flag);
141  void SmearVertexXY(Bool_t flag);
142  void SmearGausVertexXY(Bool_t flag);
143 
144  TObjArray *GetListOfGenerators() { return fGenList; }
145 
147  void SetEvent(FairMCEventHeader *event) { fEvent = event; };
148 
151 
154  void DoTracking(Bool_t doTracking = kTRUE) { fdoTracking = doTracking; }
155 
156  Int_t GetTotPrimary() { return fTotPrim; }
157 
158  protected:
163 
165  Double_t fBeamX0;
167  Double_t fBeamY0;
169  Double_t fBeamSigmaX;
171  Double_t fBeamSigmaY;
172 
174  Double_t fBeamAngleX0;
176  Double_t fBeamAngleY0;
178  Double_t fBeamAngleX;
180  Double_t fBeamAngleY;
186  TVector3 fBeamDirection;
187 
189  Double_t fPhiMin;
191  Double_t fPhiMax;
193  Double_t fPhi;
194 
196  Double_t *fTargetZ;
197 
198  Int_t fNrTargets;
200  Double_t fTargetDz;
201 
203  TVector3 fVertex;
204 
206  Int_t fNTracks;
207 
217  Bool_t fBeamAngle;
219  Bool_t fEventPlane;
220 
223 
224  TObjArray *fGenList;
226  TIterator *fListIter;
227 
229 
230  Bool_t fdoTracking;
231 
234 
235  static Int_t fTotPrim;
236 
239  Int_t fEventNr;
240 
250  virtual void MakeVertex();
251 
259  virtual void MakeBeamAngle();
260 
268  void MakeEventPlane();
269 
271 };
272 
273 #endif
void SmearVertexXY(Bool_t flag)
void AddGenerator(FairGenerator *generator)
void DoTracking(Bool_t doTracking=kTRUE)
void SetTarget(Double_t targetZ, Double_t targetDz)
ClassDef(FairPrimaryGenerator, 5)
virtual FairPrimaryGenerator * ClonePrimaryGenerator() const
void SetEventPlane(Double_t phiMin, Double_t phiMax)
void SetEventNr(Int_t evtNr)
void SmearGausVertexZ(Bool_t flag)
FairPrimaryGenerator & operator=(const FairPrimaryGenerator &)
void SetBeamAngle(Double_t beamAngleX0, Double_t beamAngleY0, Double_t beamAngleSigmaX, Double_t beamAngleSigmaY)
virtual void AddTrack(Int_t pdgid, Double_t px, Double_t py, Double_t pz, Double_t vx, Double_t vy, Double_t vz, Int_t parent=-1, Bool_t wanttracking=true, Double_t e=-9e9, Double_t tof=0., Double_t weight=0., TMCProcess proc=kPPrimary)
void SetBeam(Double_t beamX0, Double_t beamY0, Double_t beamSigmaX, Double_t beamSigmaY)
virtual Bool_t GenerateEvent(FairGenericStack *pStack)
void SetEvent(FairMCEventHeader *event)
void SmearVertexZ(Bool_t flag)
FairMCEventHeader * GetEvent()
FairGenericStack * fStack
FairMCEventHeader * fEvent
void SetMultTarget(Int_t nroftargets, Double_t *targetZ, Double_t targetDz)
TObjArray * GetListOfGenerators()
void SmearGausVertexXY(Bool_t flag)