KaliVeda  1.12/06
Heavy-Ion Analysis Toolkit
KVGemini.cpp
Go to the documentation of this file.
1 //Created by KVClassFactory on Fri Jul 25 10:05:03 2014
2 //Author: John Frankland,,,
3 #include "TTree.h"
4 
5 #include "KVGemini.h"
6 #include "KVSimEvent.h"
7 #include "KVSimNucleus.h"
8 
9 #include "CNucleus.h"
10 #include "CYrast.h"
11 
13 
14 
15 
16 //CYrast * KVGemini::yrast;
17 
18 
19 
26 KVGemini::KVGemini() : KVBase("gemini++", "Calculate statistical decay of excited nuclei")
27 {
28  // Default constructor
29  // enhanceIMF=true: call CWeight::setWeightIMF() before any decay
30  // (enhance the probabilty of IMF emission)
31  // in this case the event weights returned by
32 
33  // yrast = CYrast::instance(); //!< gives fission barriers and rotational energies
34 
35 }
36 
37 
38 
41 
43 {
44  // Destructor
45 }
46 
47 
48 
65 
66 void KVGemini::DecaySingleNucleus(KVSimNucleus& toDecay, KVSimEvent* decayProducts, bool addRotationalEnergy)
67 {
68  // Calculate decay products of excited nucleus toDecay
69  // Takes into account Z, A, E*, v and angular momentum of nucleus.
70  // Adds all decay products to decayProducts event:
71  // - call decayProducts->Clear() before calling this method
72  // if you just want to keep the products of a single nucleus
73  //
74  //
75  // \param[in] addRotationalEnergy enable or not the addition of the Yrast rotational energy to the excitation energy E*.
76  // Up to now, two cases are known:
77  // - HIPSE and DIT model => the excitation energy consists just in the thermal energy, the addition of the rotational energy is needed
78  // This is because Gemini systematically subtracts the Yrast energy from the given excitation energy of each nucleus
79  // in order to calculate the thermal excitation energy.
80  //
81  // If there is a problem with the decay of the nucleus,
82  // we throw an exception of type gemini_bad_decay
83 
84  CNucleus CN(toDecay.GetZ(), toDecay.GetA());
85  try {
86  // Adding rotationnal energy following K.Mazurek and M. Ciemala suggestion
87  // so if and J>0, E* is never 0 in CNucleus but follow the yrast line
88  // Could induce some problems for very exotic nuclei since yrast->getYrast()
89  // is defined only for a defined range of isotopes per element.
90 
91  // Info("DecaySingleNucleus", "Decaying: Z=%d A=%d E*=%g S=%g", toDecay.GetZ(), toDecay.GetA(), toDecay.GetExcitEnergy(), toDecay.GetAngMom().Mag());
92  if (addRotationalEnergy) CN.setCompoundNucleus(toDecay.GetExcitEnergy() + CYrast::instance()->getYrast(toDecay.GetZ(), toDecay.GetA(), toDecay.GetAngMom().Mag()), toDecay.GetAngMom().Mag());
93  else CN.setCompoundNucleus(toDecay.GetExcitEnergy(), toDecay.GetAngMom().Mag());
94  // set velocity
95  CN.setVelocityCartesian(toDecay.GetVelocity().X(), toDecay.GetVelocity().Y(), toDecay.GetVelocity().Z());
96 
97  // set angle of spin axis
98  CAngle ang(toDecay.GetAngMom().Theta(), toDecay.GetAngMom().Phi());
99  CN.setSpinAxis(ang);
100 
101  // Info("DecaySingleNucleus", "again Decaying: Z=%d A=%d E*=%g S=%g", toDecay.GetZ(), toDecay.GetA(), toDecay.GetExcitEnergy(), toDecay.GetAngMom().Mag());
102 
103  CN.decay();
104  // Info("DecaySingleNucleus", "decay done...");
105 
106  }
107  catch (std::exception& e) {
108  Info("DecaySingleNucleus", "Caught std::exception: %s", e.what());
109  Info("DecaySingleNucleus", "While decaying: Z=%d A=%d E*=%g S=%g", toDecay.GetZ(), toDecay.GetA(), toDecay.GetExcitEnergy(), toDecay.GetAngMom().Mag());
110  CN.reset();
111  throw gemini_bad_decay();
112  }
113  catch (...) {
114  Info("DecaySingleNucleus", "Caught unknown exception");
115  Info("DecaySingleNucleus", "While decaying: Z=%d A=%d E*=%g S=%g", toDecay.GetZ(), toDecay.GetA(), toDecay.GetExcitEnergy(), toDecay.GetAngMom().Mag());
116  CN.reset();
117  throw gemini_bad_decay();
118  }
119 
120  if (CN.abortEvent) { //problem with decay?
121  // throw exception on gemini++ error
122  CN.reset();
123  Info("DecaySingleNucleus", "Bad Gemini decay (CNucleus::abortEvent=true)");
124  Info("DecaySingleNucleus", "While decaying: Z=%d A=%d E*=%g S=%g", toDecay.GetZ(), toDecay.GetA(), toDecay.GetExcitEnergy(), toDecay.GetAngMom().Mag());
125  throw gemini_bad_decay();
126  }
127 
128  CNucleus* nuc = CN.getProducts(0);
129 
130  while (nuc) {
131 
132  KVNucleus* N = decayProducts->AddParticle();
133  N->SetParameter("GEMINI_PARENT_INDEX", part_index); //set parameter with index of parent nucleus in 'hot' event
134  N->SetZ(nuc->iZ);
135  N->SetA(nuc->iA);
136  N->SetExcitEnergy(nuc->fEx);
137  float* v = nuc->getVelocityVector();
138  TVector3 vv(v);
139  N->SetVelocity(vv);
140  nuc = CN.getProducts();
141 
142  }
143 
144  CN.reset();//<-- essential!
145 }
146 
147 
148 
160 
161 void KVGemini::DecayEvent(const KVSimEvent* hot, KVSimEvent* cold, bool addRotationalEnergy)
162 {
163  // Perform statistical decay of all nuclei in 'hot' event
164  // Takes into account Z, A, E*, v and spin of nuclei.
165  // Adds all decay products to 'cold' output event (the event is first reset, i.e. KVEvent::Clear()
166  // method is called, removing all particles AND parameters from the event).
167  //
168  // If there is a problem with the decay of any of the nuclei in the event,
169  // we throw an exception of type gemini_bad_decay
170  //
171  // \param[in] addRotationalEnergy defines whether or not to add Yrast rotational energy to E* of each nucleus (\sa DecaySingleNucleus() )
172  // \note Do not set any parameters for the `cold` event before calling this method, they will be cleared before decay
173 
174  cold->Clear();
175  KVSimNucleus* hotnuc;
176  part_index = 1;
177  while ((hotnuc = (KVSimNucleus*)const_cast<KVSimEvent*>(hot)->GetNextParticle())) {
178  try {
179  DecaySingleNucleus(*hotnuc, cold, addRotationalEnergy);
180  ++part_index;
181  }
182  catch (...) {
183  // rethrow any exceptions
184  throw;
185  }
186  }
187 }
188 
189 
190 
195 
196 void KVGemini::FillTreeWithEvents(KVSimNucleus& toDecay, bool addRotationalEnergy, Int_t nDecays, TTree* theTree, TString branchname)
197 {
198  // Perform nDecays decays of nucleus toDecay and write the events
199  // containing all decay products of each decay in theTree
200  // If no branchname is given, branch is called "gemini"
201 
202  if (branchname == "") branchname = "gemini";
203  KVSimEvent* decayProducts = new KVSimEvent;
204  KVEvent::MakeEventBranch(theTree, branchname, "KVSimEvent", decayProducts);
205 
206  while (nDecays--) {
207  decayProducts->Clear();
208  try {
209  DecaySingleNucleus(toDecay, decayProducts, addRotationalEnergy);
210  }
211  catch (exception& e) {
212  continue;
213  }
214  theTree->Fill();
215  std::cout << "\xd" << "Gemini++ processing, " << nDecays << " decays left ..." << std::flush;
216  }
217  std::cout << std::endl;
218 
219 }
220 
221 
222 
226 
228 {
229  // Maximum angular momentum with non-zero fission barrier (Sierk model)
230  // Only for 19<=Z<=111
231  if (z < 19 || z > 111) {
232  Error("GetMaxSpinWithFissionBarrier", "Only use for 19<=Z<=111");
233  return -1;
234  }
235  float amin = 1.2 * z + 0.01 * pow(z, 2);
236  float amax = 5.8 * z - 0.024 * pow(z, 2);
237 
238  if (a < amin || a > amax) {
239  Error("GetMaxSpinWithFissionBarrier", "Only valid for Z=%d with %d<=A<=%d",
240  z, (int)amin, (int)amax);
241  return -1;
242  }
243  return CYrast::instance()->getJmaxSierk(z, a);
244 }
245 
246 
247 
250 
252 {
253  // Return Rotating Liquid Drop Model fission barrier for given spin in hbar units
254  return CYrast::instance()->getBarrierFissionRLDM(z, a, J);
255 }
256 
257 
260 
262 {
263  // Return Sierk fission barrier for zero angular momentum
264  if (GetMaxSpinWithFissionBarrier(z, a) < 0) return -1;
265  return CYrast::instance()->getBarrierFissionSierk(0.);
266 }
267 
268 
int Int_t
ClassImp(KVPartitionList) void KVPartitionList
Initialisation.
#define e(i)
float Float_t
#define N
double pow(double, double)
Base class for KaliVeda framework.
Definition: KVBase.h:135
virtual void Clear(Option_t *opt="")
Definition: KVEvent.cpp:200
static void MakeEventBranch(TTree *tree, const TString &branchname, const TString &classname, T &event, Int_t bufsize=10000000)
Definition: KVEvent.h:468
KVNucleus * AddParticle()
Definition: KVEvent.cpp:166
Interface to GEMINI++.
Definition: KVGemini.h:36
void FillTreeWithEvents(KVSimNucleus &, bool, Int_t, TTree *, TString branchname="")
Definition: KVGemini.cpp:196
void DecaySingleNucleus(KVSimNucleus &, KVSimEvent *, bool)
Definition: KVGemini.cpp:66
Float_t GetMaxSpinWithFissionBarrier(int, int)
Definition: KVGemini.cpp:227
int part_index
Definition: KVGemini.h:38
Float_t GetFissionBarrierSierk(int z, int a)
Return Sierk fission barrier for zero angular momentum.
Definition: KVGemini.cpp:261
Float_t GetFissionBarrierRLDM(int z, int a, float J)
Return Rotating Liquid Drop Model fission barrier for given spin in hbar units.
Definition: KVGemini.cpp:251
virtual ~KVGemini()
Destructor.
Definition: KVGemini.cpp:42
void DecayEvent(const KVSimEvent *, KVSimEvent *, bool)
Definition: KVGemini.cpp:161
Description of properties and kinematics of atomic nuclei.
Definition: KVNucleus.h:125
Double_t GetExcitEnergy() const
Definition: KVNucleus.h:285
Int_t GetA() const
Definition: KVNucleus.cpp:799
Int_t GetZ() const
Return the number of proton / atomic number.
Definition: KVNucleus.cpp:770
TVector3 GetVelocity() const
returns velocity vector in cm/ns units
Container class for simulated nuclei, KVSimNucleus.
Definition: KVSimEvent.h:15
Nucleus in a simulated event.
Definition: KVSimNucleus.h:20
const TVector3 * GetAngMom() const
Definition: KVSimNucleus.h:92
virtual void Error(const char *method, const char *msgfmt,...) const
virtual void Info(const char *method, const char *msgfmt,...) const
virtual Int_t Fill()
Double_t Z() const
Double_t Phi() const
Double_t Y() const
Double_t X() const
Double_t Mag() const
Double_t Theta() const
Exception(s) thrown by KVGemini.
Definition: KVGemini.h:59
v
auto * a