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
G4DNAPTBElasticModel.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// Authors: S. Meylan and C. Villagrasa (IRSN, France)
27// Models come from
28// M. Bug et al, Rad. Phys and Chem. 130, 459-479 (2017)
29//
30
34#include "G4SystemOfUnits.hh"
36#include "G4Proton.hh"
37
39 const G4String& nam)
40 : G4VDNAModel(nam, applyToMaterial)
41{
42 fKillBelowEnergy = 10*eV; // will be override by the limits defined for each material
43
44 verboseLevel= 0;
45 // Verbosity scale:
46 // 0 = nothing
47 // 1 = warning for energy non-conservation
48 // 2 = details of energy budget
49 // 3 = calculation of cross sections, file openings, sampling of atoms
50 // 4 = entering in methods
51
52 if( verboseLevel>0 )
53 {
54 G4cout << "PTB Elastic model is constructed " << G4endl;
55 }
56}
57
58//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
59
61{
62
63}
64
65//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
66
68 const G4DataVector& /*cuts*/, G4ParticleChangeForGamma*)
69{
70 if (verboseLevel > 3)
71 G4cout << "Calling G4DNAPTBElasticModel::Initialise()" << G4endl;
72
73 G4double scaleFactor = 1e-16*cm*cm;
74
76
77 //*******************************************************
78 // Cross section data
79 //*******************************************************
80
81 if(particle == electronDef)
82 {
83 G4String particleName = particle->GetParticleName();
84
85 // MPietrzak, adding paths for N2
87 particleName,
88 "dna/sigma_elastic_e-_PTB_N2",
89 "dna/sigmadiff_cumulated_elastic_e-_PTB_N2",
90 scaleFactor);
91 SetLowELimit("N2", particleName, 10*eV);
92 SetHighELimit("N2", particleName, 1.02*MeV);
93 // MPietrzak
94
96 particleName,
97 "dna/sigma_elastic_e-_PTB_THF",
98 "dna/sigmadiff_cumulated_elastic_e-_PTB_THF",
99 scaleFactor);
100 SetLowELimit("THF", particleName, 10*eV);
101 SetHighELimit("THF", particleName, 1*keV);
102
104 particleName,
105 "dna/sigma_elastic_e-_PTB_PY",
106 "dna/sigmadiff_cumulated_elastic_e-_PTB_PY",
107 scaleFactor);
108 SetLowELimit("PY", particleName, 10*eV);
109 SetHighELimit("PY", particleName, 1*keV);
110
112 particleName,
113 "dna/sigma_elastic_e-_PTB_PU",
114 "dna/sigmadiff_cumulated_elastic_e-_PTB_PU",
115 scaleFactor);
116 SetLowELimit("PU", particleName, 10*eV);
117 SetHighELimit("PU", particleName, 1*keV);
118
120 particleName,
121 "dna/sigma_elastic_e-_PTB_TMP",
122 "dna/sigmadiff_cumulated_elastic_e-_PTB_TMP",
123 scaleFactor);
124 SetLowELimit("TMP", particleName, 10*eV);
125 SetHighELimit("TMP", particleName, 1*keV);
126
127 AddCrossSectionData("G4_WATER",
128 particleName,
129 "dna/sigma_elastic_e_champion",
130 "dna/sigmadiff_cumulated_elastic_e_champion",
131 scaleFactor);
132 SetLowELimit("G4_WATER", particleName, 10*eV);
133 SetHighELimit("G4_WATER", particleName, 1*keV);
134
135 // DNA materials
136 //
137 AddCrossSectionData("backbone_THF",
138 particleName,
139 "dna/sigma_elastic_e-_PTB_THF",
140 "dna/sigmadiff_cumulated_elastic_e-_PTB_THF",
141 scaleFactor*33./30);
142 SetLowELimit("backbone_THF", particleName, 10*eV);
143 SetHighELimit("backbone_THF", particleName, 1*keV);
144
145 AddCrossSectionData("cytosine_PY",
146 particleName,
147 "dna/sigma_elastic_e-_PTB_PY",
148 "dna/sigmadiff_cumulated_elastic_e-_PTB_PY",
149 scaleFactor*42./30);
150 SetLowELimit("cytosine_PY", particleName, 10*eV);
151 SetHighELimit("cytosine_PY", particleName, 1*keV);
152
153 AddCrossSectionData("thymine_PY",
154 particleName,
155 "dna/sigma_elastic_e-_PTB_PY",
156 "dna/sigmadiff_cumulated_elastic_e-_PTB_PY",
157 scaleFactor*48./30);
158 SetLowELimit("thymine_PY", particleName, 10*eV);
159 SetHighELimit("thymine_PY", particleName, 1*keV);
160
161 AddCrossSectionData("adenine_PU",
162 particleName,
163 "dna/sigma_elastic_e-_PTB_PU",
164 "dna/sigmadiff_cumulated_elastic_e-_PTB_PU",
165 scaleFactor*50./44);
166 SetLowELimit("adenine_PU", particleName, 10*eV);
167 SetHighELimit("adenine_PU", particleName, 1*keV);
168
169 AddCrossSectionData("guanine_PU",
170 particleName,
171 "dna/sigma_elastic_e-_PTB_PU",
172 "dna/sigmadiff_cumulated_elastic_e-_PTB_PU",
173 scaleFactor*56./44);
174 SetLowELimit("guanine_PU", particleName, 10*eV);
175 SetHighELimit("guanine_PU", particleName, 1*keV);
176
177 AddCrossSectionData("backbone_TMP",
178 particleName,
179 "dna/sigma_elastic_e-_PTB_TMP",
180 "dna/sigmadiff_cumulated_elastic_e-_PTB_TMP",
181 scaleFactor*33./50);
182 SetLowELimit("backbone_TMP", particleName, 10*eV);
183 SetHighELimit("backbone_TMP", particleName, 1*keV);
184 }
185
186 //*******************************************************
187 // Load the data
188 //*******************************************************
189
191
192 //*******************************************************
193 // Verbose output
194 //*******************************************************
195
196 if (verboseLevel > 2)
197 G4cout << "Loaded cross section files for PTB Elastic model" << G4endl;
198
199 if( verboseLevel>0 )
200 {
201 G4cout << "PTB Elastic model is initialized " << G4endl;
202 }
203}
204
205//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
206
207void G4DNAPTBElasticModel::ReadDiffCSFile(const G4String& materialName,
208 const G4String& particleName,
209 const G4String& file,
210 const G4double)
211{
212 // Method to read and save the information contained within the differential cross section files.
213 // This method is not yet standard.
214
215 // get the path of the G4LEDATA data folder
216 const char* path = G4FindDataDir("G4LEDATA");
217 // if it is not found then quit and print error message
218 if(!path)
219 {
220 G4Exception("G4DNAPTBElasticModel::ReadAllDiffCSFiles","em0006",
221 FatalException,"G4LEDATA environment variable not set.");
222 return;
223 }
224
225 // build the fullFileName path of the data file
226 std::ostringstream fullFileName;
227 fullFileName << path <<"/"<< file<<".dat";
228
229 // open the data file
230 std::ifstream diffCrossSection (fullFileName.str().c_str());
231 // error if file is not there
232 std::stringstream endPath;
233 if (!diffCrossSection)
234 {
235 endPath << "Missing data file: "<<file;
236 G4Exception("G4DNAPTBElasticModel::Initialise","em0003",
237 FatalException, endPath.str().c_str());
238 }
239
240 tValuesVec[materialName][particleName].push_back(0.);
241
242 G4String line;
243
244 // read the file line by line until we reach the end of file point
245 while(std::getline(diffCrossSection, line))
246 {
247 // check if the line is comment or empty
248 //
249 std::istringstream testIss(line);
250 G4String test;
251 testIss >> test;
252 // check first caracter to determine if following information is data or comments
253 if(test=="#")
254 {
255 // skip the line by beginning a new while loop.
256 continue;
257 }
258 // check if line is empty
259 else if(line.empty())
260 {
261 // skip the line by beginning a new while loop.
262 continue;
263 }
264 //
265 // end of the check
266
267 // transform the line into a iss
268 std::istringstream iss(line);
269
270 // Variables to be filled by the input file
271 double tDummy;
272 double eDummy;
273
274 // fill the variables with the content of the line
275 iss>>tDummy>>eDummy;
276
277 // SI : mandatory Vecm initialization
278
279 // Fill two vectors contained in maps of types:
280 // [materialName][particleName]=vector
281 // [materialName][particleName][T]=vector
282 // to list all the incident energies (tValues) and all the output energies (eValues) within the file
283 //
284 // Check if we already have the current T value in the vector.
285 // If not then add it
286 if (tDummy != tValuesVec[materialName][particleName].back())
287 {
288 // Add the current T value
289 tValuesVec[materialName][particleName].push_back(tDummy);
290
291 // Make it correspond to a default zero E value
292 eValuesVect[materialName][particleName][tDummy].push_back(0.);
293 }
294
295 // Put the differential cross section value of the input file within the diffCrossSectionData map
296 iss>>diffCrossSectionData[materialName][particleName][tDummy][eDummy];
297
298 // If the current E value (eDummy) is different from the one already registered in the eVector then add it to the vector
299 if (eDummy != eValuesVect[materialName][particleName][tDummy].back()) eValuesVect[materialName][particleName][tDummy].push_back(eDummy);
300 }
301}
302
303//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
304
306 const G4String& materialName,
307 const G4ParticleDefinition* p,
308 G4double ekin,
309 G4double /*emin*/,
310 G4double /*emax*/)
311{
312 if (verboseLevel > 3)
313 G4cout << "Calling CrossSectionPerVolume() of G4DNAPTBElasticModel" << G4endl;
314
315 // Get the name of the current particle
316 const G4String& particleName = p->GetParticleName();
317
318 // set killBelowEnergy value for current material
319 fKillBelowEnergy = GetLowELimit(materialName, particleName);
320
321 // initialise the return value (cross section) to zero
322 G4double sigma(0);
323
324 // check if we are below the high energy limit
325 if (ekin < GetHighELimit(materialName, particleName) )
326 {
327 // This is used to kill the particle if its kinetic energy is below fKillBelowEnergy.
328 // If the energy is lower then we return a maximum cross section and thus the SampleSecondaries method will be called for sure.
329 // SampleSecondaries will remove the particle from the simulation.
330 //
331 //SI : XS must not be zero otherwise sampling of secondaries method ignored
332 if (ekin < fKillBelowEnergy) return DBL_MAX;
333
334 // Get the tables with the cross section data
335 TableMapData* tableData = GetTableData();
336
337 // Retrieve the cross section value
338 sigma = (*tableData)[materialName][particleName]->FindValue(ekin);
339 }
340
341 if (verboseLevel > 2)
342 {
343 G4cout << "__________________________________" << G4endl;
344 G4cout << "°°° G4DNAPTBElasticModel - XS INFO START" << G4endl;
345 G4cout << "°°° Kinetic energy(eV)=" << ekin/eV << " particle : " << particleName << G4endl;
346 G4cout << "°°° Cross section per molecule (cm^2)=" << sigma/cm/cm << G4endl;
347 G4cout << "°°° G4DNAPTBElasticModel - XS INFO END" << G4endl;
348 }
349
350 // Return the cross section
351 return sigma;
352}
353
354//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
355
356void G4DNAPTBElasticModel::SampleSecondaries(std::vector<G4DynamicParticle*>* /*fvect*/,
357 const G4MaterialCutsCouple* /*couple*/,
358 const G4String& materialName,
359 const G4DynamicParticle* aDynamicElectron,
360 G4ParticleChangeForGamma* particleChangeForGamma,
361 G4double /*tmin*/,
362 G4double /*tmax*/)
363{
364 if (verboseLevel > 3)
365 G4cout << "Calling SampleSecondaries() of G4DNAPTBElasticModel" << G4endl;
366
367 G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
368
369 const G4String& particleName = aDynamicElectron->GetParticleDefinition()->GetParticleName();
370
371 // set killBelowEnergy value for material
372 fKillBelowEnergy = GetLowELimit(materialName, particleName);
373
374 // If the particle (electron here) energy is below the kill limit then we remove it from the simulation
375 if (electronEnergy0 < fKillBelowEnergy)
376 {
377 particleChangeForGamma->SetProposedKineticEnergy(0.);
378 particleChangeForGamma->ProposeTrackStatus(fStopAndKill);
379 particleChangeForGamma->ProposeLocalEnergyDeposit(electronEnergy0);
380 }
381 // If we are above the kill limite and below the high limit then we proceed
382 else if (electronEnergy0>= fKillBelowEnergy && electronEnergy0 < GetHighELimit(materialName, particleName) )
383 {
384 // Random sampling of the cosTheta
385 G4double cosTheta = RandomizeCosTheta(electronEnergy0, materialName);
386
387 // Random sampling of phi
388 G4double phi = 2. * pi * G4UniformRand();
389
390 G4ThreeVector zVers = aDynamicElectron->GetMomentumDirection();
391 G4ThreeVector xVers = zVers.orthogonal();
392 G4ThreeVector yVers = zVers.cross(xVers);
393
394 G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
395 G4double yDir = xDir;
396 xDir *= std::cos(phi);
397 yDir *= std::sin(phi);
398
399 // Particle direction after ModelInterface
400 G4ThreeVector zPrikeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
401
402 // Give the new direction
403 particleChangeForGamma->ProposeMomentumDirection(zPrikeVers.unit()) ;
404
405 // Update the energy which does not change here
406 particleChangeForGamma->SetProposedKineticEnergy(electronEnergy0);
407 }
408}
409
410//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
411
412G4double G4DNAPTBElasticModel::Theta
413(G4ParticleDefinition * particleDefinition, G4double k, G4double integrDiff, const G4String& materialName)
414{
415 G4double theta = 0.;
416 G4double valueT1 = 0;
417 G4double valueT2 = 0;
418 G4double valueE21 = 0;
419 G4double valueE22 = 0;
420 G4double valueE12 = 0;
421 G4double valueE11 = 0;
422 G4double xs11 = 0;
423 G4double xs12 = 0;
424 G4double xs21 = 0;
425 G4double xs22 = 0;
426 G4String particleName = particleDefinition->GetParticleName();
427
428 if (particleDefinition == G4Electron::ElectronDefinition())
429 {
430 std::vector<double>::iterator t2 = std::upper_bound(tValuesVec[materialName][particleName].begin(),tValuesVec[materialName][particleName].end(), k);
431 std::vector<double>::iterator t1 = t2-1;
432
433 std::vector<double>::iterator e12 = std::upper_bound(eValuesVect[materialName][particleName][(*t1)].begin(),eValuesVect[materialName][particleName][(*t1)].end(), integrDiff);
434 std::vector<double>::iterator e11 = e12-1;
435
436 std::vector<double>::iterator e22 = std::upper_bound(eValuesVect[materialName][particleName][(*t2)].begin(),eValuesVect[materialName][particleName][(*t2)].end(), integrDiff);
437 std::vector<double>::iterator e21 = e22-1;
438
439 valueT1 =*t1;
440 valueT2 =*t2;
441 valueE21 =*e21;
442 valueE22 =*e22;
443 valueE12 =*e12;
444 valueE11 =*e11;
445
446 xs11 = diffCrossSectionData[materialName][particleName][valueT1][valueE11];
447 xs12 = diffCrossSectionData[materialName][particleName][valueT1][valueE12];
448 xs21 = diffCrossSectionData[materialName][particleName][valueT2][valueE21];
449 xs22 = diffCrossSectionData[materialName][particleName][valueT2][valueE22];
450 }
451
452 if (xs11==0 && xs12==0 && xs21==0 && xs22==0) return (0.);
453
454 theta = QuadInterpolator ( valueE11, valueE12,
455 valueE21, valueE22,
456 xs11, xs12,
457 xs21, xs22,
458 valueT1, valueT2,
459 k, integrDiff );
460
461 return theta;
462}
463
464//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
465
466G4double G4DNAPTBElasticModel::LinLogInterpolate(G4double e1,
467 G4double e2,
468 G4double e,
469 G4double xs1,
470 G4double xs2)
471{
472 G4double d1 = std::log(xs1);
473 G4double d2 = std::log(xs2);
474 G4double value = std::exp(d1 + (d2 - d1)*(e - e1)/ (e2 - e1));
475 return value;
476}
477
478//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
479
480G4double G4DNAPTBElasticModel::LinLinInterpolate(G4double e1,
481 G4double e2,
482 G4double e,
483 G4double xs1,
484 G4double xs2)
485{
486 G4double d1 = xs1;
487 G4double d2 = xs2;
488 G4double value = (d1 + (d2 - d1)*(e - e1)/ (e2 - e1));
489 return value;
490}
491
492//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
493
494G4double G4DNAPTBElasticModel::LogLogInterpolate(G4double e1,
495 G4double e2,
496 G4double e,
497 G4double xs1,
498 G4double xs2)
499{
500 G4double a = (std::log10(xs2)-std::log10(xs1)) / (std::log10(e2)-std::log10(e1));
501 G4double b = std::log10(xs2) - a*std::log10(e2);
502 G4double sigma = a*std::log10(e) + b;
503 G4double value = (std::pow(10.,sigma));
504 return value;
505}
506
507//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
508
509G4double G4DNAPTBElasticModel::QuadInterpolator(G4double e11, G4double e12,
510 G4double e21, G4double e22,
511 G4double xs11, G4double xs12,
512 G4double xs21, G4double xs22,
513 G4double t1, G4double t2,
514 G4double t, G4double e)
515{
516 // Log-Log
517 /*
518 G4double interpolatedvalue1 = LogLogInterpolate(e11, e12, e, xs11, xs12);
519 G4double interpolatedvalue2 = LogLogInterpolate(e21, e22, e, xs21, xs22);
520 G4double value = LogLogInterpolate(t1, t2, t, interpolatedvalue1, interpolatedvalue2);
521
522
523 // Lin-Log
524 G4double interpolatedvalue1 = LinLogInterpolate(e11, e12, e, xs11, xs12);
525 G4double interpolatedvalue2 = LinLogInterpolate(e21, e22, e, xs21, xs22);
526 G4double value = LinLogInterpolate(t1, t2, t, interpolatedvalue1, interpolatedvalue2);
527*/
528
529 // Lin-Lin
530 G4double interpolatedvalue1 = LinLinInterpolate(e11, e12, e, xs11, xs12);
531 G4double interpolatedvalue2 = LinLinInterpolate(e21, e22, e, xs21, xs22);
532 G4double value = LinLinInterpolate(t1, t2, t, interpolatedvalue1, interpolatedvalue2);
533
534 return value;
535}
536
537//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
538
539G4double G4DNAPTBElasticModel::RandomizeCosTheta(G4double k, const G4String& materialName)
540{
541 G4double integrdiff=0;
542 G4double uniformRand=G4UniformRand();
543 integrdiff = uniformRand;
544
545 G4double theta=0.;
546 G4double cosTheta=0.;
547 theta = Theta(G4Electron::ElectronDefinition(),k/eV,integrdiff, materialName);
548
549 cosTheta= std::cos(theta*pi/180);
550
551 return cosTheta;
552}
553
554
555
556
const char * G4FindDataDir(const char *)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
@ fStopAndKill
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector unit() const
Hep3Vector orthogonal() const
Hep3Vector cross(const Hep3Vector &) const
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4String &materialName, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
CrossSectionPerVolume This method is mandatory for any model class. It finds and return the cross sec...
G4DNAPTBElasticModel(const G4String &applyToMaterial="all", const G4ParticleDefinition *p=0, const G4String &nam="DNAPTBElasticModel")
G4DNAPTBElasticModel Constructor.
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4String &materialName, const G4DynamicParticle *, G4ParticleChangeForGamma *particleChangeForGamma, G4double tmin, G4double tmax)
SampleSecondaries Method called after CrossSectionPerVolume if the process is the one which is select...
virtual ~G4DNAPTBElasticModel()
~G4DNAPTBElasticModel Destructor
virtual void Initialise(const G4ParticleDefinition *particle, const G4DataVector &, G4ParticleChangeForGamma *fpChangeForGamme=nullptr)
Initialise Mandatory method for every model class. The material/particle for which the model can be u...
const G4ThreeVector & GetMomentumDirection() const
const G4ParticleDefinition * GetParticleDefinition() const
G4double GetKineticEnergy() const
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:88
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
const G4String & GetParticleName() const
The G4VDNAModel class.
Definition: G4VDNAModel.hh:50
TableMapData * GetTableData()
GetTableData.
Definition: G4VDNAModel.hh:193
void SetHighELimit(const G4String &material, const G4String &particle, G4double lim)
SetHighEnergyLimit.
Definition: G4VDNAModel.hh:169
std::map< G4String, std::map< G4String, G4DNACrossSectionDataSet *, std::less< G4String > > > TableMapData
Definition: G4VDNAModel.hh:183
void AddCrossSectionData(G4String materialName, G4String particleName, G4String fileCS, G4String fileDiffCS, G4double scaleFactor)
AddCrossSectionData Method used during the initialization of the model class to add a new material....
Definition: G4VDNAModel.cc:58
void SetLowELimit(const G4String &material, const G4String &particle, G4double lim)
SetLowEnergyLimit.
Definition: G4VDNAModel.hh:177
G4double GetHighELimit(const G4String &material, const G4String &particle)
GetHighEnergyLimit.
Definition: G4VDNAModel.hh:153
void LoadCrossSectionData(const G4String &particleName)
LoadCrossSectionData Method to loop on all the registered materials in the model and load the corresp...
Definition: G4VDNAModel.cc:75
G4double GetLowELimit(const G4String &material, const G4String &particle)
GetLowEnergyLimit.
Definition: G4VDNAModel.hh:161
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
#define DBL_MAX
Definition: templates.hh:62