GammaWare  Head Version for release 0.9
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
BaseGEM.h
Go to the documentation of this file.
1 /***************************************************************************
2  * Copyright (C) 2004 by Olivier Stezowski *
3  * stezow(AT)ipnl.in2p3.fr *
4  * *
5  * This program is free software; you can redistribute it and/or modify *
6  * it under the terms of the GNU General Public License as published by *
7  * the Free Software Foundation; either version 2 of the License, or *
8  * (at your option) any later version. *
9  * *
10  * This program is distributed in the hope that it will be useful, *
11  * but WITHOUT ANY WARRANTY; without even the implied warranty of *
12  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
13  * GNU General Public License for more details. *
14  * *
15  * You should have received a copy of the GNU General Public License *
16  * along with this program; if not, write to the *
17  * Free Software Foundation, Inc., *
18  * 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA. *
19  ***************************************************************************/
20 
23 #ifndef GW_BASEGEM_H
24 #define GW_BASEGEM_H
25 
26 #include <RandObj.h>
27 #include <LevelScheme.h>
28 #include <Link.h>
29 
30 class TLorentzVector;
31 
32 namespace Gw {
33 
107  class BaseGEM : public TObject
108  {
109  public:
110  enum EDirection { kDown, kUp } ;
111 
112  private:
113  LevelScheme *fLevelScheme; // level scheme on which is performed the Monte-Carlo
114 
115  private:
116  TObjArray fFeeding; // list to store all the constructed feeding links
117  TObjArray fRandFeeding; // list to store the random selectors for all the feeding links
118  TObjArray fRandDown; // list to store the random selectors for all the links of the level scheme
119  TObjArray fAD; // list of functions used for angular distributions
120 
121  private:
122  RandObj fRand0; // random selector to determine the starting point (one of the feeding links)
123 
124  private:
125  TList fTmp; // temporary list that owns objects to be deleted by BaseGEM
126  TObject fEoC; // EoC (End of Cascade). The random cascade stops when this particular link is reached.
127 
128  protected:
129  Int_t fDirection; // Simulated cascade up or down.
130 
131  protected:
132  Double_t fMinEnergyFactor; // minimum intensity of a link with a strength equal to 0
133  Int_t fADType; // 0 means isotropic, 1 means angular distribution from pure m-substate, 2 with a gaussian-distribution driven by sigma
134  Double_t fCutLifeTime; // under such life time, angular distribution is isotropic
135  Double_t fSigma; // sigma of the gaussian for orientation paramters. 0 by default
136 
137  protected:
139  virtual void InitRandom();
141  virtual void DoAngularDistribution(Int_t which_gamma, TLorentzVector &, Bool_t forceiso = false);
142 
143  public:
144  BaseGEM();
145  virtual ~BaseGEM();
146 
148  const TObjArray &GetFeedings() const
149  {
150  return fFeeding;
151  }
152 
154  Bool_t SetParameter(const char *name, Int_t value);
155  Bool_t SetParameter(const char *name, Double_t value);
156 
158 
161  virtual void SetDirection(EDirection d = kDown)
162  {
163  fDirection = d;
164  }
165 
167 
172  //Float_t SetMinFactor(Float_t newmin) { Float_t tmp = fMinEnergyFactor ; fMinEnergyFactor = newmin; return tmp; }
173 
175  const LevelScheme *GetLS() const
176  {
177  return fLevelScheme;
178  }
179 
181 
189  virtual Int_t Import(const Char_t *, Option_t *opt="152DY");
190 
192 
195  virtual void Clear(Option_t *opt = "all");
196 
198 
209  virtual Int_t DoCascade(TSeqCollection &cascade, Option_t *opt = "");
211 
217  virtual Int_t DoCascade(TSeqCollection &cascade, TSeqCollection &directions, Option_t *opt = "");
218 
220 
223  virtual Int_t DoCascade(TSeqCollection &cas, Int_t from, Option_t *opt = "");
224  virtual Int_t DoCascade(TSeqCollection &cas, TSeqCollection &directions, Int_t from, Option_t *opt = "");
225 
226 
228  // virtual void FillRandom(TH1 *h, Int_t ntimes = 5000);
229 
231  ClassDef(BaseGEM,0); // Base for a Monte-Carlo on a Level Scheme
232  };
233 
234 
235 }
236 #endif
237 
A level Scheme.
Definition: LevelScheme.h:82
Double_t fCutLifeTime
Definition: BaseGEM.h:134
virtual void Clear(Option_t *opt="all")
Clear everything.
Definition: BaseGEM.cpp:182
Class to get randomly cascades of gammas on the basis of a level scheme.
Definition: BaseGEM.h:107
Int_t fADType
Definition: BaseGEM.h:133
virtual ~BaseGEM()
Definition: BaseGEM.cpp:172
A class to select randomly an object in a TObjArray of objects.
Definition: RandObj.h:55
UInt_t value[MaxValue]
Definition: ReadDaqAlone.C:29
virtual void DoAngularDistribution(Int_t which_gamma, TLorentzVector &, Bool_t forceiso=false)
Definition: BaseGEM.cpp:536
header file for a RandObj
Double_t fMinEnergyFactor
Definition: BaseGEM.h:132
Int_t fDirection
Definition: BaseGEM.h:129
Double_t fSigma
Definition: BaseGEM.h:135
const LevelScheme * GetLS() const
To set the intensity of the lowest link of the level scheme.
Definition: BaseGEM.h:175
header file for a LevelScheme
const TObjArray & GetFeedings() const
To get all the entry point in this level scheme with their intensities.
Definition: BaseGEM.h:148
virtual Int_t DoCascade(TSeqCollection &cascade, Option_t *opt="")
cascade simulation
Definition: BaseGEM.cpp:523
virtual void SetDirection(EDirection d=kDown)
To determine the order of the simulate cascade.
Definition: BaseGEM.h:161
virtual Int_t Import(const Char_t *, Option_t *opt="152DY")
Load the level scheme and init the Monte-Carlo.
Definition: BaseGEM.cpp:194
Bool_t SetParameter(const char *name, Int_t value)
to set some parameters that modify the way the simulation is performed. return false is the parameter...
Definition: BaseGEM.cpp:227
ClassDef(BaseGEM, 0)
Fill the histogram with a random distribution corresponding to the first selected Link...
virtual void InitRandom()
to be called to init the random generator
Definition: BaseGEM.cpp:251