Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4DNABornExcitationModel2.cc
Go to the documentation of this file.
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//
27
29#include "G4SystemOfUnits.hh"
32#include "G4PhysicsTable.hh"
33#include "G4PhysicsVector.hh"
34#include "G4UnitsTable.hh"
35#include <map>
36
37//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
38
39using namespace std;
40
41//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
42
44 const G4String& nam) :
45G4VEmModel(nam), isInitialised(false), fTableData(0)
46{
47 fpMolWaterDensity = 0;
48 fHighEnergy = 0;
49 fLowEnergy = 0;
50 fParticleDefinition = 0;
51
52 verboseLevel = 0;
53 // Verbosity scale:
54 // 0 = nothing
55 // 1 = warning for energy non-conservation
56 // 2 = details of energy budget
57 // 3 = calculation of cross sections, file openings, sampling of atoms
58 // 4 = entering in methods
59
60 if (verboseLevel > 0)
61 {
62 G4cout << "Born excitation model is constructed " << G4endl;
63 }
65 fLastBinCallForFinalXS = 0;
66 fTotalXS = 0;
67 fTableData = 0;
68
69 // Selection of stationary mode
70
71 statCode = false;
72}
73
74//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
75
77{
78 // Cross section
79 if (fTableData)
80 delete fTableData;
81}
82
83//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
84
86 const G4DataVector& /*cuts*/)
87{
88
89 if (verboseLevel > 3)
90 {
91 G4cout << "Calling G4DNABornExcitationModel2::Initialise()" << G4endl;
92 }
93
94 if(fParticleDefinition != 0 && fParticleDefinition != particle)
95 {
96 G4Exception("G4DNABornExcitationModel2::Initialise","em0001",
97 FatalException,"Model already initialized for another particle type.");
98 }
99
100 fParticleDefinition = particle;
101
102 std::ostringstream fullFileName;
103 const char* path = G4FindDataDir("G4LEDATA");
104
105 if(G4String(path) == "")
106 {
107 G4Exception("G4DNABornExcitationModel2::Initialise","G4LEDATA-CHECK",
108 FatalException, "G4LEDATA not defined in environment variables");
109 }
110
111 fullFileName << path;
112
113 if(particle->GetParticleName() == "e-")
114 {
115 fullFileName << "/dna/bornExcitation-e.dat";
116 fLowEnergy = 9*eV;
117 fHighEnergy = 1*MeV;
118 }
119 else if(particle->GetParticleName() == "proton")
120 {
121 fullFileName << "/dna/bornExcitation-p.dat";
122 fLowEnergy = 500. * keV;
123 fHighEnergy = 100. * MeV;
124 }
125
126 SetLowEnergyLimit(fLowEnergy);
127 SetHighEnergyLimit(fHighEnergy);
128
129 //G4double scaleFactor = (1.e-22 / 3.343) * m*m;
130
131 fTableData = new G4PhysicsTable();
132 fTableData->RetrievePhysicsTable(fullFileName.str().c_str(), true);
133 /*
134 for(std::size_t level = 0; level<fTableData->size(); ++level)
135 {
136 //(*fTableData)(level)->ScaleVector(1,scaleFactor);
137 }
138 */
139 std::size_t finalBin_i = 2000;
140 G4double E_min = fLowEnergy;
141 G4double E_max = fHighEnergy;
142 fTotalXS = new G4PhysicsLogVector(E_min, E_max, finalBin_i, true);
143 G4double energy;
144 G4double finalXS;
145
146 for(std::size_t energy_i = 0; energy_i < finalBin_i; ++energy_i)
147 {
148 energy = fTotalXS->Energy(energy_i);
149 finalXS = 0;
150
151 for(std::size_t level = 0; level<fTableData->size(); ++level)
152 {
153 finalXS += (*fTableData)(level)->Value(energy);
154 }
155 fTotalXS->PutValue(energy_i, finalXS);
156 //G4cout << "energy = " << energy << " " << fTotalXS->Value(energy)
157 // << " " << energy_i << " " << finalXS << G4endl;
158 }
159
160 // for(energy = LowEnergyLimit() ; energy < HighEnergyLimit() ; energy += 1*pow(10,log10(energy)))
161 // {
162 // G4cout << "energy = " << energy << " " << fTotalXS->Value(energy) << G4endl;
163 // }
164
165 if( verboseLevel>0 )
166 {
167 G4cout << "Born excitation model is initialized " << G4endl
168 << "Energy range: "
169 << LowEnergyLimit() / eV << " eV - "
170 << HighEnergyLimit() / keV << " keV for "
171 << particle->GetParticleName()
172 << G4endl;
173 }
174
175 // Initialize water density pointer
177
178 if (isInitialised)
179 { return;}
181 isInitialised = true;
182}
183
184//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
185
187 const G4ParticleDefinition* particleDefinition,
188 G4double ekin,
189 G4double,
190 G4double)
191{
192 if (verboseLevel > 3)
193 {
194 G4cout << "Calling CrossSectionPerVolume() of G4DNABornExcitationModel2"
195 << G4endl;
196 }
197
198 if(particleDefinition != fParticleDefinition) return 0;
199
200 // Calculate total cross section for model
201
202 G4double sigma=0;
203
204 G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
205
206 if (ekin >= fLowEnergy && ekin <= fHighEnergy)
207 {
208 sigma = fTotalXS->Value(ekin, fLastBinCallForFinalXS);
209
210 // for(std::size_t i = 0; i < 5; ++i)
211 // sigma += (*fTableData)[i]->Value(ekin);
212
213 if(sigma == 0)
214 {
215 G4cerr << "PROBLEM SIGMA = 0 at " << G4BestUnit(ekin, "Energy")<< G4endl;
216 }
217 }
218
219 if (verboseLevel > 2)
220 {
221 G4cout << "__________________________________" << G4endl;
222 G4cout << "G4DNABornExcitationModel2 - XS INFO START" << G4endl;
223 G4cout << "Kinetic energy(eV)=" << ekin/eV << " particle : " << particleDefinition->GetParticleName() << G4endl;
224 G4cout << "Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
225 G4cout << "Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
226 G4cout << "G4DNABornExcitationModel2 - XS INFO END" << G4endl;
227 }
228
229 return sigma*waterDensity;
230}
231
232//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
233
234void G4DNABornExcitationModel2::SampleSecondaries(std::vector<G4DynamicParticle*>* /*fvect*/,
235 const G4MaterialCutsCouple* /*couple*/,
236 const G4DynamicParticle* aDynamicParticle,
237 G4double,
238 G4double)
239{
240
241 if (verboseLevel > 3)
242 {
243 G4cout << "Calling SampleSecondaries() of G4DNABornExcitationModel2"
244 << G4endl;
245 }
246
247 G4double k = aDynamicParticle->GetKineticEnergy();
248
249 G4int level = RandomSelect(k);
250 G4double excitationEnergy = waterStructure.ExcitationEnergy(level);
251 G4double newEnergy = k - excitationEnergy;
252
253 if (newEnergy > 0)
254 {
256
257 if (!statCode) fParticleChangeForGamma->SetProposedKineticEnergy(newEnergy);
259
261 }
262
263 const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
265 level,
266 theIncomingTrack);
267}
268
269//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
270
272 G4int level,
273 const G4ParticleDefinition* particle,
274 G4double kineticEnergy)
275{
276 if (fParticleDefinition != particle)
277 {
278 G4Exception("G4DNABornExcitationModel2::GetPartialCrossSection",
279 "bornParticleType",
281 "Model initialized for another particle type.");
282 }
283
284 return (*fTableData)(level)->Value(kineticEnergy);
285}
286
287//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
288
289G4int G4DNABornExcitationModel2::RandomSelect(G4double k)
290{
291 const std::size_t n(fTableData->size());
292 std::size_t i(n);
293
294 G4double value = fTotalXS->Value(k, fLastBinCallForFinalXS);
295
296 value *= G4UniformRand();
297 i = n;
298
299 G4double partialXS;
300
301 while (i > 0)
302 {
303 i--;
304
305 partialXS = (*fTableData)(i)->Value(k);
306 if (partialXS > value)
307 {
308 return (G4int)i;
309 }
310 value -= partialXS;
311 }
312
313 return 0;
314}
315
@ eExcitedMolecule
const char * G4FindDataDir(const char *)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
G4ParticleChangeForGamma * fParticleChangeForGamma
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
G4DNABornExcitationModel2(const G4ParticleDefinition *p=0, const G4String &nam="DNABornExcitationModel")
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &= *(new G4DataVector()))
virtual G4double GetPartialCrossSection(const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
static G4DNAChemistryManager * Instance()
void CreateWaterMolecule(ElectronicModification, G4int, const G4Track *)
const std::vector< G4double > * GetNumMolPerVolTableFor(const G4Material *) const
Retrieve a table of molecular densities (number of molecules per unit volume) in the G4 unit system f...
static G4DNAMolecularMaterial * Instance()
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
size_t GetIndex() const
Definition: G4Material.hh:255
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:691
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
const G4String & GetParticleName() const
G4bool RetrievePhysicsTable(const G4String &filename, G4bool ascii=false, G4bool spline=false)
void PutValue(const std::size_t index, const G4double value)
G4double Energy(const std::size_t index) const
G4double Value(const G4double energy, std::size_t &lastidx) const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:746
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:124
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
virtual G4double Value(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:349
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:753
const G4Track * GetCurrentTrack() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)