CbmRoot
CbmUnigenGenerator.h
Go to the documentation of this file.
1 
6 #ifndef CBMUNIGENGENERATOR
7 #define CBMUNIGENGENERATOR
8 
9 
10 #include "FairGenerator.h"
11 #include "TMath.h"
12 #include "TString.h"
13 #include "TVector3.h"
14 #include <map>
15 
16 class TFile;
17 class TTree;
18 class UEvent;
19 class UParticle;
20 class FairIon;
21 class FairPrimaryGenerator;
22 
23 
44 class CbmUnigenGenerator : public FairGenerator {
45 
46 public:
48  enum EMode {
53  };
54 
55 
62  CbmUnigenGenerator(const char* fileName = "", EMode mode = kStandard);
63 
64 
70  CbmUnigenGenerator(const char* fileName, EMode mode, Double_t phi);
71 
72 
74  virtual ~CbmUnigenGenerator();
75 
76 
87  virtual Bool_t ReadEvent(FairPrimaryGenerator* primGen);
88 
89 
90 private:
91  TString fFileName;
93  Double_t fPhi;
94  Bool_t fIsInit;
95  TFile* fFile;
96  TTree* fTree;
97  Int_t fCurrentEntry;
99  Int_t fNofPrimaries;
100  Int_t fNofEvents;
101  Double_t fBetaCM;
102  Double_t fGammaCM;
103  std::map<TString, FairIon*> fIonMap;
104 
105  // Constants for decimal decomposition of ion PDG.
106  // For ions the PDG code is +-10LZZZAAAI, with L = number of Lambdas,
107  // ZZZ = charge (number of protons), AAA = mass (sum of numbers
108  // of Lambdas, protons and neutrons, I = isomer level
109  static const Int_t kPdgLambda = 10000000;
110  static const Int_t kPdgCharge = 10000;
111  static const Int_t kPdgMass = 10;
112 
113 
119  void AddPrimary(FairPrimaryGenerator* primGen,
120  Int_t pdgCode,
121  const TVector3& momentum);
122 
123 
125  void CloseInput();
126 
127 
134  Int_t GetIonCharge(Int_t pdgCode) const {
135  return (pdgCode % kPdgLambda) / kPdgCharge;
136  }
137 
138 
145  Int_t GetIonLambdas(Int_t pdgCode) const {
146  return (pdgCode % (10 * kPdgLambda)) / kPdgLambda;
147  }
148 
149 
156  Int_t GetIonMass(Int_t pdgCode) const {
157  return (pdgCode % kPdgCharge) / kPdgMass;
158  }
159 
160 
164  Bool_t GetNextEntry();
165 
166 
173  Bool_t Init();
174 
175 
186  void ProcessIon(FairPrimaryGenerator* primGen,
187  Int_t pdgCode,
188  const TVector3& momentum);
189 
190 
199  Int_t RegisterIons();
200 
201 
204 
205 
208 
209 
211 };
212 
213 #endif
CbmUnigenGenerator::EMode
EMode
Mode enumerator.
Definition: CbmUnigenGenerator.h:48
CbmUnigenGenerator::ProcessIon
void ProcessIon(FairPrimaryGenerator *primGen, Int_t pdgCode, const TVector3 &momentum)
Treat a composite particle (ion)
Definition: CbmUnigenGenerator.cxx:216
CbmUnigenGenerator::GetIonLambdas
Int_t GetIonLambdas(Int_t pdgCode) const
Number of Lambdas in an ion.
Definition: CbmUnigenGenerator.h:145
CbmUnigenGenerator::fBetaCM
Double_t fBetaCM
CM velocity in the lab frame.
Definition: CbmUnigenGenerator.h:101
CbmUnigenGenerator::AddPrimary
void AddPrimary(FairPrimaryGenerator *primGen, Int_t pdgCode, const TVector3 &momentum)
Add a primary particle to the event generator.
Definition: CbmUnigenGenerator.cxx:89
UParticle
Definition: UParticle.h:10
CbmUnigenGenerator::fEvent
UEvent * fEvent
Current input event.
Definition: CbmUnigenGenerator.h:98
CbmUnigenGenerator::Init
Bool_t Init()
Initialisation.
Definition: CbmUnigenGenerator.cxx:138
CbmUnigenGenerator::fIonMap
std::map< TString, FairIon * > fIonMap
Map from ion name to FairIon.
Definition: CbmUnigenGenerator.h:103
CbmUnigenGenerator::fFile
TFile * fFile
Input ROOT file.
Definition: CbmUnigenGenerator.h:95
CbmUnigenGenerator
Generates input to transport simulation from files in Unigen format.
Definition: CbmUnigenGenerator.h:44
CbmUnigenGenerator::fNofEvents
Int_t fNofEvents
Number of processed events.
Definition: CbmUnigenGenerator.h:100
CbmUnigenGenerator::kStandard
@ kStandard
Rotate events to zero event plane (default)
Definition: CbmUnigenGenerator.h:49
CbmUnigenGenerator::operator=
CbmUnigenGenerator & operator=(const CbmUnigenGenerator &)=delete
Assignment operator forbidden.
CbmUnigenGenerator::kRotateFixed
@ kRotateFixed
Rotate events around z by a fixed angle.
Definition: CbmUnigenGenerator.h:51
UEvent
Definition: UEvent.h:12
CbmUnigenGenerator::ReadEvent
virtual Bool_t ReadEvent(FairPrimaryGenerator *primGen)
Read one event from the input file.
Definition: CbmUnigenGenerator.cxx:249
CbmUnigenGenerator::fNofPrimaries
Int_t fNofPrimaries
Number of primaries registered in current event.
Definition: CbmUnigenGenerator.h:99
CbmUnigenGenerator::fFileName
TString fFileName
Input file name.
Definition: CbmUnigenGenerator.h:91
CbmUnigenGenerator::kPdgMass
static const Int_t kPdgMass
Decomposition of ion PDG code.
Definition: CbmUnigenGenerator.h:111
CbmUnigenGenerator::kReuseEvents
@ kReuseEvents
Reuse events if more are requested than present.
Definition: CbmUnigenGenerator.h:52
CbmUnigenGenerator::RegisterIons
Int_t RegisterIons()
Register ions to the simulation.
Definition: CbmUnigenGenerator.cxx:319
CbmUnigenGenerator::kPdgLambda
static const Int_t kPdgLambda
Decomposition of ion PDG code.
Definition: CbmUnigenGenerator.h:109
CbmUnigenGenerator::fMode
EMode fMode
Rotation mode.
Definition: CbmUnigenGenerator.h:92
CbmUnigenGenerator::fPhi
Double_t fPhi
Event plane rotation angle.
Definition: CbmUnigenGenerator.h:93
CbmUnigenGenerator::kNoRotation
@ kNoRotation
No event rotation.
Definition: CbmUnigenGenerator.h:50
CbmUnigenGenerator::fGammaCM
Double_t fGammaCM
Gamma factor of CM in lab frame.
Definition: CbmUnigenGenerator.h:102
CbmUnigenGenerator::GetNextEntry
Bool_t GetNextEntry()
Get next entry from input tree.
Definition: CbmUnigenGenerator.cxx:110
CbmUnigenGenerator::fTree
TTree * fTree
Input ROOT tree.
Definition: CbmUnigenGenerator.h:96
CbmUnigenGenerator::CbmUnigenGenerator
CbmUnigenGenerator(const CbmUnigenGenerator &)=delete
Copy constructor forbidden.
CbmUnigenGenerator::CbmUnigenGenerator
CbmUnigenGenerator(const char *fileName="", EMode mode=kStandard)
Default constructor.
Definition: CbmUnigenGenerator.cxx:27
CbmUnigenGenerator::kPdgCharge
static const Int_t kPdgCharge
Decomposition of ion PDG code.
Definition: CbmUnigenGenerator.h:110
CbmUnigenGenerator::fCurrentEntry
Int_t fCurrentEntry
Current entry number.
Definition: CbmUnigenGenerator.h:97
CbmUnigenGenerator::CloseInput
void CloseInput()
Close the input file.
Definition: CbmUnigenGenerator.cxx:100
CbmUnigenGenerator::fIsInit
Bool_t fIsInit
Flag whether generator is initialised.
Definition: CbmUnigenGenerator.h:94
CbmUnigenGenerator::~CbmUnigenGenerator
virtual ~CbmUnigenGenerator()
Destructor.
Definition: CbmUnigenGenerator.cxx:80
CbmUnigenGenerator::GetIonCharge
Int_t GetIonCharge(Int_t pdgCode) const
Charge number of an ion.
Definition: CbmUnigenGenerator.h:134
CbmUnigenGenerator::GetIonMass
Int_t GetIonMass(Int_t pdgCode) const
Mass number of an ion.
Definition: CbmUnigenGenerator.h:156
CbmUnigenGenerator::ClassDef
ClassDef(CbmUnigenGenerator, 4)