KaliVeda
Toolkit for HIC analysis
KVGemini Class Reference

Detailed Description

Interface to GEMINI++.

You can use this class as an interface to Gemini++ to perform statistical decay of excited nuclei. In case of problems with the decay (CHucleus::abortEvent==1) we throw a standard exception of type gemini_bad_decay. This should be caught in your code when using DecaySingleNucleus() or DecayEvent():

KVGemini gem;
KVSimEvent *hot, *cold;
try {
gem.DecayEvent(hot,cold);
}
catch (gemini_bad_decay& e) {
}
Interface to GEMINI++.
Definition: KVGemini.h:37
void DecayEvent(const KVSimEvent *, KVSimEvent *, bool=true)
Definition: KVGemini.cpp:173
Container class for simulated nuclei, KVSimNucleus.
Definition: KVSimEvent.h:22
Exception(s) thrown by KVGemini.
Definition: KVGemini.h:60

Definition at line 37 of file KVGemini.h.

#include <KVGemini.h>

Inheritance diagram for KVGemini:

Public Member Functions

 KVGemini ()
 static CYrast* yrast; More...
 
virtual ~KVGemini ()
 Destructor. More...
 
void DecayEvent (const KVSimEvent *, KVSimEvent *, bool=true)
 
void DecaySingleNucleus (KVSimNucleus &, KVSimEvent *, bool=true)
 
void FillTreeWithEvents (KVSimNucleus &, bool, Int_t, TTree *, TString branchname="")
 
Float_t GetFissionBarrierRLDM (int z, int a, float J)
 Return Rotating Liquid Drop Model fission barrier for given spin in hbar units. More...
 
Float_t GetFissionBarrierSierk (int z, int a)
 Return Sierk fission barrier for zero angular momentum. More...
 
Float_t GetMaxSpinWithFissionBarrier (int, int)
 
- Public Member Functions inherited from KVBase
 KVBase ()
 Default constructor. More...
 
 KVBase (const Char_t *name, const Char_t *title="")
 Ctor for object with given name and type. More...
 
 KVBase (const KVBase &)
 copy ctor More...
 
virtual ~ KVBase ()
 
virtual void Clear (Option_t *opt="")
 Clear object properties : name, type/title, number, label. More...
 
virtual void Copy (TObject &) const
 Make a copy of this object. More...
 
const Char_t * GetLabel () const
 
UInt_t GetNumber () const
 
UInt_t GetNumberOfObjects () const
 
virtual TObject * GetObject () const
 
virtual const Char_t * GetType () const
 
Bool_t HasLabel () const
 
virtual Bool_t IsCalled (const Char_t *name) const
 
Bool_t IsLabelled (const Char_t *l) const
 
virtual Bool_t IsType (const Char_t *typ) const
 
virtual void List ()
 
KVBaseoperator= (const KVBase &)
 copy assignment operator More...
 
virtual void Print (Option_t *option="") const
 
Double_t ProtectedGetX (const TF1 *func, Double_t val, int &status, Double_t xmin=0.0, Double_t xmax=0.0) const
 
void SetLabel (const Char_t *lab)
 
virtual void SetNumber (UInt_t num)
 
virtual void SetType (const Char_t *str)
 

Private Attributes

int part_index
 

Additional Inherited Members

- Public Types inherited from KVBase
enum  EKaliVedaBits { kIsKaliVedaObject = BIT(23) }
 
- Static Public Member Functions inherited from KVBase
static Bool_t AreEqual (Double_t x, Double_t y, Long64_t maxdif=1)
 Comparison between two 64-bit floating-point values. More...
 
static void BackupFileWithDate (const Char_t *path)
 
static void CombineFiles (const Char_t *file1, const Char_t *file2, const Char_t *newfilename, Bool_t keep=kTRUE)
 
static void Deprecated (const char *method, const char *advice)
 
static Bool_t FindClassSourceFiles (const Char_t *class_name, KVString &imp_file, KVString &dec_file, const Char_t *dir_name=".")
 
static Bool_t FindExecutable (TString &exec, const Char_t *path="$(PATH)")
 
static const Char_t * FindFile (const Char_t *search, TString &wfil)
 
static const Char_t * GetBINDIRFilePath (const Char_t *namefile="")
 
static const Char_t * GetDATABASEFilePath ()
 
static const Char_t * GetDATADIRFilePath (const Char_t *namefile="")
 
static Bool_t GetDataSetEnv (const Char_t *dataset, const Char_t *type, Bool_t defval)
 
static const Char_t * GetDataSetEnv (const Char_t *dataset, const Char_t *type, const Char_t *defval)
 
static Double_t GetDataSetEnv (const Char_t *dataset, const Char_t *type, Double_t defval)
 
static const Char_t * GetETCDIRFilePath (const Char_t *namefile="")
 
static const Char_t * GetExampleFilePath (const Char_t *library, const Char_t *namefile)
 Return full path to example file for given library (="KVMultiDet", "BackTrack", etc.) More...
 
static const Char_t * GetINCDIRFilePath (const Char_t *namefile="")
 
static const Char_t * GetKVBuildDate ()
 Returns KaliVeda build date. More...
 
static const Char_t * GetKVBuildDir ()
 Returns top-level directory used for build. More...
 
static const Char_t * GetKVBuildTime ()
 Returns KaliVeda build time. More...
 
static const Char_t * GetKVBuildType ()
 Returns KaliVeda build type (cmake build: Release, Debug, RelWithDebInfo, ...) More...
 
static const Char_t * GetKVBuildUser ()
 Returns username of person who performed build. More...
 
static const Char_t * GetKVSourceDir ()
 Returns top-level directory of source tree used for build. More...
 
