/*************************************************************************** * Copyright (C) 2004 by Olivier Stezowski * * stezow(AT)ipnl.in2p3.fr * * * * This program is free software; you can redistribute it and/or modify * * it under the terms of the GNU General Public License as published by * * the Free Software Foundation; either version 2 of the License, or * * (at your option) any later version. * * * * This program is distributed in the hope that it will be useful, * * but WITHOUT ANY WARRANTY; without even the implied warranty of * * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the * * GNU General Public License for more details. * * * * You should have received a copy of the GNU General Public License * * along with this program; if not, write to the * * Free Software Foundation, Inc., * * 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA. * ***************************************************************************/ /** \file BaseGEM.h header file for a BaseGEM (gamma-rays generator) */ #ifndef GW_BASEGEM_H #define GW_BASEGEM_H #include #include #include class TLorentzVector; namespace Gw { /** Class to get randomly cascades of gammas on the basis of a level scheme A BaseGEM wraps a LevelScheme in order to produce randomly cascades of gamma-rays on the basis of the level scheme. BE SURE OF YOUR LEVEL SCHEME BEFORE RUNNING IT !! \n The first version of the BaseGEM class (Release 0.2 and 0.4, svn revision < 188 ) takes into account the lateral feeding ONLY to stop the random process. A side effect of such a 'simple' treatment has been found by Gabriel and gives spectra with relative intensities between peaks that are not what is expected from the level scheme. To properly simulate cascades of gamma-rays, a new implementation has been developed. It DOES NOT change anything if you use BaseGEM (public interface), except the order of the links in the produced cascade which now by default goes from the top to the bottom contrary to the previous implementation. If your program needs absolutely the previous behavior, it is possible just by calling SetDirection(kUp). If you have developed a class that inherits from BaseGEM, there are important modifications of the protected interface. If you have troubles, please asked (agata-AT-ipnl.in2p3.fr) to find solutions to help you resolving them. In the new implementation, once the level scheme has been initialised with Import(), the whole level scheme is processed to determine the 'entry points' i.e. the lateral feedings with an intensity > 0. A list of Gw::Link is produced and is available through GetFeedings(). The random procedure consists then in: -# a random selection of one of the feeding link on the basis of the calculated intensities (not added to the cascade) -# from this starting point, the cascade down is determined. -# if kUp has been set, the cascade is inverted. Here is the illustration with a 'factice' level scheme and the produced spectra (sum, and two gated to show the correlation are reproduced). \image html simu_ags.gif \image html simu1_sp1.gif \deprecated (Release 0.2 and 0.4, svn revision < 188 ) The first step is to randomly select a gamma-ray on the basis on the relative intensities of all gamma-rays in the level scheme. One cannot starts from the bottom of the level scheme for several reasons: - there may be several gammas feeding the ground state - there may be several "ground state" corresponding to levels that no decay by gamma - there may be some isomeric states .. in this case we can simulate a lost by giving a smaller intensity below the isomeric state - there may be strongly converted gamma rays at the beginning of the level scheme Once a gamma-ray has been selected, the path down and up in the level scheme is generated based on the relative intensities of the links feeding (desexciting) the initial (final) link's level. This generator takes into account lateral feeding from the continuum or the lost of intensity due to isomeric or converted gammas. \n For each level, a lateral feeding is calculated as the difference between all the intensity that arrives on a level compared to all the intensity that leaves a level. It could be then from 0 (no lateral feeding) to 100 % (the highest discret gamma in the cascade).\n If in the random process up in the level scheme a lateral feeding is found, the process random stops.\n For each level, a probability to decay to something else that a gamma in the level scheme is calculated as the difference between all the intensity that leave a level compare to all the intensity that arrive on a level. It goes then from 0 (the next gamma belongs to the level scheme) to 100% (last level, decay by electron, long life time isomeric state ...). The Monte-Carlo process down stops when one of this situation is reached. \todo initialisation with a existing level scheme by copying it (copy constructor for level scheme) \todo Check out what happens when a link is not bind to two levels \remarks{ for the developers: fRandUp (fRandDown) is a TObjArray. With fRand0, all gamma-rays of the level scheme are associated to a unique integer corresponding to their slot position in the underlying TObjArray. Concerning fRandUp (fRandDown), each slot number i corresponds to a RandObj in which is stored the list of links needed to go up (down) after the link labelled number i in fRand0. Besides storing a list of links, in each slot of fRandUp (fRandDown) is stored a collection of TObject and their uniqueID (see ROOT TObject) are set to establish a correspondance with a link in the fRand0 TObjArray. Those objects are stored in a TList fTmp that owms tehm and thus will delete them. } \author Olivier Stezowski */ class BaseGEM : public TObject { public: enum EDirection { kDown, kUp } ; private: LevelScheme *fLevelScheme; // level scheme on which is performed the Monte-Carlo private: TObjArray fFeeding; // list to store all the constructed feeding links TObjArray fRandFeeding; // list to store the random selectors for all the feeding links TObjArray fRandDown; // list to store the random selectors for all the links of the level scheme TObjArray fAD; // list of functions used for angular distributions private: RandObj fRand0; // random selector to determine the starting point (one of the feeding links) private: TList fTmp; // temporary list that owns objects to be deleted by BaseGEM TObject fEoC; // EoC (End of Cascade). The random cascade stops when this particular link is reached. protected: Int_t fDirection; // Simulated cascade up or down. protected: Double_t fMinEnergyFactor; // minimum intensity of a link with a strength equal to 0 Int_t fADType; // 0 means isotropic, 1 means angular distribution from pure m-substate, 2 with a gaussian-distribution driven by sigma Double_t fCutLifeTime; // under such life time, angular distribution is isotropic Double_t fSigma; // sigma of the gaussian for orientation paramters. 0 by default protected: //! to be called to init the random generator virtual void InitRandom(); //! virtual void DoAngularDistribution(Int_t which_gamma, TLorentzVector &, Bool_t forceiso = false); public: BaseGEM(); BaseGEM(LevelScheme *); virtual ~BaseGEM(); //! To get all the entry point in this level scheme with their intensities const TObjArray &GetFeedings() const { return fFeeding; } //! to set some parameters that modify the way the simulation is performed. return false is the parameter is not known Bool_t SetParameter(const char *name, Int_t value); Bool_t SetParameter(const char *name, Double_t value); //! To determine the order of the simulate cascade /*! Default is kDown which means the first GammaLink in the cascade corresponds to the entry point. */ virtual void SetDirection(EDirection d = kDown) { fDirection = d; } //! To set the intensity of the lowest link of the level scheme /*! In level schemes, the intensity of a gamma-ray might be 0 ... which is a problem for the generator. By default, all links found with no intensity when Import is called are set with an intensity equal to the lowest found / MinFactor with MinFactor = 10. You can modifiy this factor with this method. */ //Float_t SetMinFactor(Float_t newmin) { Float_t tmp = fMinEnergyFactor ; fMinEnergyFactor = newmin; return tmp; } //! returns the level scheme embedded in this class const LevelScheme *GetLS() const { return fLevelScheme; } //! Load the level scheme and init the Monte-Carlo /*! The level scheme is loaded and, in case of success, the InitRandom method is called to initialised the random generator. If you creat a class that inherits from BaseGEM, in principle, you just need to implement correctly the InitRandom method. An example could be found in GEM It allows to load ags / ens / root level schemes - for ROOT first arg is the name of the file, the second one is the name of the level in the root file */ virtual Int_t Import(const Char_t *, Option_t *opt="152DY"); //! Clear everything /*! Default ("all"): clear the random generator and the level scheme. All other options ("") will reset the generator but not the level scheme */ virtual void Clear(Option_t *opt = "all"); //! cascade simulation /*! Fill the collection with the Links of the randomly determined path. The 'entry point' is determined randomly. If you would like to fix the enty point, you should use DoCascade(TSeqCollection &cas, Int_t from, Option_t *opt = "") The cascade is cleared before running the Monte-Carlo unless the option contains add. \n The sequence is ordered from the top to the bottom. Thus the first link gives the entry point in the level scheme. The order could be modified with SetDirection() It returns the number of links added to the cascade */ virtual Int_t DoCascade(TSeqCollection &cascade, Option_t *opt = ""); //! cascade simulation with angular distributions /*! Same that DoCascade(TSeqCollection &cascade, Option_t *opt = "") but simulate also the emission direction the directions collection contains TLorentzVector s with the spacial part being used for directions and t potentially for life time By default isotropic direction. This could be modified to take into account multipolarities through the Options [SetParameter] of BaseGEM */ virtual Int_t DoCascade(TSeqCollection &cascade, TSeqCollection &directions, Option_t *opt = ""); //! cascade simulation with the 'entry point' given by the user /*! GetFeedings() gives you the list of all the feeding links. from corresponds to the slot in this array of the feeding link you would like to start with. */ virtual Int_t DoCascade(TSeqCollection &cas, Int_t from, Option_t *opt = ""); virtual Int_t DoCascade(TSeqCollection &cas, TSeqCollection &directions, Int_t from, Option_t *opt = ""); //! Fill the histogram with a random distribution corresponding to the first selected Link // virtual void FillRandom(TH1 *h, Int_t ntimes = 5000); //! rootcint dictionary ClassDef(BaseGEM,0); // Base for a Monte-Carlo on a Level Scheme }; } #endif