G4EmModelManager.hh

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 // $Id$
00027 //
00028 // -------------------------------------------------------------------
00029 //
00030 // GEANT4 Class header file
00031 //
00032 // File name:     G4EmModelManager
00033 //
00034 // Author:        Vladimir Ivanchenko
00035 //
00036 // Creation date: 07.05.2002
00037 //
00038 // Modifications:
00039 //
00040 // 03-12-02 V.Ivanchenko fix a bug in model selection
00041 // 20-01-03 Migrade to cut per region (V.Ivanchenko)
00042 // 27-01-03 Make models region aware (V.Ivanchenko)
00043 // 13-02-03 The set of models is defined for region (V.Ivanchenko)
00044 // 26-03-03 Add GetDEDXDispersion (V.Ivanchenko)
00045 // 13-04-03 Add startFromNull (V.Ivanchenko)
00046 // 13-05-03 Add calculation of precise range (V.Ivanchenko)
00047 // 21-07-03 Add UpdateEmModel method (V.Ivanchenko)
00048 // 03-11-03 Substitute STL vector for G4RegionModels (V.Ivanchenko)
00049 // 11-04-05 Remove access to fluctuation models (V.Ivanchenko)
00050 // 10-01-06 PreciseRange -> CSDARange (V.Ivantchenko)
00051 // 20-01-06 Introduce G4EmTableType and reducing number of methods (VI)
00052 // 13-05-06 Add GetModel by index method (VI)
00053 // 15-03-07 Add maxCutInRange (V.Ivanchenko)
00054 // 08-04-08 Simplify Select method for only one G4RegionModel (VI)
00055 // 03-08-09 Removed unused members and simplify model search if only one
00056 //          model is used (VI)
00057 // 14-07-11 Use pointer to the vector of cuts and not local copy (VI)
00058 //
00059 // Class Description:
00060 //
00061 // It is the unified energy loss process it calculates the continuous
00062 // energy loss for charged particles using a set of Energy Loss
00063 // models valid for different energy regions. There are a possibility
00064 // to create and access to dE/dx and range tables, or to calculate
00065 // that information on fly.
00066 
00067 // -------------------------------------------------------------------
00068 //
00069 
00070 
00071 #ifndef G4EmModelManager_h
00072 #define G4EmModelManager_h 1
00073 
00074 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00075 
00076 #include "globals.hh"
00077 #include "G4DataVector.hh"
00078 #include "G4EmTableType.hh"
00079 #include "G4EmProcessSubType.hh"
00080 #include "G4Region.hh"
00081 
00082 class G4RegionModels
00083 {
00084 
00085 friend class G4EmModelManager;
00086 
00087 public:
00088 
00089 private:
00090 
00091   G4RegionModels(G4int nMod, std::vector<G4int>& indx, 
00092                  G4DataVector& lowE, const G4Region* reg);
00093 
00094   ~G4RegionModels();
00095 
00096   inline G4int SelectIndex(G4double e) const {
00097     G4int idx = 0;
00098     if (nModelsForRegion>1) {
00099       idx = nModelsForRegion;
00100       do {--idx;} while (idx > 0 && e <= lowKineticEnergy[idx]);
00101     }
00102     return theListOfModelIndexes[idx];
00103   };
00104 
00105   inline G4int ModelIndex(G4int n) const {
00106     return theListOfModelIndexes[n];
00107   };
00108 
00109   inline G4int NumberOfModels() const {
00110     return nModelsForRegion;
00111   };
00112 
00113   inline G4double LowEdgeEnergy(G4int n) const {
00114     return lowKineticEnergy[n];
00115   };
00116 
00117   inline const G4Region* Region() const {
00118     return theRegion;
00119   };
00120 
00121   G4RegionModels(G4RegionModels &);
00122   G4RegionModels & operator=(const G4RegionModels &right);
00123 
00124   const G4Region*    theRegion;
00125   G4int              nModelsForRegion;
00126   G4int*             theListOfModelIndexes;
00127   G4double*          lowKineticEnergy;
00128 
00129 };
00130 
00131 #include "G4VEmModel.hh"
00132 #include "G4VEmFluctuationModel.hh"
00133 #include "G4DynamicParticle.hh"
00134 
00135 class G4Region;
00136 class G4ParticleDefinition;
00137 class G4PhysicsVector;
00138 class G4MaterialCutsCouple;
00139 
00140 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00141 
00142 class G4EmModelManager
00143 {
00144 public:
00145 
00146   G4EmModelManager();
00147 
00148   ~G4EmModelManager();
00149 
00150   void Clear();
00151 
00152   const G4DataVector* Initialise(const G4ParticleDefinition*,
00153                                  const G4ParticleDefinition*,
00154                                        G4double,
00155                                        G4int);
00156 
00157   void FillDEDXVector(G4PhysicsVector*, const G4MaterialCutsCouple*, 
00158                       G4EmTableType t = fRestricted);
00159 
00160   void FillLambdaVector(G4PhysicsVector*, const G4MaterialCutsCouple*, 
00161                         G4bool startFromNull = true, G4EmTableType t = fRestricted);
00162 
00163   G4VEmModel* GetModel(G4int, G4bool ver = false);
00164 
00165   void AddEmModel(G4int, G4VEmModel*, G4VEmFluctuationModel*, const G4Region*);
00166   
00167   void UpdateEmModel(const G4String&, G4double, G4double);
00168 
00169   void DumpModelList(G4int verb);
00170 
00171   inline G4VEmModel* SelectModel(G4double& energy, size_t& index);
00172 
00173   inline const G4DataVector* Cuts() const;
00174 
00175   inline const G4DataVector* SubCutoff() const;
00176 
00177   inline void SetFluoFlag(G4bool val);
00178 
00179   inline G4int NumberOfModels() const;
00180 
00181 private:
00182 
00183   inline G4double ComputeDEDX(G4VEmModel* model,
00184                               const G4MaterialCutsCouple*,
00185                               G4double kinEnergy,
00186                               G4double cutEnergy,
00187                               G4double minEnergy);
00188 
00189   // hide  assignment operator
00190 
00191   G4EmModelManager(G4EmModelManager &);
00192   G4EmModelManager & operator=(const G4EmModelManager &right);
00193 
00194 // =====================================================================
00195 
00196 private:
00197 
00198   const G4DataVector*          theCuts;
00199   G4DataVector*                theSubCuts;
00200 
00201   std::vector<G4VEmModel*>                models;
00202   std::vector<G4VEmFluctuationModel*>     flucModels;
00203   std::vector<const G4Region*>            regions;
00204   std::vector<G4int>                      orderOfModels;
00205   std::vector<G4int>                      isUsed;
00206 
00207   G4int                       nEmModels;
00208   G4int                       nRegions;
00209 
00210   std::vector<G4int>            idxOfRegionModels;
00211   std::vector<G4RegionModels*>  setOfRegionModels;
00212 
00213   G4double                    maxSubCutInRange;
00214 
00215   const G4ParticleDefinition* particle;
00216 
00217   G4int                       verboseLevel;
00218   G4bool                      severalModels;
00219   G4bool                      fluoFlag;
00220 
00221   // may be changed in run time
00222   G4RegionModels*             currRegionModel;
00223   G4VEmModel*                 currModel;
00224 };
00225 
00226 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00227 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00228 
00229 inline G4VEmModel* G4EmModelManager::SelectModel(G4double& kinEnergy, 
00230                                                  size_t& index)
00231 {
00232   if(severalModels) {
00233     if(nRegions > 1) {
00234       currRegionModel = setOfRegionModels[idxOfRegionModels[index]];
00235     }
00236     currModel = models[currRegionModel->SelectIndex(kinEnergy)];
00237   }
00238   return currModel;
00239 }
00240 
00241 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00242 
00243 inline const G4DataVector* G4EmModelManager::Cuts() const
00244 {
00245   return theCuts;
00246 }
00247 
00248 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00249 
00250 inline const G4DataVector* G4EmModelManager::SubCutoff() const
00251 {
00252   return theSubCuts;
00253 }
00254 
00255 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00256 
00257 inline void G4EmModelManager::SetFluoFlag(G4bool val)
00258 {
00259   fluoFlag = val;
00260 }
00261 
00262 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00263 
00264 inline G4int G4EmModelManager::NumberOfModels() const
00265 {
00266   return nEmModels;
00267 }
00268 
00269 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00270 
00271 inline G4double 
00272 G4EmModelManager::ComputeDEDX(G4VEmModel* model,
00273                               const G4MaterialCutsCouple* couple,
00274                               G4double e,
00275                               G4double cut,
00276                               G4double emin)
00277 {
00278   G4double dedx = 0.0;
00279   if(model && cut > emin) {
00280     dedx = model->ComputeDEDX(couple,particle,e,cut); 
00281     if(emin > 0.0) {dedx -= model->ComputeDEDX(couple,particle,e,emin);} 
00282   }
00283   return dedx;
00284 }
00285 
00286 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00287 
00288 #endif
00289 

Generated on Mon May 27 17:48:09 2013 for Geant4 by  doxygen 1.4.7