KaliVeda  1.13/01
Heavy-Ion Analysis Toolkit
KVGenPhaseSpace.cpp
Go to the documentation of this file.
1 //Created by KVClassFactory on Fri Apr 17 10:19:02 2015
2 //Author: John Frankland,,,
3 
4 #include "KVGenPhaseSpace.h"
5 #include "TGenPhaseSpace.h"
6 
8 
9 
10 
11 
12 
15 void KVGenPhaseSpace::init()
16 {
17  // default initialisations
18 
19  fMCSampler = 0;
20  fEvent = 0;
21  fOK = kFALSE;
22 }
23 
24 
25 
31 
33 {
34  // initialise the TGenPhaseSpace for the decay channel
35  // We store the name of the generator in the list of parameters of the
36  // KVEvent object pointed to by fEvent:
37  // fEvent->GetParameters()->SetValue("PHASESPACE_GENERATOR","TGenPhaseSpace");
38 
41  fEvent->GetParameters()->SetValue("PHASESPACE_GENERATOR", "TGenPhaseSpace");
42 }
43 
44 
45 
48 
50 {
51  // Default constructor
52  init();
53 }
54 
55 
56 
57 
60 
61 KVGenPhaseSpace::KVGenPhaseSpace(const Char_t* name, const Char_t* title) : KVBase(name, title)
62 {
63  // Constructor with name and title
64  init();
65 }
66 
67 
68 
71 
73 {
74  // Destructor
76 }
77 
78 
79 
80 
88 
90 {
91  // This method copies the current state of 'this' object into 'obj'
92  // You should add here any member variables, for example:
93  // (supposing a member variable KVGenPhaseSpace::fToto)
94  // CastedObj.fToto = fToto;
95  // or
96  // CastedObj.SetToto( GetToto() );
97 
98  KVBase::Copy(obj);
99  //KVGenPhaseSpace& CastedObj = (KVGenPhaseSpace&)obj;
100 }
101 
102 
103 
110 
112 {
113  // Define the break-up channel, i.e. the initial compound nucleus with its
114  // kinematics & excitation energy, and a list of nuclei to be produced
115  // by the decay. Nuclei produced by decay may be excited.
116  //
117  // Returns kTRUE if decay is energetically allowed, kFALSE if not.
118 
119  fEvent = e;
120  fCompound = CN;
121  fMult = e->GetMult();
122 
124  if (!fOK) return fOK;
125 
126  fMasses.clear();
127  fMasses.reserve(fMult);
129 
131 
132  return fOK;
133 }
134 
135 
136 
147 
149 {
150  // Generate 1 event for this break-up channel.
151  // The weight associated to the event is returned.
152  // We also store the value in the list of parameters of the
153  // KVEvent object pointed to by fEvent:
154  // fEvent->GetParameters()->SetValue("PHASESPACE_WEIGHT", weight);
155  //
156  // NB: in order to get e.g. correct kinetic energy distributions for particles,
157  // this weight MUST be used as e.g. the filling weight for any histograms
158  // (see example in Class Description)
159 
160  if (!fOK) {
161  Warning("Generate", "Generator is not initialised correctly.");
162  return 0.;
163  }
164  Double_t weight = ((TGenPhaseSpace*)fMCSampler)->Generate();
165  fEvent->GetParameters()->SetValue("PHASESPACE_WEIGHT", weight);
166  for (int i = 0; i < fMult; i++) {
167  fEvent->GetParticle(i + 1)->Set4Mom(*(((TGenPhaseSpace*)fMCSampler)->GetDecay(i)));
168  }
169  return weight;
170 }
171 
172 
173 
179 
181 {
182  // We check that:
183  // - mass and charge are conserved
184  // - the excitation energy of the compound is greater than the channel Q-value
185  // (taking into account any excitation energy of nuclei in exit channel)
186 
187  Int_t ztot = (Int_t)fEvent->GetSum("GetZ");
188  Int_t atot = (Int_t)fEvent->GetSum("GetA");
189  if (ztot != fCompound.GetZ() || atot != fCompound.GetA()) {
190  Warning("CheckBreakUpChannel", "Compound has (Z,A)=(%d,%d), break-up channel has (Z,A)=(%d,%d)",
191  fCompound.GetZ(), fCompound.GetA(), ztot, atot);
192  return kFALSE;
193  }
194 
195  Double_t exmin = -fEvent->GetChannelQValue();
196  fEtot = fCompound.GetExcitEnergy() - exmin;
197  if (fEtot < 0.) {
198  Warning("CheckBreakUpChannel", "Excitation energy of compound must be > %.2lf MeV", exmin);
199  return kFALSE;
200  }
201 
202  return kTRUE;
203 }
204 
205 
int Int_t
ClassImp(KVPartitionList) void KVPartitionList
Initialisation.
#define SafeDelete(p)
#define e(i)
char Char_t
const Bool_t kFALSE
bool Bool_t
double Double_t
const Bool_t kTRUE
Base class for KaliVeda framework.
Definition: KVBase.h:141
virtual void Copy(TObject &) const
Make a copy of this object.
Definition: KVBase.cpp:394
Abstract base class container for multi-particle events.
Definition: KVEvent.h:66
KVNameValueList * GetParameters() const
Definition: KVEvent.h:203
virtual Double_t GetSum(const Char_t *, Option_t *="")=0
virtual KVParticle * GetParticle(Int_t npart) const =0
virtual Double_t GetChannelQValue() const =0
virtual void GetMasses(std::vector< Double_t > &)=0
virtual Int_t GetMult(Option_t *opt="") const
Definition: KVEvent.h:141
Generate momenta for an event using microcanonical phase space sampling.
Bool_t CheckBreakUpChannel()
Bool_t SetBreakUpChannel(const KVNucleus &CN, KVEvent *e)
Bool_t fOK
ready to calculate or not
KVNucleus fCompound
initial nucleus which undergoes break-up
virtual ~KVGenPhaseSpace()
Destructor.
KVEvent * fEvent
break-up channel
TObject * fMCSampler
Monte-Carlo phase space sampler.
void init()
default initialisations
std::vector< Double_t > fMasses
masses of nuclei in break-up channel
Double_t fEtot
available kinetic energy for decay
virtual Double_t Generate()
KVGenPhaseSpace()
Default constructor.
void Copy(TObject &obj) const
virtual void InitialiseMCSampler()
Int_t fMult
multiplicity of channel
void SetValue(const Char_t *name, value_type value)
Description of properties and kinematics of atomic nuclei.
Definition: KVNucleus.h:125
Double_t GetExcitEnergy() const
Definition: KVNucleus.h:282
Int_t GetA() const
Definition: KVNucleus.cpp:799
Int_t GetZ() const
Return the number of proton / atomic number.
Definition: KVNucleus.cpp:770
void Set4Mom(const TLorentzVector &p)
Definition: KVParticle.h:591
Bool_t SetDecay(TLorentzVector &P, Int_t nt, const Double_t *mass, Option_t *opt="")
virtual void Warning(const char *method, const char *msgfmt,...) const