Geant4 Cross Reference |
1 // 2 // ******************************************************************** 3 // * License and Disclaimer * 4 // * * 5 // * The Geant4 software is copyright of the Copyright Holders of * 6 // * the Geant4 Collaboration. It is provided under the terms and * 7 // * conditions of the Geant4 Software License, included in the file * 8 // * LICENSE and available at http://cern.ch/geant4/license . These * 9 // * include a list of copyright holders. * 10 // * * 11 // * Neither the authors of this software system, nor their employing * 12 // * institutes,nor the agencies providing financial support for this * 13 // * work make any representation or warranty, express or implied, * 14 // * regarding this software system or assume any liability for its * 15 // * use. Please see the license in the file LICENSE and URL above * 16 // * for the full disclaimer and the limitation of liability. * 17 // * * 18 // * This code implementation is the result of the scientific and * 19 // * technical work of the GEANT4 collaboration. * 20 // * By using, copying, modifying or distributing the software (or * 21 // * any work based on the software) you agree to acknowledge its * 22 // * use in resulting scientific publications, and indicate your * 23 // * acceptance of all terms of the Geant4 Software license. * 24 // ******************************************************************** 25 // 26 // Author: Sebastien Incerti 27 // 22 January 2012 28 // on base of G4LivermoreGammaConversionModel (original version) 29 // and G4LivermoreRayleighModel (MT version) 30 // 31 // Modifications: Zhuxin Li@CENBG 32 // 11 March 2020 33 // derives from G4PairProductionRelModel 34 // ------------------------------------------------------------------- 35 36 #include "G4LivermoreGammaConversionModel.hh" 37 38 #include "G4AutoLock.hh" 39 #include "G4Electron.hh" 40 #include "G4EmParameters.hh" 41 #include "G4Exp.hh" 42 #include "G4ParticleChangeForGamma.hh" 43 #include "G4PhysicalConstants.hh" 44 #include "G4PhysicsFreeVector.hh" 45 #include "G4SystemOfUnits.hh" 46 47 namespace 48 { 49 G4Mutex LivermoreGammaConversionModelMutex = G4MUTEX_INITIALIZER; 50 } 51 52 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 53 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 54 55 G4PhysicsFreeVector* G4LivermoreGammaConversionModel::data[] = {nullptr}; 56 G4String G4LivermoreGammaConversionModel::gDataDirectory = ""; 57 58 G4LivermoreGammaConversionModel::G4LivermoreGammaConversionModel(const G4ParticleDefinition* p, 59 const G4String& nam) 60 : G4PairProductionRelModel(p, nam) 61 { 62 fParticleChange = nullptr; 63 lowEnergyLimit = 2. * CLHEP::electron_mass_c2; 64 verboseLevel = 0; 65 // Verbosity scale for debugging purposes: 66 // 0 = nothing 67 // 1 = calculation of cross sections, file openings... 68 // 2 = entering in methods 69 if (verboseLevel > 0) { 70 G4cout << "G4LivermoreGammaConversionModel is constructed " << G4endl; 71 } 72 } 73 74 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 75 76 G4LivermoreGammaConversionModel::~G4LivermoreGammaConversionModel() 77 { 78 if (IsMaster()) { 79 for (G4int i = 0; i <= maxZ; ++i) { 80 if (data[i]) { 81 delete data[i]; 82 data[i] = nullptr; 83 } 84 } 85 } 86 } 87 88 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 89 90 void G4LivermoreGammaConversionModel::Initialise(const G4ParticleDefinition* particle, 91 const G4DataVector& cuts) 92 { 93 G4PairProductionRelModel::Initialise(particle, cuts); 94 if (verboseLevel > 1) { 95 G4cout << "Calling Initialise() of G4LivermoreGammaConversionModel." << G4endl 96 << "Energy range: " << LowEnergyLimit() / MeV << " MeV - " << HighEnergyLimit() / GeV 97 << " GeV isMater: " << IsMaster() << G4endl; 98 } 99 100 if (IsMaster()) { 101 // Initialise element selector 102 InitialiseElementSelectors(particle, cuts); 103 104 // Access to elements 105 const G4ElementTable* elemTable = G4Element::GetElementTable(); 106 std::size_t numElems = (*elemTable).size(); 107 for (std::size_t ie = 0; ie < numElems; ++ie) { 108 const G4Element* elem = (*elemTable)[ie]; 109 const G4int Z = std::min(maxZ, elem->GetZasInt()); 110 if (data[Z] == nullptr) { 111 ReadData(Z); 112 } 113 } 114 } 115 if (isInitialised) { 116 return; 117 } 118 fParticleChange = GetParticleChangeForGamma(); 119 isInitialised = true; 120 } 121 122 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 123 124 const G4String& G4LivermoreGammaConversionModel::FindDirectoryPath() 125 { 126 // no check in this method - environment variable is check by utility 127 if (gDataDirectory.empty()) { 128 auto param = G4EmParameters::Instance(); 129 std::ostringstream ost; 130 if (param->LivermoreDataDir() == "livermore") { 131 ost << param->GetDirLEDATA() << "/livermore/pair/"; 132 useSpline = true; 133 } 134 else { 135 ost << param->GetDirLEDATA() << "/epics2017/pair/"; 136 } 137 gDataDirectory = ost.str(); 138 } 139 return gDataDirectory; 140 } 141 142 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 143 144 void G4LivermoreGammaConversionModel::ReadData(const G4int Z) 145 { 146 if (verboseLevel > 1) { 147 G4cout << "Calling ReadData() of G4LivermoreGammaConversionModel" << G4endl; 148 } 149 150 if (data[Z] != nullptr) { 151 return; 152 } 153 154 std::ostringstream ost; 155 ost << FindDirectoryPath() << "pp-cs-" << Z << ".dat"; 156 157 data[Z] = new G4PhysicsFreeVector(useSpline); 158 159 std::ifstream fin(ost.str().c_str()); 160 161 if (!fin.is_open()) { 162 G4ExceptionDescription ed; 163 ed << "G4LivermoreGammaConversionModel data file <" << ost.str().c_str() << "> is not opened!" 164 << G4endl; 165 G4Exception("G4LivermoreGammaConversionModel::ReadData()", "em0003", FatalException, ed, 166 "G4LEDATA version should be G4EMLOW8.0 or later."); 167 return; 168 } 169 else { 170 if (verboseLevel > 1) { 171 G4cout << "File " << ost.str() << " is opened by G4LivermoreGammaConversionModel" << G4endl; 172 } 173 174 data[Z]->Retrieve(fin, true); 175 } 176 // Activation of spline interpolation 177 if (useSpline) data[Z]->FillSecondDerivatives(); 178 } 179 180 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 181 182 G4double 183 G4LivermoreGammaConversionModel::ComputeCrossSectionPerAtom(const G4ParticleDefinition* particle, 184 G4double GammaEnergy, G4double Z, 185 G4double, G4double, G4double) 186 { 187 if (verboseLevel > 1) { 188 G4cout << "G4LivermoreGammaConversionModel::ComputeCrossSectionPerAtom() Z= " << Z << G4endl; 189 } 190 191 if (GammaEnergy < lowEnergyLimit) { 192 return 0.0; 193 } 194 195 G4double xs = 0.0; 196 197 G4int intZ = std::max(1, std::min(G4lrint(Z), maxZ)); 198 199 G4PhysicsFreeVector* pv = data[intZ]; 200 201 // if element was not initialised 202 // do initialisation safely for MT mode 203 if (pv == nullptr) { 204 InitialiseForElement(particle, intZ); 205 pv = data[intZ]; 206 if (pv == nullptr) { 207 return xs; 208 } 209 } 210 // x-section is taken from the table 211 xs = pv->Value(GammaEnergy); 212 213 if (verboseLevel > 0) { 214 G4cout << "*** Gamma conversion xs for Z=" << Z << " at energy E(MeV)=" << GammaEnergy / MeV 215 << " cs=" << xs / millibarn << " mb" << G4endl; 216 } 217 return xs; 218 } 219 220 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 221 222 void G4LivermoreGammaConversionModel::InitialiseForElement(const G4ParticleDefinition*, G4int Z) 223 { 224 G4AutoLock l(&LivermoreGammaConversionModelMutex); 225 if (data[Z] == nullptr) { 226 ReadData(Z); 227 } 228 l.unlock(); 229 } 230 231 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 232