Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/lowenergy/src/G4LivermoreRayleighModel.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 //         31 March 2012
 28 //         on base of G4LivermoreRayleighModel
 29 //
 30 
 31 #include "G4LivermoreRayleighModel.hh"
 32 
 33 #include "G4AutoLock.hh"
 34 #include "G4EmParameters.hh"
 35 #include "G4RayleighAngularGenerator.hh"
 36 #include "G4SystemOfUnits.hh"
 37 
 38 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 39 
 40 namespace
 41 {
 42 G4Mutex LivermoreRayleighModelMutex = G4MUTEX_INITIALIZER;
 43 }
 44 
 45 G4PhysicsFreeVector* G4LivermoreRayleighModel::dataCS[] = {nullptr};
 46 G4String G4LivermoreRayleighModel::gDataDirectory = "";
 47 
 48 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 49 
 50 G4LivermoreRayleighModel::G4LivermoreRayleighModel() : G4VEmModel("LivermoreRayleigh")
 51 {
 52   fParticleChange = nullptr;
 53   lowEnergyLimit = 10 * CLHEP::eV;
 54 
 55   SetAngularDistribution(new G4RayleighAngularGenerator());
 56 
 57   verboseLevel = 0;
 58   // Verbosity scale for debugging purposes:
 59   // 0 = nothing
 60   // 1 = calculation of cross sections, file openings...
 61   // 2 = entering in methods
 62 
 63   if (verboseLevel > 0) {
 64     G4cout << "G4LivermoreRayleighModel is constructed " << G4endl;
 65   }
 66 }
 67 
 68 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 69 
 70 G4LivermoreRayleighModel::~G4LivermoreRayleighModel()
 71 {
 72   if (IsMaster()) {
 73     for (G4int i = 0; i <= maxZ; ++i) {
 74       if (nullptr != dataCS[i]) {
 75         delete dataCS[i];
 76         dataCS[i] = nullptr;
 77       }
 78     }
 79   }
 80 }
 81 
 82 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 83 
 84 void G4LivermoreRayleighModel::Initialise(const G4ParticleDefinition* particle,
 85                                           const G4DataVector& cuts)
 86 {
 87   if (verboseLevel > 1) {
 88     G4cout << "Calling Initialise() of G4LivermoreRayleighModel." << G4endl
 89            << "Energy range: " << LowEnergyLimit() / eV << " eV - " << HighEnergyLimit() / GeV
 90            << " GeV" << G4endl;
 91   }
 92 
 93   if (IsMaster()) {
 94     // Initialise element selector
 95     InitialiseElementSelectors(particle, cuts);
 96 
 97     // Access to elements
 98     const G4ElementTable* elemTable = G4Element::GetElementTable();
 99     std::size_t numElems = (*elemTable).size();
100     for (std::size_t ie = 0; ie < numElems; ++ie) {
101       const G4Element* elem = (*elemTable)[ie];
102       const G4int Z = std::min(maxZ, elem->GetZasInt());
103       if (dataCS[Z] == nullptr) {
104         ReadData(Z);
105       }
106     }
107   }
108   if (isInitialised) {
109     return;
110   }
111   fParticleChange = GetParticleChangeForGamma();
112   isInitialised = true;
113 }
114 
115 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
116 
117 void G4LivermoreRayleighModel::InitialiseLocal(const G4ParticleDefinition*, G4VEmModel* masterModel)
118 {
119   SetElementSelectors(masterModel->GetElementSelectors());
120 }
121 
122 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
123 
124 const G4String& G4LivermoreRayleighModel::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/rayl/";
132     }
133     else {
134       ost << param->GetDirLEDATA() << "/epics2017/rayl/";
135     }
136     gDataDirectory = ost.str();
137   }
138   return gDataDirectory;
139 }
140 
141 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
142 
143 void G4LivermoreRayleighModel::ReadData(const G4int ZZ)
144 {
145   if (verboseLevel > 1) {
146     G4cout << "Calling ReadData() of G4LivermoreRayleighModel for Z=" << ZZ << G4endl;
147   }
148   const G4int Z = std::min(ZZ, maxZ);
149 
150   if (nullptr != dataCS[Z]) {
151     return;
152   }
153 
154   dataCS[Z] = new G4PhysicsFreeVector();
155 
156   std::ostringstream ostCS;
157   ostCS << FindDirectoryPath() << "re-cs-" << Z << ".dat";
158 
159   std::ifstream finCS(ostCS.str().c_str());
160 
161   if (!finCS.is_open()) {
162     G4ExceptionDescription ed;
163     ed << "G4LivermoreRayleighModel data file <" << ostCS.str().c_str() << "> is not opened!"
164        << G4endl;
165     G4Exception("G4LivermoreRayleighModel::ReadData()", "em0003", FatalException, ed,
166                 "G4LEDATA version should be G4EMLOW8.0 or later.");
167     return;
168   }
169   else {
170     if (verboseLevel > 3) {
171       G4cout << "File " << ostCS.str() << " is opened by G4LivermoreRayleighModel" << G4endl;
172     }
173     dataCS[Z]->Retrieve(finCS, true);
174   }
175 }
176 
177 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
178 
179 G4double G4LivermoreRayleighModel::ComputeCrossSectionPerAtom(const G4ParticleDefinition*,
180                                                               G4double GammaEnergy, G4double Z,
181                                                               G4double, G4double, G4double)
182 {
183   if (verboseLevel > 1) {
184     G4cout << "G4LivermoreRayleighModel::ComputeCrossSectionPerAtom()" << G4endl;
185   }
186 
187   if (GammaEnergy < lowEnergyLimit) {
188     return 0.0;
189   }
190 
191   G4double xs = 0.0;
192   G4int intZ = G4lrint(Z);
193   if (intZ < 1 || intZ > maxZ) {
194     return xs;
195   }
196 
197   G4PhysicsFreeVector* pv = dataCS[intZ];
198 
199   // if element was not initialised
200   // do initialisation safely for MT mode
201   if (nullptr == pv) {
202     InitialiseForElement(nullptr, intZ);
203     pv = dataCS[intZ];
204     if (nullptr == pv) {
205       return xs;
206     }
207   }
208 
209   auto n = G4int(pv->GetVectorLength() - 1);
210   G4double e = GammaEnergy / MeV;
211   if (e >= pv->Energy(n)) {
212     xs = (*pv)[n] / (e * e);
213   }
214   else if (e >= pv->Energy(0)) {
215     xs = pv->Value(e) / (e * e);
216   }
217 
218   if (verboseLevel > 1) {
219     G4cout << "****** DEBUG: tcs value for Z=" << Z << " at energy (MeV)=" << e << G4endl;
220     G4cout << "  cs (Geant4 internal unit)=" << xs << G4endl;
221     G4cout << "    -> first E*E*cs value in CS data file (iu) =" << (*pv)[0] << G4endl;
222     G4cout << "    -> last  E*E*cs value in CS data file (iu) =" << (*pv)[n] << G4endl;
223     G4cout << "*********************************************************" << G4endl;
224   }
225   return xs;
226 }
227 
228 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
229 
230 void G4LivermoreRayleighModel::SampleSecondaries(std::vector<G4DynamicParticle*>*,
231                                                  const G4MaterialCutsCouple* couple,
232                                                  const G4DynamicParticle* aDynamicGamma, G4double,
233                                                  G4double)
234 {
235   if (verboseLevel > 1) {
236     G4cout << "Calling SampleSecondaries() of G4LivermoreRayleighModel" << G4endl;
237   }
238   G4double photonEnergy0 = aDynamicGamma->GetKineticEnergy();
239 
240   // Select randomly one element in the current material
241   const G4ParticleDefinition* particle = aDynamicGamma->GetDefinition();
242   const G4Element* elm = SelectRandomAtom(couple, particle, photonEnergy0);
243   G4int Z = elm->GetZasInt();
244 
245   // Sample the angle of the scattered photon
246   G4ThreeVector photonDirection = GetAngularDistribution()->SampleDirection(
247     aDynamicGamma, photonEnergy0, Z, couple->GetMaterial());
248   fParticleChange->ProposeMomentumDirection(photonDirection);
249 }
250 
251 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
252 
253 void G4LivermoreRayleighModel::InitialiseForElement(const G4ParticleDefinition*, G4int Z)
254 {
255   if (nullptr != dataCS[Z]) {
256     return;
257   }
258   G4AutoLock l(&LivermoreRayleighModelMutex);
259   if (nullptr == dataCS[Z]) {
260     ReadData(Z);
261   }
262   l.unlock();
263 }
264 
265 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
266