Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNABornExcitationModel.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// $Id$
27//
28
30#include "G4SystemOfUnits.hh"
33
34//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
35
36using namespace std;
37
38//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
39
41 const G4String& nam)
42 :G4VEmModel(nam),isInitialised(false)
43{
44 // nistwater = G4NistManager::Instance()->FindOrBuildMaterial("G4_WATER");
45 fpMolWaterDensity = 0;
46
47 verboseLevel= 0;
48 // Verbosity scale:
49 // 0 = nothing
50 // 1 = warning for energy non-conservation
51 // 2 = details of energy budget
52 // 3 = calculation of cross sections, file openings, sampling of atoms
53 // 4 = entering in methods
54
55 if( verboseLevel>0 )
56 {
57 G4cout << "Born excitation model is constructed " << G4endl;
58 }
60}
61
62//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
63
65{
66 // Cross section
67
68 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
69 for (pos = tableData.begin(); pos != tableData.end(); ++pos)
70 {
71 G4DNACrossSectionDataSet* table = pos->second;
72 delete table;
73 }
74
75}
76
77//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
78
80 const G4DataVector& /*cuts*/)
81{
82
83 if (verboseLevel > 3)
84 G4cout << "Calling G4DNABornExcitationModel::Initialise()" << G4endl;
85
86 G4String fileElectron("dna/sigma_excitation_e_born");
87 G4String fileProton("dna/sigma_excitation_p_born");
88
91
92 G4String electron;
93 G4String proton;
94
95 G4double scaleFactor = (1.e-22 / 3.343) * m*m;
96
97 // *** ELECTRON
98
99 electron = electronDef->GetParticleName();
100
101 tableFile[electron] = fileElectron;
102
103 lowEnergyLimit[electron] = 9. * eV;
104 highEnergyLimit[electron] = 1. * MeV;
105
106 // Cross section
107
109 tableE->LoadData(fileElectron);
110
111 tableData[electron] = tableE;
112
113 // *** PROTON
114
115 proton = protonDef->GetParticleName();
116
117 tableFile[proton] = fileProton;
118
119 lowEnergyLimit[proton] = 500. * keV;
120 highEnergyLimit[proton] = 100. * MeV;
121
122 // Cross section
123
125 tableP->LoadData(fileProton);
126
127 tableData[proton] = tableP;
128
129 //
130
131 if (particle==electronDef)
132 {
133 SetLowEnergyLimit(lowEnergyLimit[electron]);
134 SetHighEnergyLimit(highEnergyLimit[electron]);
135 }
136
137 if (particle==protonDef)
138 {
139 SetLowEnergyLimit(lowEnergyLimit[proton]);
140 SetHighEnergyLimit(highEnergyLimit[proton]);
141 }
142
143 if( verboseLevel>0 )
144 {
145 G4cout << "Born excitation model is initialized " << G4endl
146 << "Energy range: "
147 << LowEnergyLimit() / eV << " eV - "
148 << HighEnergyLimit() / keV << " keV for "
149 << particle->GetParticleName()
150 << G4endl;
151 }
152
153 // Initialize water density pointer
155
156 if (isInitialised) { return; }
158 isInitialised = true;
159}
160
161//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
162
164 const G4ParticleDefinition* particleDefinition,
165 G4double ekin,
166 G4double,
167 G4double)
168{
169 if (verboseLevel > 3)
170 G4cout << "Calling CrossSectionPerVolume() of G4DNABornExcitationModel" << G4endl;
171
172 if (
173 particleDefinition != G4Proton::ProtonDefinition()
174 &&
175 particleDefinition != G4Electron::ElectronDefinition()
176 )
177
178 return 0;
179
180 // Calculate total cross section for model
181
182 G4double lowLim = 0;
183 G4double highLim = 0;
184 G4double sigma=0;
185
186 G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
187
188 if(waterDensity!= 0.0)
189 // if (material == nistwater || material->GetBaseMaterial() == nistwater)
190 {
191 const G4String& particleName = particleDefinition->GetParticleName();
192
193 std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
194 pos1 = lowEnergyLimit.find(particleName);
195 if (pos1 != lowEnergyLimit.end())
196 {
197 lowLim = pos1->second;
198 }
199
200 std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
201 pos2 = highEnergyLimit.find(particleName);
202 if (pos2 != highEnergyLimit.end())
203 {
204 highLim = pos2->second;
205 }
206
207 if (ekin >= lowLim && ekin < highLim)
208 {
209 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
210 pos = tableData.find(particleName);
211
212 if (pos != tableData.end())
213 {
214 G4DNACrossSectionDataSet* table = pos->second;
215 if (table != 0)
216 {
217 sigma = table->FindValue(ekin);
218 }
219 }
220 else
221 {
222 G4Exception("G4DNABornExcitationModel::CrossSectionPerVolume","em0002",
223 FatalException,"Model not applicable to particle type.");
224 }
225 }
226
227 if (verboseLevel > 2)
228 {
229 G4cout << "__________________________________" << G4endl;
230 G4cout << "°°° G4DNABornExcitationModel - XS INFO START" << G4endl;
231 G4cout << "°°° Kinetic energy(eV)=" << ekin/eV << " particle : " << particleName << G4endl;
232 G4cout << "°°° Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
233 G4cout << "°°° Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
234 // G4cout << " - Cross section per water molecule (cm^-1)=" << sigma*material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl;
235 G4cout << "°°° G4DNABornExcitationModel - XS INFO END" << G4endl;
236 }
237
238 } // if (waterMaterial)
239
240 return sigma*waterDensity;
241 // return sigma*material->GetAtomicNumDensityVector()[1];
242}
243
244//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
245
246void G4DNABornExcitationModel::SampleSecondaries(std::vector<G4DynamicParticle*>* /*fvect*/,
247 const G4MaterialCutsCouple* /*couple*/,
248 const G4DynamicParticle* aDynamicParticle,
249 G4double,
250 G4double)
251{
252
253 if (verboseLevel > 3)
254 G4cout << "Calling SampleSecondaries() of G4DNABornExcitationModel" << G4endl;
255
256 G4double k = aDynamicParticle->GetKineticEnergy();
257
258 const G4String& particleName = aDynamicParticle->GetDefinition()->GetParticleName();
259
260 G4int level = RandomSelect(k,particleName);
261 G4double excitationEnergy = waterStructure.ExcitationEnergy(level);
262 G4double newEnergy = k - excitationEnergy;
263
264 if (newEnergy > 0)
265 {
269 }
270
271 const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
273 level,
274 theIncomingTrack);
275}
276
277//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
278
279G4int G4DNABornExcitationModel::RandomSelect(G4double k, const G4String& particle)
280{
281 G4int level = 0;
282
283 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
284 pos = tableData.find(particle);
285
286 if (pos != tableData.end())
287 {
288 G4DNACrossSectionDataSet* table = pos->second;
289
290 if (table != 0)
291 {
292 G4double* valuesBuffer = new G4double[table->NumberOfComponents()];
293 const size_t n(table->NumberOfComponents());
294 size_t i(n);
295 G4double value = 0.;
296
297 while (i>0)
298 {
299 i--;
300 valuesBuffer[i] = table->GetComponent(i)->FindValue(k);
301 value += valuesBuffer[i];
302 }
303
304 value *= G4UniformRand();
305
306 i = n;
307
308 while (i > 0)
309 {
310 i--;
311
312 if (valuesBuffer[i] > value)
313 {
314 delete[] valuesBuffer;
315 return i;
316 }
317 value -= valuesBuffer[i];
318 }
319
320 if (valuesBuffer) delete[] valuesBuffer;
321
322 }
323 }
324 else
325 {
326 G4Exception("G4DNABornExcitationModel::RandomSelect","em0002",
327 FatalException,"Model not applicable to particle type.");
328 }
329 return level;
330}
331
332
@ eExcitedMolecule
@ FatalException
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
G4DNABornExcitationModel(const G4ParticleDefinition *p=0, const G4String &nam="DNABornExcitationModel")
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &= *(new G4DataVector()))
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)
static G4DNAChemistryManager * Instance()
void CreateWaterMolecule(ElectronicModification, G4int, const G4Track *)
virtual G4double FindValue(G4double e, G4int componentId=0) const
virtual size_t NumberOfComponents(void) const
virtual const G4VEMDataSet * GetComponent(G4int componentId) const
virtual G4bool LoadData(const G4String &argFileName)
static G4DNAMolecularMaterial * Instance()
const std::vector< double > * GetNumMolPerVolTableFor(const G4Material *) const
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
size_t GetIndex() const
Definition: G4Material.hh:261
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:576
const G4Track * GetCurrentTrack() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
const G4String & GetParticleName() const
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88
virtual G4double FindValue(G4double x, G4int componentId=0) const =0
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:585
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:109
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:529
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:522
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:592
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41