static const Char_t * GetKVVersion ()
 Returns KaliVeda version string. More...
 
static const Char_t * GetLIBDIRFilePath (const Char_t *namefile="")
 
static const Char_t * GetListOfPlugins (const Char_t *base)
 
static const Char_t * GetListOfPluginURIs (const Char_t *base)
 
static const Char_t * GetPluginURI (const Char_t *base, const Char_t *plugin)
 
static void GetTempFileName (TString &base)
 
static const Char_t * GetTEMPLATEDIRFilePath (const Char_t *namefile="")
 
static const Char_t * GetWORKDIRFilePath (const Char_t *namefile="")
 
static const Char_t * gitBranch ()
 Returns git branch of sources. More...
 
static const Char_t * gitCommit ()
 Returns last git commit of sources. More...
 
static void InitEnvironment ()
 
static bool is_gnuinstall ()
 
static Bool_t IsThisAPlugin (const TString &uri, TString &base)
 
static TPluginHandler * LoadPlugin (const Char_t *base, const Char_t *uri="0")
 
static Bool_t OpenContextMenu (const char *method, TObject *obj, const char *alt_method_name="")
 
static void OpenTempFile (TString &base, std::ofstream &fp)
 
static void PrintSplashScreen ()
 Prints welcome message and infos on version etc. More...
 
static Bool_t SearchAndOpenKVFile (const Char_t *name, KVSQLite::database &dbfile, const Char_t *kvsubdir="")
 
static Bool_t SearchAndOpenKVFile (const Char_t *name, std::ifstream &file, const Char_t *kvsubdir="", KVLockfile *locks=0)
 
static Bool_t SearchAndOpenKVFile (const Char_t *name, std::ofstream &file, const Char_t *kvsubdir="", KVLockfile *locks=0)
 
static Bool_t SearchKVFile (const Char_t *name, TString &fullpath, const Char_t *kvsubdir="")
 
static const Char_t * WorkingDirectory ()
 

Constructor & Destructor Documentation

◆ KVGemini()

KVGemini::KVGemini ( )

static CYrast* yrast;

Default constructor enhanceIMF=true: call CWeight::setWeightIMF() before any decay (enhance the probabilty of IMF emission) in this case the event weights returned by yrast = CYrast::instance(); //!< gives fission barriers and rotational energies

Definition at line 26 of file KVGemini.cpp.

◆ ~KVGemini()

KVGemini::~KVGemini ( )
virtual

Destructor.

Definition at line 42 of file KVGemini.cpp.

Member Function Documentation

◆ DecayEvent()

void KVGemini::DecayEvent ( const KVSimEvent hot,
KVSimEvent cold,
bool  addRotationalEnergy = true 
)

Perform statistical decay of all nuclei in 'hot' event Takes into account Z, A, E*, v and spin of nuclei. Adds all decay products to 'cold' output event (the event is first reset, i.e. KVEvent::Clear() method is called, removing all particles AND parameters from the event).

If there is a problem with the decay of any of the nuclei in the event, we throw an exception of type gemini_bad_decay

Parameters
[in]addRotationalEnergydefines whether or not to add Yrast rotational energy to E* of each nucleus
See also
DecaySingleNucleus()
Note
Do not set any parameters for the cold event before calling this method, they will be cleared before decay

Definition at line 173 of file KVGemini.cpp.

◆ DecaySingleNucleus()

void KVGemini::DecaySingleNucleus ( KVSimNucleus toDecay,
KVSimEvent decayProducts,
bool  addRotationalEnergy = true 
)

Calculate decay products of excited nucleus toDecay Takes into account Z, A, E*, v and angular momentum of nucleus. Adds all decay products to decayProducts event:

  • call decayProducts->Clear() before calling this method if you just want to keep the products of a single nucleus
Parameters
[in]addRotationalEnergy

enable or not the addition of the Yrast rotational energy to the excitation energy E*.

Up to now, two cases are known:

  • HIPSE and DIT model => the excitation energy consists just in the thermal energy, the addition of the rotational energy is needed

This is done because Gemini systematically subtracts the Yrast energy from the given excitation energy of each nucleus in order to calculate the thermal excitation energy.

If there is a problem with the decay of the nucleus, we throw an exception of type gemini_bad_decay

Definition at line 71 of file KVGemini.cpp.

◆ FillTreeWithEvents()

void KVGemini::FillTreeWithEvents ( KVSimNucleus toDecay,
bool  addRotationalEnergy,
Int_t  nDecays,
TTree *  theTree,
TString  branchname = "" 
)

Perform nDecays decays of nucleus toDecay and write the events containing all decay products of each decay in theTree If no branchname is given, branch is called "gemini"

Definition at line 209 of file KVGemini.cpp.

◆ GetFissionBarrierRLDM()

Float_t KVGemini::GetFissionBarrierRLDM ( int  z,
int  a,
float  J 
)

Return Rotating Liquid Drop Model fission barrier for given spin in hbar units.

Definition at line 264 of file KVGemini.cpp.

◆ GetFissionBarrierSierk()

Float_t KVGemini::GetFissionBarrierSierk ( int  z,
int  a 
)

Return Sierk fission barrier for zero angular momentum.

Definition at line 274 of file KVGemini.cpp.

◆ GetMaxSpinWithFissionBarrier()

Float_t KVGemini::GetMaxSpinWithFissionBarrier ( int  z,
int  a 
)

Maximum angular momentum with non-zero fission barrier (Sierk model) Only for 19<=Z<=111

Definition at line 240 of file KVGemini.cpp.

Member Data Documentation

◆ part_index

int KVGemini::part_index
private

used for labelling decay products

Definition at line 39 of file KVGemini.h.