|
||||
Warning, file /include/Geant4/G4SBBremTable.hh was not indexed or was modified since last indexation (in which case cross-reference links may be missing, inaccurate or erroneous).
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 // ---------------------------------------------------------------------------- 0028 // 0029 // GEANT4 Class header file 0030 // 0031 // 0032 // File name: G4SBBremTable 0033 // 0034 // Author: Mihaly Novak 0035 // 0036 // Creation date: 15.07.2018 0037 // 0038 // Modifications: 0039 // 0040 // Class description: 0041 // 0042 // Utility class to handle sampling tables for the Seltzer-Berger scalled brems- 0043 // strahlung differential cross sections. It makes possible fast (significantly 0044 // faster than the rejection) sampling of the emitted photon energy in case of 0045 // interactions. An object from this class is supposed to be a member of the 0046 // Seltzer-Berger model for e-/e+ bremsstrahlung photon emission model. Note, 0047 // that one object from this class can handle both e- and e+ cases (containes 0048 // e+ correction in the SampleEnergyTransfer method only). 0049 // 0050 // ---------------------------------------------------------------------------- 0051 0052 #ifndef G4SBBremTable_h 0053 #define G4SBBremTable_h 1 0054 0055 #include "globals.hh" 0056 #include "G4String.hh" 0057 0058 #include <vector> 0059 0060 // forward declar 0061 class G4MaterialCutsCouple; 0062 0063 class G4SBBremTable { 0064 0065 public: 0066 // CTR/DTR 0067 G4SBBremTable(); 0068 0069 ~G4SBBremTable(); 0070 0071 // loads and init sampling tables: lowe/highe are the low/high energy usage 0072 // limits of the corresponding Seltzerberger-model. 0073 void Initialize(const G4double lowe, const G4double highe); 0074 0075 // clean away all sampling tables and makes ready for re-initialisation 0076 void ClearSamplingTables(); 0077 0078 // run-time method to sample energy transferred to the emitted photon 0079 double SampleEnergyTransfer(const G4double eekin, const G4double leekin, 0080 const G4double gcut , const G4double dielSupConst, 0081 const G4int izet , const G4int matCutIndx, 0082 const bool iselectron); 0083 0084 // used only for development: print out table related information 0085 // void Dump(); 0086 0087 private: 0088 0089 void BuildSamplingTables(); 0090 0091 void InitSamplingTables(); 0092 0093 void LoadSTGrid(); 0094 0095 void LoadSamplingTables(G4int iz); 0096 0097 void ReadCompressedFile(const G4String &fname, std::istringstream &iss); 0098 0099 private: 0100 0101 // Sampling-Table point: describes one [E_i],[kappa_j] point 0102 struct STPoint { 0103 G4double fCum; // value of the cumulative function 0104 G4double fParA; // rational function approximation based interp. parameter 0105 G4double fParB; // rational function approximation based interp. parameter 0106 }; 0107 0108 // Sampling-Table: describes one [E_j] e- energy point i.e. one Table 0109 struct STable { 0110 // cumulative values for the kappa-cuts: kappa_cut_i=E_gamma_cut_i/E_el_j 0111 std::vector<G4double> fCumCutValues; 0112 // as many STPoint-s as kappa values 0113 std::vector<STPoint> fSTable; 0114 }; 0115 0116 // Sampling-Tables for a given Z: 0117 // describes all tables (i.e. for all e- energies) for a given element (Z) 0118 struct SamplingTablePerZ { 0119 SamplingTablePerZ() : fNumGammaCuts(0), fMinElEnergyIndx(-1), fMaxElEnergyIndx(-1) {} 0120 size_t fNumGammaCuts; // number of gamma-cut for this 0121 G4int fMinElEnergyIndx; // max(i) such E_i <= E for all E 0122 G4int fMaxElEnergyIndx; // min(i) such E_i >= E for all E 0123 std::vector<STable*> fTablesPerEnergy; // as many table as e-ekin grid point 0124 //the different gamma-cut values that are defined for this element(Z) and ln 0125 std::vector<G4double> fGammaECuts; 0126 std::vector<G4double> fLogGammaECuts; 0127 // the couple index element stores the corresponding (sorted) gamma-cut index 0128 std::vector<size_t> fMatCutIndxToGamCutIndx; 0129 // temporary vector to store some indecis during initialisation 0130 std::vector< std::vector<size_t> > fGamCutIndxToMatCutIndx; 0131 }; 0132 0133 // simple linear search: most of the time faster than anything in our case 0134 G4int LinSearch(const std::vector<STPoint>& vect, 0135 const G4int size, 0136 const G4double val); 0137 0138 private: 0139 0140 // pre-prepared sampling tables are available: 0141 G4int fMaxZet; // max Z number 0142 G4int fNumElEnergy; // # e- kine (E_k) per Z 0143 G4int fNumKappa; // # red. photon eners per E_k 0144 0145 // min/max electron kinetic energy usage limits 0146 G4double fUsedLowEenergy; 0147 G4double fUsedHighEenergy; 0148 G4double fLogMinElEnergy; 0149 G4double fILDeltaElEnergy; 0150 0151 // e- kinetic energy and reduced photon energy grids and tehir logarithms 0152 std::vector<G4double> fElEnergyVect; 0153 std::vector<G4double> fLElEnergyVect; 0154 std::vector<G4double> fKappaVect; 0155 std::vector<G4double> fLKappaVect; 0156 0157 // container to store samplingtables per Z (size is fMaxZet+1) 0158 std::vector<SamplingTablePerZ*> fSBSamplingTables; 0159 0160 }; 0161 0162 #endif
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |