|
||||
File indexing completed on 2025-01-18 09:58:54
0001 // 0002 // ******************************************************************** 0003 // * License and Disclaimer * 0004 // * * 0005 // * The Geant4 software is copyright of the Copyright Holders of * 0006 // * the Geant4 Collaboration. It is provided under the terms and * 0007 // * conditions of the Geant4 Software License, included in the file * 0008 // * LICENSE and available at http://cern.ch/geant4/license . These * 0009 // * include a list of copyright holders. * 0010 // * * 0011 // * Neither the authors of this software system, nor their employing * 0012 // * institutes,nor the agencies providing financial support for this * 0013 // * work make any representation or warranty, express or implied, * 0014 // * regarding this software system or assume any liability for its * 0015 // * use. Please see the license in the file LICENSE and URL above * 0016 // * for the full disclaimer and the limitation of liability. * 0017 // * * 0018 // * This code implementation is the result of the scientific and * 0019 // * technical work of the GEANT4 collaboration. * 0020 // * By using, copying, modifying or distributing the software (or * 0021 // * any work based on the software) you agree to acknowledge its * 0022 // * use in resulting scientific publications, and indicate your * 0023 // * acceptance of all terms of the Geant4 Software license. * 0024 // ******************************************************************** 0025 // 0026 // 0027 // Author: Luciano Pandola 0028 // 0029 // History: 0030 // ----------- 0031 // 20 Oct 2010 L. Pandola 1st implementation. 0032 // 02 May 2011 L. Pandola Remove dependency on CLHEP::HepMatrix 0033 // 24 May 2011 L. Pandola Renamed (make default Penelope) 0034 // 03 Oct 2013 L. Pandola Migration to MT 0035 // 07 Oct 2013 L. Pandola Add verbosity and ismaster flag for the 0036 // master-only methods 0037 // 30 Oct 2013 L. Pandola Add a G4Cache of temporary vectors as 0038 // private member. 0039 // 0040 // ------------------------------------------------------------------- 0041 // 0042 // Class description: 0043 // Helper class for the calculation of final state (energy and angular 0044 // distribution) for bremsstrahlung, Penelope Model, version 2008 0045 // ------------------------------------------------------------------- 0046 0047 #ifndef G4PENELOPEBREMSSTRAHLUNGFS_HH 0048 #define G4PENELOPEBREMSSTRAHLUNGFS_HH 1 0049 0050 #include "globals.hh" 0051 #include "G4DataVector.hh" 0052 #include "G4Cache.hh" 0053 #include <map> 0054 0055 class G4PhysicsFreeVector; 0056 class G4PhysicsLogVector; 0057 class G4Material; 0058 class G4PhysicsTable; 0059 0060 class G4PenelopeBremsstrahlungFS 0061 { 0062 public: 0063 //! Only master models are supposed to create instances 0064 explicit G4PenelopeBremsstrahlungFS(G4int verbosity=0); 0065 ~G4PenelopeBremsstrahlungFS(); 0066 0067 //! 0068 //! Master and workers (do not touch tables) 0069 //! All of them are const 0070 //! 0071 G4double GetEffectiveZSquared(const G4Material* mat) const; 0072 size_t GetNBinsX() const {return fNBinsX;}; 0073 G4double GetMomentumIntegral(G4double* y, 0074 G4double up,G4int momOrder) const; 0075 const G4PhysicsTable* GetScaledXSTable(const G4Material*, 0076 const G4double cut) const; 0077 G4double SampleGammaEnergy(G4double energy, 0078 const G4Material*, const G4double cut) const; 0079 0080 //! Reserved for the master model: they build and handle tables 0081 void ClearTables(G4bool isMaster=true); 0082 void BuildScaledXSTable(const G4Material* material,G4double cut, 0083 G4bool isMaster); 0084 0085 void SetVerbosity(G4int ver){fVerbosity=ver;}; 0086 G4int GetVerbosity(){return fVerbosity;}; 0087 0088 G4PenelopeBremsstrahlungFS & operator=(const G4PenelopeBremsstrahlungFS &right) = delete; 0089 G4PenelopeBremsstrahlungFS(const G4PenelopeBremsstrahlungFS&) = delete; 0090 0091 private: 0092 void ReadDataFile(G4int Z); 0093 //Tables for energy sampling 0094 void InitializeEnergySampling(const G4Material*,G4double cut); 0095 0096 //Differential cross section tables 0097 //Table contains G4PhysicsVectors of log(XS(E,x)) vs. log(E) 0098 //for a grid of 32 values in x (= reduced photon energy) 0099 std::map< std::pair<const G4Material*,G4double> , 0100 G4PhysicsTable*> *fReducedXSTable; 0101 0102 std::map<const G4Material*,G4double> *fEffectiveZSq; 0103 0104 //Map of element data vs. Z. 0105 //This is conceptually an array [57][33] with 57 energy 0106 //points and 32 points in x. The 33-th column gives the total XS vs. E. 0107 //It is implemented as a one-dimensional array of dimension 0108 //57*33=1881 elements. data[e][x] --> theElementData[e*(nBinsX+1)+x] 0109 std::map<G4int,G4DataVector*> *fElementData; 0110 0111 //Table contains G4PhysicsVectors of integral_XS(E,x) vs. x for a grid of 0112 //57 values in energy 0113 std::map< std::pair<const G4Material*,G4double> , 0114 G4PhysicsTable*> *fSamplingTable; 0115 std::map< std::pair<const G4Material*,G4double> , 0116 G4PhysicsFreeVector* > *fPBcut; 0117 0118 //temporary vector. Used as member variable to avoid to book/release the 0119 //memory on the fly. This vector is over-written at every call of 0120 //SampleGammaEnergy(). It is thread-local (each thread has its own) 0121 //and managed by G4Cache 0122 G4Cache<G4PhysicsFreeVector*> fCache; 0123 0124 //Element data table 0125 static const size_t fNBinsE = 57; 0126 static const size_t fNBinsX = 32; 0127 //x and E grids used in the data tables 0128 G4double theXGrid[fNBinsX]; 0129 G4double theEGrid[fNBinsE]; 0130 0131 G4int fVerbosity; 0132 0133 }; 0134 0135 #endif 0136
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |