Geant4-11
G4DNAModelInterface.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// Contact authors: S. Meylan, C. Villagrasa
28//
29// email: sylvain.meylan@symalgo-tech.com, carmen.villagrasa@irsn.fr
30
33#include "G4SystemOfUnits.hh"
35
36
38 : G4VEmModel(nam), fName(nam), fpParticleChangeForGamma(0), fSampledMat("")
39{
40
41}
42
43//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
44
46{
47 // Loop on all the registered models to properly delete them (free the memory)
48 for(unsigned int i=0, ie = fRegisteredModels.size(); i<ie; ++i)
49 {
50 if(fRegisteredModels.at(i) != nullptr) delete fRegisteredModels.at(i);
51 }
52}
53
54//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
55
57 const G4DataVector& cuts)
58{
59 // Those two statements are necessary to override the energy limits set in the G4DNAProcesses (ionisation, elastic, etc...).
60 // Indeed, with the ModelInterface system, the model define themselves their energy limits per material and particle.
61 // Therefore, such a limit should not be in the G4DNAProcess classes.
62 //
65
67
68 // Loop on all the registered models to initialise them
69 for(unsigned int i=0, ie = fRegisteredModels.size(); i<ie; ++i)
70 {
71 fRegisteredModels.at(i)->Initialise(particle, cuts, fpParticleChangeForGamma);
72 }
73
74
75 // Build the [material][particle]=Models table
76 // used to retrieve the model corresponding to the current material/particle couple
78
80}
81
82//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
83
85 const G4ParticleDefinition* p,
86 G4double ekin,
87 G4double emin,
89{
90 // Method to return the crossSection * nbMoleculePerUnitVolume to the process class.
91 // Process class then calculates the path.
92 // The cross section is calculated in the registered model(s) and this class just call the method
93 // Two cases are handled here: normal material and composite material.
94 //
95 // Idea:
96 // *** Simple material ***
97 // Ask for the cross section of the chosen model.
98 // Multiply it by the number of medium molecules per volume unit.
99 // Return the value.
100 // *** Composite material ***
101 // Ask for the cross section of the chosen model for each component.
102 // Apply a factor to each cross section and sum the results. The factor is the molecule number of component per composite volume unit.
103 // The total cross section is returned.
104
105 // To reset the sampledMat variable.
106 // Can be used by user to retrieve current component
107 fSampledMat = "";
108
109 // This is the value to be sum up and to be returned at then end
110 G4double crossSectionTimesNbMolPerVol (0);
111
112 // Reset the map saving the material and the cumulated corresponding cross section
113 // Used in SampleSecondaries if the interaction is selected for the step and if the material is a composite
114 fMaterialCS.clear();
115
116 // This is the value to be used by SampleSecondaries
117 fCSsumTot = 0;
118
119 // *****************************
120 // Material is not a composite
121 // *****************************
122 //
123 if(material->GetMatComponents().empty())
124 {
125 // Get the material name
126 const G4String& materialName = material->GetName();
127
128 // Use the table to get the model
129 G4VDNAModel* model = GetDNAModel(materialName, p->GetParticleName(), ekin);
130
131 // Get the nunber of molecules per volume unit for that material
132 G4double nbOfMoleculePerVolumeUnit = GetNumMoleculePerVolumeUnitForMaterial(material);
133
134 // Calculate the cross section times the number of molecules
135 if(model != 0)
136 crossSectionTimesNbMolPerVol = nbOfMoleculePerVolumeUnit * model->CrossSectionPerVolume(material, materialName, p, ekin, emin, emax);
137 else // no model was selected, we are out of the energy ranges
138 crossSectionTimesNbMolPerVol = 0.;
139 }
140
141 // ********************************
142 // Material is a composite
143 // ********************************
144 //
145 else
146 {
147 // Copy the map in a local variable
148 // Otherwise we get segmentation fault and iterator pointing to nowhere: do not know why...
149 // Maybe MatComponents map is overrided by something somewhere ?
150 std::map<G4Material*, G4double> componentsMap = material->GetMatComponents();
151
152 // Retrieve the iterator
153 std::map<G4Material*, G4double>::const_iterator it = componentsMap.begin();
154
155 // Get the size
156 unsigned int componentNumber = componentsMap.size();
157
158 // Loop on all the components
159 //for(it = material->GetMatComponents().begin(); it!=material->GetMatComponents().end();++it)
160 for(unsigned int i=0; i<componentNumber; ++i)
161 {
162 // Get the current component
163 G4Material* component = it->first;
164
165 // Get the current component mass fraction
166 //G4double massFraction = it->second;
167
168 // Get the number of component molecules in a volume unit of composite material
169 G4double nbMoleculeOfComponentInCompositeMat = GetNumMolPerVolUnitForComponentInComposite(component, material);
170
171 // Get the current component name
172 const G4String componentName = component->GetName();
173
174 // Retrieve the model corresponding to the current component (ie material)
175 G4VDNAModel* model = GetDNAModel(componentName, p->GetParticleName(), ekin);
176
177 // Add the component part of the cross section to the cross section variable.
178 // The component cross section is multiplied by the total molecule number in the composite scaled by the mass fraction.
179 if(model != 0)
180 crossSectionTimesNbMolPerVol =
181 nbMoleculeOfComponentInCompositeMat * model->CrossSectionPerVolume(component, componentName, p, ekin, emin, emax);
182 else // no model was selected, we are out of the energy ranges
183 crossSectionTimesNbMolPerVol = 0.;
184
185 // Save the component name and its calculated crossSectionTimesNbMolPerVol
186 // To be used by sampling secondaries if the interaction is selected for the step
187 fMaterialCS[componentName] = crossSectionTimesNbMolPerVol;
188
189 // Save the component name and its calculated crossSectionTimesNbMolPerVol
190 // To be used by sampling secondaries if the interaction is selected for the step
191 fCSsumTot += crossSectionTimesNbMolPerVol;
192
193 // Move forward the iterator
194 ++it;
195 }
196
197 crossSectionTimesNbMolPerVol = fCSsumTot;
198
199 }
200
201 // return the cross section times the number of molecules
202 // the path of the interaction will be calculated using that value
203 return crossSectionTimesNbMolPerVol;
204}
205
206//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
207
208void G4DNAModelInterface::SampleSecondaries(std::vector<G4DynamicParticle*>* fVect,
209 const G4MaterialCutsCouple* couple,
210 const G4DynamicParticle* aDynamicParticle,
211 G4double tmin,
212 G4double tmax)
213{
214 // To call the sampleSecondaries method of the registered model(s)
215 // In the case of composite material, we need to choose a component to call the method from.
216 // To do so we use a random sampling on the crossSectionTimesNbMolPerVol used in CrossSectionPerVolume method.
217 // If we enter that method it means the corresponding interaction (and process) has been chosen for the current step.
218
219 G4String materialName;
220
221 // *******************************
222 // Material is not a composite
223 // *******************************
224 //
225 if(couple->GetMaterial()->GetMatComponents().empty())
226 {
227 materialName = couple->GetMaterial()->GetName();
228 }
229
230 // ****************************
231 // Material is a composite
232 // ****************************
233 //
234 else
235 {
236 // Material is a composite
237 // We need to select a component
238
239 // We select a random number between 0 and fCSSumTot
241
242 G4double cumulCS (0);
243
244 G4bool result = false;
245
246 // We loop on each component cumulated cross section
247 //
248 // Retrieve the iterators
249 std::map<const G4String , G4double>::const_iterator it = fMaterialCS.begin();
250 std::map<const G4String , G4double>::const_iterator ite = fMaterialCS.end();
251 // While this is true we do not have found our component.
252 while(rand>cumulCS)
253 {
254 // Check if the sampling is ok
255 if(it==ite)
256 {
257 G4Exception("G4DNAModelManager::SampleSecondaries","em0006",
259 "The random component selection has failed: we ran into the end of the map without having a selected component");
260 return; // to make some compilers happy
261 }
262
263 // Set the cumulated value for the iteration
264 cumulCS += it->second;
265
266 // Check if we have reach the material to be selected
267 // The DBL_MAX is here to take into account a return DBL_MAX in CSPerVol for the elastic model
268 // to force elastic sampleSecondaries where the particle can be killed.
269 // Used when paticle energy is lower than limit.
270 if(rand<cumulCS || cumulCS >= DBL_MAX)
271 {
272 // we have our selected material
273 materialName = it->first;
274 result = true;
275 break;
276 }
277
278 // make the iterator move forward
279 ++it;
280 }
281
282 // Check that we get a result
283 if(!result)
284 {
285 // it is possible to end up here if the return DBL_MAX of CSPerVol in the elastic model is not taken into account
286
287 G4Exception("G4DNAModelManager::SampleSecondaries","em0006",
289 "The random component selection has failed: while loop ended without a selected component.");
290 return; // to make some compilers happy
291 }
292
293 }
294
295 // **************************************
296 // Call the SampleSecondaries method
297 // **************************************
298
299 // Rename material if modified NIST material
300 // This is needed when material is obtained from G4MaterialCutsCouple
301 if(materialName.find("_MODIFIED")!=G4String::npos)
302 {
303 materialName = materialName.substr(0,materialName.size()-9);
304 }
305
306 fSampledMat = materialName;
307
308 G4VDNAModel* model = GetDNAModel(materialName,
309 aDynamicParticle->GetParticleDefinition()->GetParticleName(),
310 aDynamicParticle->GetKineticEnergy() );
311 //fMaterialParticleModelTable[materialName][aDynamicParticle->GetDefinition()->GetParticleName()][0];
312
313 model->SampleSecondaries(fVect, couple, materialName, aDynamicParticle, fpParticleChangeForGamma, tmin, tmax);
314}
315
317{
318 fRegisteredModels.push_back(model);
319}
320
322{
323 G4DNADummyModel* dummyWrapper = new G4DNADummyModel("G4_WATER", particle, model->GetName(), model);
324
325 RegisterModel(dummyWrapper);
326}
327
328//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
329
331{
332 // Method to build a map: [material][particle] = Model*.
333 // The map is used to retrieve the correct model for the current particle/material couple.
334
335 // Get the current particle name
336 const G4String& pName = p->GetParticleName();
337
338 // Retrieve the iterator
339 G4MaterialTable::iterator it;
340
341 // Loop on all materials registered in the simulation
342 for(it = G4Material::GetMaterialTable()->begin(); it!=G4Material::GetMaterialTable()->end(); ++it)
343 {
344 // Get the material pointer
345 G4Material* mat = *it;
346
347 // Get the map
348 std::map<G4Material*, G4double> componentMap = mat->GetMatComponents();
349
350 // Get the number of component within the composite
351 unsigned int compositeSize = componentMap.size();
352
353 // Check that the material is not a composite material
354 if(componentMap.empty())
355 {
356 // Get the material name
357 const G4String& matName = mat->GetName();
358
359 // Insert the model in the table.
360 InsertModelInTable(matName, pName);
361 }
362 // if the material is a composite material then we need to loop on all its components to register them
363 else
364 {
365 // Retrieve the component map begin iterator
366 std::map<G4Material*, G4double>::const_iterator itComp = componentMap.begin();
367
368 // Loop on all the components of the material
369 //for(itComp = mat->GetMatComponents().begin(); itComp != eitComp; ++itComp)
370 for(unsigned int k=0; k<compositeSize; ++k)
371 {
372 G4Material* component = itComp->first;
373
374// // Check that the component is not itself a composite
375// if(component->GetMatComponents().size()!=0)
376// {
377// std::ostringstream oss;
378// oss<<"Material "<<mat->GetName()<<" is a composite and its component ";
379// oss<<component->GetName()<<" is also a composite material. Building composite with other composites is not implemented yet";
380// oss<<G4endl;
381// G4Exception("G4DNAModelManager::BuildMaterialParticleModelTable","em0006",
382// FatalException, oss.str().c_str());
383// return; // to make some compilers happy
384// }
385
386 // Get the current component name
387 const G4String compName = component->GetName();
388
389 // If there is a model then insert the model corresponding to the component in the table
390 // contains a if statement to check we have not registered the material as a component or a normal material before.
391 InsertModelInTable(compName, pName);
392
393 // move forward the iterator
394 ++itComp;
395 }
396 }
397 }
398}
399
400//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
401
403{
404 // To be sure the G4DNAMolecularMaterial is initialized
406
408
409 // Loop on all the materials inside the "materialTable"
410 for(size_t i=0, ie=materialTable->size(); i<ie; i++)
411 {
412 // Current material
413 G4Material* currentMaterial = materialTable->at(i);
414
415 // Current material name
416 const G4String& currentMatName = currentMaterial->GetName();
417
418 // Will the material be used in this interface instance ?
419 // Loop on all the materials that can be dealt with in this class
420 MaterialParticleModelTable::iterator it = fMaterialParticleModelTable.begin();
421 MaterialParticleModelTable::iterator ite = fMaterialParticleModelTable.end();
422 for(; it != ite; it++)
423 {
424 const G4String& materialName = it->first;
425
426 if(materialName == currentMatName)
427 {
428 const std::vector<double>* numMolPerVolForMat = G4DNAMolecularMaterial::Instance()->GetNumMolPerVolTableFor(currentMaterial);
429 fMaterialMolPerVol[materialName] = numMolPerVolForMat;
430 }
431 }
432 }
433}
434
435//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
436
438{
439 // To insert the model(s) in the table Material Particule -> Model(s)
440
441 // First, we need to check if the current material has already been inserted in the table.
442 // This is possible because of the composite material. We could add a component M1 and then try to add the independant M1 material.
443 // This case must be avoided. Checking if M1 is already in the table is the way to avoid it.
444 //
445 // Chech if the current material and particle are already in the table.
446 // If they are: do nothing.
447 // If they are not: add the model(s)
448 //
449 // Check for the material
451 {
452 // Check for the particle
453 if(fMaterialParticleModelTable[matName].find(pName) == fMaterialParticleModelTable[matName].end())
454 {
455 G4int modelNbForMaterial (0);
456
457 // Loop on all models registered in the simulation to check:
458 // 1- if they can be applied to the current material
459 // 2- if they can be applied to the current particle
460 for(unsigned int i=0, ie=fRegisteredModels.size(); i<ie; ++i)
461 {
462 // check if the model is correct for material and particle (previous 1 and 2)
463 if(fRegisteredModels[i]->IsParticleExistingInModelForMaterial(pName, matName))
464 {
465 // if yes then add the model in the map
466 fMaterialParticleModelTable[matName][pName].push_back(fRegisteredModels[i]);
467
468 // and add one to the "there is a model" material flag
469 ++modelNbForMaterial;
470 }
471 }
472
473 // The model(s) applicable to the currently selected material should be in the map.
474 // We check if there are several models for the material.
475 if(modelNbForMaterial>1)
476 {
477 // If there are several models for a given material and particle couple it could be
478 // because of the energy ranges. We will check if the energy ranges are coherent.
479
480 // Get the models (vector)
481 std::vector<G4VDNAModel*>& models = fMaterialParticleModelTable[matName][pName];
482
483 // Declare a map to sort the limits (G4double) and a model "id" (G4int).
484 // The model id is created on the fly here.
485 // The idea is to fill a map with [limit] = counter. This map will be auto-sorted
486 // and we will check by iterating on it that the counter order is maintained.
487
488 // Delcare the map
489 std::map<G4double, G4int, std::less<G4double> > sortMap;
490
491 G4double smallDiff = 0.01 *eV;
492
493 // Loop on all the model for the current couple
494 // and fill a map with [lim] = modelNumber
495 for(unsigned int ii=0, em=models.size(); ii<em; ++ii)
496 {
497 G4double lowLim = models[ii]->GetLowELimit(matName, pName);
498 G4double highLim = models[ii]->GetHighELimit(matName, pName);
499
500 if(sortMap.find(lowLim) != sortMap.end() )
501 {
502 lowLim += smallDiff;
503 }
504
505 sortMap[lowLim] = ii;
506
507 if(sortMap.find(highLim) != sortMap.end() )
508 {
509 highLim -= smallDiff;
510 }
511
512 sortMap[highLim] = ii;
513 }
514
515 // The map has been created and ordered at this point.
516 // We will check the map order.
517
518 // Loop on the sortMap with iterator and check the order is correct.
519 std::map<G4double, G4int>::iterator it = sortMap.begin();
520
521 // First energy limit value
522 G4double dummyLim = it->first - smallDiff;
523
524 // Loop on all the models again.
525 // The goal is to check if for each limit pairs we have the same model number
526 // and that the upper and lower limit are consistent.
527 for(unsigned int ii=0, eii=models.size(); ii<eii; ++ii)
528 {
529 G4double lim1 = it->first - smallDiff;
530 G4int count1 = it->second;
531
532 // Iterate
533 ++it;
534
535 G4double lim2 = it->first + smallDiff;
536 G4int count2 = it->second;
537
538 // Iterate
539 ++it;
540
541 // Check model number and energy limit consistency
542 // std::abs(dummyLim - lim1) > 1.*eV because we cannot do (dummyLim != lim1)
543 // without experimenting precision loss. Therefore, the std::abs(...) > tolerance is the usual way of avoiding
544 // the issue.
545 if( (count1 != count2) || ( std::abs(dummyLim - lim1) > 1.*eV ) )
546 {
547 // Error
548
549 std::ostringstream oss;
550 oss<<"The material "<<matName<<" and the particle "<<pName;
551 oss<<" have several models registered for the "<<fName<<" interaction and their energy ranges ";
552 oss<<"do not match. \nEnergy ranges: \n";
553
554 for(int iii=0, eiii=models.size(); iii<eiii; ++iii)
555 {
556 oss<<models[iii]->GetName()<<"\n";
557 oss<<"low: "<<models[iii]->GetLowELimit(matName, pName)/eV<<" eV \n";
558 oss<<"high: "<<models[iii]->GetHighELimit(matName, pName)/eV<<" eV \n";
559 }
560
561 G4Exception("G4DNAModelManager::InsertModelInTable","em0006",
562 FatalException, oss.str().c_str());
563 return; // to make some compilers happy
564 }
565
566 dummyLim = lim2;
567 }
568
569 // If we are here then everything was ok.
570 }
571 // no model for the material case
572 else if(modelNbForMaterial==0)
573 {
574// std::ostringstream oss;
575// oss<<"The material "<<matName<<" and the particle "<<pName;
576// oss<<" does not have any model registered for the "<<fName<<" interaction. ";
577
578// G4Exception("G4DNAModelManager::InsertModelInTable","em0006",
579// FatalException, oss.str().c_str());
580// return; // to make some compilers happy
581 }
582 }
583 }
584}
585
587{
588 // Output pointer
589 G4VDNAModel* model = 0;
590
591 // Get a reference to all the models for the couple (material and particle)
592 std::vector<G4VDNAModel*>& models = fMaterialParticleModelTable[material][particle];
593
594 // We must choose one of the model(s) accordingly to the particle energy and the model energy range(s)
595
596 //G4bool isOneModelSelected = false;
597
598 // Loop on all the models within the models vector and check if ekin is within the energy range.
599 for(int i=0, ie=models.size(); i<ie; ++i)
600 {
601 // ekin is in the energy range: we select the model and stop the loop.
602 if( ekin >= models[i]->GetLowELimit(material, particle)
603 && ekin < models[i]->GetHighELimit(material, particle) )
604 {
605 // Select the model
606 model = models[i];
607
608 // Boolean flag
609 //isOneModelSelected = true;
610
611 // Quit the for loop
612 break;
613 }
614
615 // ekin is not in the energy range: we continue the loop.
616 }
617
618// // If no model was selected then fatal error
619// if(!isOneModelSelected)
620// {
621// G4String msg = "No model has ";
622// msg += ekin/eV;
623// msg += " eV in its energy range. Therefore nothing was selected.";
624
625// G4Exception("G4DNAModelManager::GetDNAModel","em0006",
626// FatalException,
627// msg);
628// }
629
630 // Return a pointer to the selected model
631 return model;
632}
633
634//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
635
637{
638 return fMaterialMolPerVol[mat->GetName()]->at(mat->GetIndex() );
639}
640
641//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
642
644{
645 return fMaterialMolPerVol[component->GetName() ]->at(composite->GetIndex() );
646}
static const G4double emax
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::vector< G4Material * > G4MaterialTable
static constexpr double eV
Definition: G4SIunits.hh:201
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4UniformRand()
Definition: Randomize.hh:52
virtual ~G4DNAModelInterface()
~G4DNAModelManager Destructor
std::map< G4String, const std::vector< double > * > fMaterialMolPerVol
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
CrossSectionPerVolume Method called by the process and used to call the CrossSectionPerVolume method ...
std::vector< G4VDNAModel * > fRegisteredModels
vector containing all the registered models
G4String fSampledMat
for the user to retrieve selected material/component
G4double GetNumMoleculePerVolumeUnitForMaterial(const G4Material *mat)
virtual void Initialise(const G4ParticleDefinition *particle, const G4DataVector &cuts)
Initialise Initialise method to call all the initialise methods of the registered models.
void BuildMaterialParticleModelTable(const G4ParticleDefinition *p)
BuildMaterialParticleModelTable Method used to build a map allowing the code to quickly retrieve the ...
std::map< const G4String, G4double > fMaterialCS
map used to share information between CrossSectionPerVolume and SampleSecondaries
const G4String fName
name of the interaction
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *fVect, const G4MaterialCutsCouple *couple, const G4DynamicParticle *aDynamicElectron, G4double tmin, G4double tmax)
SampleSecondaries Used to call the SampleSecondaries method of the registered models....
G4double GetNumMolPerVolUnitForComponentInComposite(const G4Material *component, const G4Material *composite)
G4VDNAModel * GetDNAModel(const G4String &material, const G4String &particle, G4double ekin)
GetDNAModel.
void InsertModelInTable(const G4String &matName, const G4String &pName)
InsertModelInTable Used to put a model in the table after performing some checks.
void RegisterModel(G4VDNAModel *model)
RegisterModel Method used to associate a model with the interaction.
G4double fCSsumTot
value which contains the sum of all the component cross sections in case of a composite material
G4DNAModelInterface(const G4String &nam)
G4DNAModelManager Constructor.
MaterialParticleModelTable fMaterialParticleModelTable
map: [materialName][particleName] = vector of models
G4ParticleChangeForGamma * fpParticleChangeForGamma
pointer used to change the characteristics of the current particle
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 G4ParticleDefinition * GetParticleDefinition() const
G4double GetKineticEnergy() const
const G4Material * GetMaterial() const
const std::map< G4Material *, G4double > & GetMatComponents() const
Definition: G4Material.hh:233
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:672
const G4String & GetName() const
Definition: G4Material.hh:173
size_t GetIndex() const
Definition: G4Material.hh:256
const G4String & GetParticleName() const
The G4VDNAModel class.
Definition: G4VDNAModel.hh:50
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4String &materialName, const G4DynamicParticle *, G4ParticleChangeForGamma *particleChangeForGamma, G4double tmin=0, G4double tmax=DBL_MAX)=0
SampleSecondaries Each model must implement SampleSecondaries to decide if a particle will be created...
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4String &materialName, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)=0
CrossSectionPerVolume Every model must implement its own CrossSectionPerVolume method....
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:767
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:123
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:774
const G4String & GetName() const
Definition: G4VEmModel.hh:837
string material
Definition: eplot.py:19
#define DBL_MAX
Definition: templates.hh:62