Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4VDNAPTBModel Class Referenceabstract

#include <G4VDNAPTBModel.hh>

Inheritance diagram for G4VDNAPTBModel:
Collaboration diagram for G4VDNAPTBModel:

Classes

struct  MaterialData
 

Public Member Functions

 G4VDNAPTBModel (const G4String &nam, const G4String &applyToMaterial)
 
virtual ~G4VDNAPTBModel ()
 
virtual void Initialise (const G4ParticleDefinition *particle, const G4DataVector &cuts)=0
 
G4bool IsMaterialDefine (const G4String &materialName)
 
G4bool IsParticleExistingInModel (const G4String &particleName)
 
G4bool IsMaterialExistingInModelForParticle (const G4String &particleName, const G4String &materialName)
 
void SetHighELimit (const G4String &material, const G4String &particle, G4double lim)
 
void SetLowELimit (const G4String &material, const G4String &particle, G4double lim)
 
G4double GetHighELimit (const G4String &material, const G4String &particle)
 
G4double GetLowELimit (const G4String &material, const G4String &particle)
 
virtual G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)=0
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0, G4double tmax=DBL_MAX)=0
 
G4double GetHighELimit (const G4Material *material)
 
G4double GetLowELimit (const G4Material *material)
 
void SetHighELimit (const G4Material *material, G4double lim)
 
void SetLowELimit (const G4Material *material, G4double lim)
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
virtual G4double GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ComputeCrossSectionPerShell (const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ChargeSquareRatio (const G4Track &)
 
virtual G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void StartTracking (G4Track *)
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector
< G4EmElementSelector * > * 
GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
virtual G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
G4double CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4int SelectIsotopeNumber (const G4Element *)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
void SetAngularDistribution (G4VEmAngularDistribution *)
 
G4double HighEnergyLimit () const
 
G4double LowEnergyLimit () const
 
G4double HighEnergyActivationLimit () const
 
G4double LowEnergyActivationLimit () const
 
G4double PolarAngleLimit () const
 
G4double SecondaryThreshold () const
 
G4bool LPMFlag () const
 
G4bool DeexcitationFlag () const
 
G4bool ForceBuildTableFlag () const
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy)
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 

Protected Types

typedef std::map< G4String,
std::map< G4String,
G4DNACrossSectionDataSet
*, std::less< G4String > > > 
TableMapData
 

Protected Member Functions

TableMapDataGetTableData ()
 
G4DNACrossSectionDataSetGetSigmaData (const G4Material *material)
 
std::vector< G4StringBuildApplyToMatVect (const G4String &materials)
 
void ReadAndSaveCSFile (const G4String &materialName, const G4String &particleName, const G4String &file, G4double scaleFactor)
 
G4int RandomSelectShell (G4double k, const G4Material *material)
 
void AddCrossSectionData (const G4String &materialName, const G4String &particleName, const G4String &fileCS, const G4String &fileDiffCS, G4double scaleFactor)
 
void AddCrossSectionData (const G4String &materialName, const G4String &particleName, const G4String &fileCS, G4double scaleFactor)
 
void LoadCrossSectionData (const G4String &particleName)
 
virtual void ReadDiffCSFile (const G4String &materialName, const G4String &particleName, const G4String &path, const G4double scaleFactor)
 
void EnableForMaterialAndParticle (const G4String &materialName, const G4String &particleName)
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Protected Attributes

const G4String fStringOfMaterials
 
TableMapData fTableData
 
std::vector< MaterialDatafModelMaterialData
 
std::map< G4String, std::map
< G4String, G4double > > 
fLowEnergyLimits
 
std::map< G4String, std::map
< G4String, G4double > > 
fHighEnergyLimits
 
std::map< G4int,
G4DNACrossSectionDataSet * > 
fTableDataRuntime
 
std::map< G4int, G4doublefLowEnergyLimitsRuntime
 
std::map< G4int, G4doublefHighEnergyLimitsRuntime
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 

Additional Inherited Members

- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

Detailed Description

Definition at line 36 of file G4VDNAPTBModel.hh.

Member Typedef Documentation

typedef std::map<G4String, std::map<G4String,G4DNACrossSectionDataSet*, std::less<G4String> > > G4VDNAPTBModel::TableMapData
protected

Definition at line 91 of file G4VDNAPTBModel.hh.

Constructor & Destructor Documentation

G4VDNAPTBModel::G4VDNAPTBModel ( const G4String nam,
const G4String applyToMaterial 
)

Definition at line 29 of file G4VDNAPTBModel.cc.

30  : G4VEmModel(nam),
31  fStringOfMaterials(applyToMaterial)
32 {
33 
34 }
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:68
const G4String fStringOfMaterials
G4VDNAPTBModel::~G4VDNAPTBModel ( )
virtual

Definition at line 36 of file G4VDNAPTBModel.cc.

37 {
38  // Clean fTableData
39  std::map<G4String, std::map<G4String,G4DNACrossSectionDataSet*,std::less<G4String> > >::iterator posOuter;
40  std::map<G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator posInner;
41  // iterate on each particle
42  for (posOuter = fTableData.begin(); posOuter != fTableData.end(); ++posOuter)
43  {
44  // iterate on each material
45  for(posInner = posOuter->second.begin(); posInner != posOuter->second.end(); ++posInner)
46  {
47  G4DNACrossSectionDataSet* table = posInner->second;
48  if(table != 0) delete table;
49  }
50  }
51 }
TableMapData fTableData

Member Function Documentation

void G4VDNAPTBModel::AddCrossSectionData ( const G4String materialName,
const G4String particleName,
const G4String fileCS,
const G4String fileDiffCS,
G4double  scaleFactor 
)
protected

Definition at line 53 of file G4VDNAPTBModel.cc.

54 {
55  fModelMaterialData.push_back(MaterialData(materialName, particleName, fileCS, fileDiffCS, scaleFactor) );
56 }
std::vector< MaterialData > fModelMaterialData
void G4VDNAPTBModel::AddCrossSectionData ( const G4String materialName,
const G4String particleName,
const G4String fileCS,
G4double  scaleFactor 
)
protected

Definition at line 58 of file G4VDNAPTBModel.cc.

59 {
60  fModelMaterialData.push_back(MaterialData(materialName, particleName, fileCS, "", scaleFactor) );
61 }
std::vector< MaterialData > fModelMaterialData
std::vector< G4String > G4VDNAPTBModel::BuildApplyToMatVect ( const G4String materials)
protected

Definition at line 149 of file G4VDNAPTBModel.cc.

150 {
151  // output material vector
152  std::vector<G4String> materialVect;
153 
154  // if we don't find any "/" then it means we only have one "material" (could be the "all" option)
155  if(materials.find("/")==std::string::npos)
156  {
157  // we add the material to the output vector
158  materialVect.push_back(materials);
159  }
160  // if we have several materials listed in the string then we must retrieve them
161  else
162  {
163  G4String materialsNonIdentified = materials;
164 
165  while(materialsNonIdentified.find_first_of("/") != std::string::npos)
166  {
167  // we select the first material and stop at the "/" caracter
168  G4String mat = materialsNonIdentified.substr(0, materialsNonIdentified.find_first_of("/"));
169  materialVect.push_back(mat);
170 
171  // we remove the previous material from the materialsNonIdentified string
172  materialsNonIdentified = materialsNonIdentified.substr(materialsNonIdentified.find_first_of("/")+1,
173  materialsNonIdentified.size()-materialsNonIdentified.find_first_of("/"));
174  }
175 
176  // we don't find "/" anymore, it means we only have one material string left
177  // we get it
178  materialVect.push_back(materialsNonIdentified);
179  }
180 
181  return materialVect;
182 }

Here is the caller graph for this function:

virtual G4double G4VDNAPTBModel::CrossSectionPerVolume ( const G4Material material,
const G4ParticleDefinition p,
G4double  ekin,
G4double  emin,
G4double  emax 
)
pure virtual

Reimplemented from G4VEmModel.

void G4VDNAPTBModel::EnableForMaterialAndParticle ( const G4String materialName,
const G4String particleName 
)
protected

Definition at line 144 of file G4VDNAPTBModel.cc.

145 {
146  fTableData[particleName][materialName] = 0;
147 }
TableMapData fTableData
G4double G4VDNAPTBModel::GetHighELimit ( const G4String material,
const G4String particle 
)
inline

Definition at line 59 of file G4VDNAPTBModel.hh.

59 {return fHighEnergyLimits[particle][material];}
std::map< G4String, std::map< G4String, G4double > > fHighEnergyLimits
G4double G4VDNAPTBModel::GetHighELimit ( const G4Material material)
inline

Definition at line 79 of file G4VDNAPTBModel.hh.

79 {return fHighEnergyLimitsRuntime.at(material->GetIndex() );}
size_t GetIndex() const
Definition: G4Material.hh:262
std::map< G4int, G4double > fHighEnergyLimitsRuntime

Here is the call graph for this function:

G4double G4VDNAPTBModel::GetLowELimit ( const G4String material,
const G4String particle 
)
inline

Definition at line 60 of file G4VDNAPTBModel.hh.

60 {return fLowEnergyLimits[particle][material];}
std::map< G4String, std::map< G4String, G4double > > fLowEnergyLimits
G4double G4VDNAPTBModel::GetLowELimit ( const G4Material material)
inline

Definition at line 80 of file G4VDNAPTBModel.hh.

80 {return fLowEnergyLimitsRuntime.at(material->GetIndex() );}
size_t GetIndex() const
Definition: G4Material.hh:262
std::map< G4int, G4double > fLowEnergyLimitsRuntime

Here is the call graph for this function:

G4DNACrossSectionDataSet* G4VDNAPTBModel::GetSigmaData ( const G4Material material)
inlineprotected

Definition at line 140 of file G4VDNAPTBModel.hh.

140 {return fTableDataRuntime.at(material->GetIndex() );}
size_t GetIndex() const
Definition: G4Material.hh:262
std::map< G4int, G4DNACrossSectionDataSet * > fTableDataRuntime

Here is the call graph for this function:

TableMapData* G4VDNAPTBModel::GetTableData ( )
inlineprotected

Definition at line 139 of file G4VDNAPTBModel.hh.

139 {return &fTableData;}
TableMapData fTableData
virtual void G4VDNAPTBModel::Initialise ( const G4ParticleDefinition particle,
const G4DataVector cuts 
)
pure virtual

Implements G4VEmModel.

G4bool G4VDNAPTBModel::IsMaterialDefine ( const G4String materialName)

Definition at line 245 of file G4VDNAPTBModel.cc.

246 {
247  // Check if the given material is defined in the simulation
248 
249  G4bool exist (false);
250 
251  double matTableSize = G4Material::GetMaterialTable()->size();
252 
253  for(int i=0;i<matTableSize;i++)
254  {
255  if(materialName == G4Material::GetMaterialTable()->at(i)->GetName())
256  {
257  exist = true;
258  return exist;
259  }
260  }
261 
262  return exist;
263 }
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:587
bool G4bool
Definition: G4Types.hh:79
const G4String & GetName() const
Definition: G4VEmModel.hh:794

Here is the call graph for this function:

G4bool G4VDNAPTBModel::IsMaterialExistingInModelForParticle ( const G4String particleName,
const G4String materialName 
)

Definition at line 279 of file G4VDNAPTBModel.cc.

280 {
281  // To check two things:
282  // 1- is the material existing in model ?
283  // 2- if yes, is the particle defined for that material ?
284 
285  if(IsParticleExistingInModel(particleName))
286  {
287  if (fTableData[particleName].find(materialName) == fTableData[particleName].end())
288  {
289  return false;
290  }
291  else return true;
292  }
293  else return false;
294 }
G4bool IsParticleExistingInModel(const G4String &particleName)
TableMapData fTableData

Here is the call graph for this function:

Here is the caller graph for this function:

G4bool G4VDNAPTBModel::IsParticleExistingInModel ( const G4String particleName)

Definition at line 265 of file G4VDNAPTBModel.cc.

266 {
267  // Check if the given material is defined in the current model class
268 
269  if (fTableData.find(particlelName) == fTableData.end())
270  {
271  return false;
272  }
273  else
274  {
275  return true;
276  }
277 }
TableMapData fTableData

Here is the caller graph for this function:

void G4VDNAPTBModel::LoadCrossSectionData ( const G4String particleName)
protected

Definition at line 63 of file G4VDNAPTBModel.cc.

64 {
65  G4String fileCS, fileDiffCS;
66  G4String materialName, particleNameData;
67  G4double scaleFactor;
68 
69  // construct applyToMatVect with materials specified by the user
70  std::vector<G4String> applyToMatVect = BuildApplyToMatVect(fStringOfMaterials);
71 
72  // iterate on each material contained into the fStringOfMaterials variable (through applyToMatVect)
73  for(unsigned int i=0;i<applyToMatVect.size();++i)
74  {
75  // We have selected a material coming from applyToMatVect
76  // We try to find if this material correspond to a model registered material
77  // If it is, then isMatFound becomes true
78  G4bool isMatFound = false;
79 
80  // We iterate on each model registered materials to load the CS data
81  // We have to do a for loop because of the "all" option
82  // applyToMatVect[i] == "all" implies applyToMatVect.size()=1 and we want to iterate on all registered materials
83  for(unsigned int j=0, je=fModelMaterialData.size();j<je;j++)
84  {
85  materialName = fModelMaterialData[j].fMaterial;
86  particleNameData = fModelMaterialData[j].fParticle;
87 
88  if( (applyToMatVect[i] == materialName || applyToMatVect[i] == "all")
89  && particleNameData==particleName )
90  {
91  isMatFound = true;
92  fileCS = fModelMaterialData[j].fCSFile;
93  fileDiffCS = fModelMaterialData[j].fDiffCSFile;
94  scaleFactor = fModelMaterialData[j].fScaleFactor;
95 
96  ReadAndSaveCSFile(materialName, particleNameData, fileCS, scaleFactor);
97 
98  if(fileDiffCS != "") ReadDiffCSFile(materialName, particleNameData, fileDiffCS, scaleFactor);
99 
100  }
101  }
102 
103  // check if we found a correspondance, if not: fatal error
104  if(!isMatFound)
105  {
106  std::ostringstream oss;
107  oss << applyToMatVect[i] << " material was not found. It means the material specified in the UserPhysicsList is not a model material for ";
108  oss << particleName;
109  G4Exception("G4VDNAPTBModel::LoadCrossSectionData","em0003",
110  FatalException, oss.str().c_str());
111  return;
112  }
113  }
114 
115  // ************************************************
116  // Generation of the data tables used at runtime
117  // ************************************************
118 
120 
121  // Loop on all the materials registered into the table
122  for(G4int i=0, ie=table->size(); i<ie; i++)
123  {
124  G4Material* material = table->at(i);
125 
126  if(IsMaterialExistingInModelForParticle(particleName, material->GetName() ) )
127  {
128  fTableDataRuntime[i] = fTableData[particleName][materialName];
129 
130  fLowEnergyLimitsRuntime[i] = fLowEnergyLimits[particleName][materialName];
131  fHighEnergyLimitsRuntime[i] = fHighEnergyLimits[particleName][materialName];
132  }
133  }
134 }
const G4String & GetName() const
Definition: G4Material.hh:178
std::map< G4String, std::map< G4String, G4double > > fLowEnergyLimits
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:587
std::vector< G4Material * > G4MaterialTable
std::map< G4int, G4double > fLowEnergyLimitsRuntime
int G4int
Definition: G4Types.hh:78
G4bool IsMaterialExistingInModelForParticle(const G4String &particleName, const G4String &materialName)
std::vector< G4String > BuildApplyToMatVect(const G4String &materials)
bool G4bool
Definition: G4Types.hh:79
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
const G4String fStringOfMaterials
std::map< G4String, std::map< G4String, G4double > > fHighEnergyLimits
void ReadAndSaveCSFile(const G4String &materialName, const G4String &particleName, const G4String &file, G4double scaleFactor)
std::map< G4int, G4double > fHighEnergyLimitsRuntime
TableMapData fTableData
virtual void ReadDiffCSFile(const G4String &materialName, const G4String &particleName, const G4String &path, const G4double scaleFactor)
std::map< G4int, G4DNACrossSectionDataSet * > fTableDataRuntime
double G4double
Definition: G4Types.hh:76
std::vector< MaterialData > fModelMaterialData

Here is the call graph for this function:

G4int G4VDNAPTBModel::RandomSelectShell ( G4double  k,
const G4Material material 
)
protected

Definition at line 192 of file G4VDNAPTBModel.cc.

193 {
194  G4int level = 0;
195 
196  std::map<G4int, G4DNACrossSectionDataSet*>::iterator pos;
197  pos = fTableDataRuntime.find(material->GetIndex() );
198 
199  if(pos != fTableDataRuntime.end())
200  {
201  G4DNACrossSectionDataSet* table = pos->second;
202 
203  if (table != 0)
204  {
205  G4double* valuesBuffer = new G4double[table->NumberOfComponents()];
206  const size_t n(table->NumberOfComponents());
207  size_t i(n);
208  G4double value = 0.;
209 
210  while (i>0)
211  {
212  i--;
213  valuesBuffer[i] = table->GetComponent(i)->FindValue(k);
214  value += valuesBuffer[i];
215  }
216 
217  value *= G4UniformRand();
218 
219  i = n;
220 
221  while (i > 0)
222  {
223  i--;
224 
225  if (valuesBuffer[i] > value)
226  {
227  delete[] valuesBuffer;
228  return i;
229  }
230  value -= valuesBuffer[i];
231  }
232 
233  if (valuesBuffer) delete[] valuesBuffer;
234 
235  }
236  }
237  else
238  {
239  G4Exception("G4VDNAPTBModel::RandomSelectShell","em0002",
240  FatalException,"Model not applicable to particle type.");
241  }
242  return level;
243 }
virtual G4double FindValue(G4double x, G4int componentId=0) const =0
virtual const G4VEMDataSet * GetComponent(G4int componentId) const
size_t GetIndex() const
Definition: G4Material.hh:262
int G4int
Definition: G4Types.hh:78
#define G4UniformRand()
Definition: Randomize.hh:97
const XML_Char int const XML_Char * value
Definition: expat.h:331
virtual size_t NumberOfComponents(void) const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
std::map< G4int, G4DNACrossSectionDataSet * > fTableDataRuntime
double G4double
Definition: G4Types.hh:76
static const G4double pos

Here is the call graph for this function:

void G4VDNAPTBModel::ReadAndSaveCSFile ( const G4String materialName,
const G4String particleName,
const G4String file,
G4double  scaleFactor 
)
protected

Definition at line 184 of file G4VDNAPTBModel.cc.

187 {
188  fTableData[particleName][materialName] = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, eV, scaleFactor);
189  fTableData[particleName][materialName]->LoadData(file);
190 }
static constexpr double eV
Definition: G4SIunits.hh:215
TableMapData fTableData

Here is the caller graph for this function:

void G4VDNAPTBModel::ReadDiffCSFile ( const G4String materialName,
const G4String particleName,
const G4String path,
const G4double  scaleFactor 
)
protectedvirtual

Definition at line 136 of file G4VDNAPTBModel.cc.

137 {
138  G4String text("ReadDiffCSFile must be implemented in the model class using a differential cross section data file");
139 
140  G4Exception("G4VDNAPTBModel::ReadDiffCSFile","em0003",
141  FatalException, text);
142 }
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

Here is the call graph for this function:

Here is the caller graph for this function:

virtual void G4VDNAPTBModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  ,
const G4MaterialCutsCouple ,
const G4DynamicParticle ,
G4double  tmin = 0,
G4double  tmax = DBL_MAX 
)
pure virtual

Implements G4VEmModel.

void G4VDNAPTBModel::SetHighELimit ( const G4String material,
const G4String particle,
G4double  lim 
)
inline

Definition at line 56 of file G4VDNAPTBModel.hh.

56 {fHighEnergyLimits[particle][material]=lim;}
std::map< G4String, std::map< G4String, G4double > > fHighEnergyLimits
void G4VDNAPTBModel::SetHighELimit ( const G4Material material,
G4double  lim 
)
inline

Definition at line 82 of file G4VDNAPTBModel.hh.

82 {fHighEnergyLimitsRuntime[material->GetIndex()]=lim;}
size_t GetIndex() const
Definition: G4Material.hh:262
std::map< G4int, G4double > fHighEnergyLimitsRuntime

Here is the call graph for this function:

void G4VDNAPTBModel::SetLowELimit ( const G4String material,
const G4String particle,
G4double  lim 
)
inline

Definition at line 57 of file G4VDNAPTBModel.hh.

57 {fLowEnergyLimits[particle][material]=lim;}
std::map< G4String, std::map< G4String, G4double > > fLowEnergyLimits
void G4VDNAPTBModel::SetLowELimit ( const G4Material material,
G4double  lim 
)
inline

Definition at line 83 of file G4VDNAPTBModel.hh.

83 {fLowEnergyLimitsRuntime[material->GetIndex()]=lim;}
size_t GetIndex() const
Definition: G4Material.hh:262
std::map< G4int, G4double > fLowEnergyLimitsRuntime

Here is the call graph for this function:

Member Data Documentation

std::map<G4String, std::map<G4String, G4double> > G4VDNAPTBModel::fHighEnergyLimits
protected

Definition at line 121 of file G4VDNAPTBModel.hh.

std::map<G4int, G4double> G4VDNAPTBModel::fHighEnergyLimitsRuntime
protected

Definition at line 133 of file G4VDNAPTBModel.hh.

std::map<G4String, std::map<G4String, G4double> > G4VDNAPTBModel::fLowEnergyLimits
protected

Definition at line 120 of file G4VDNAPTBModel.hh.

std::map<G4int, G4double> G4VDNAPTBModel::fLowEnergyLimitsRuntime
protected

Definition at line 132 of file G4VDNAPTBModel.hh.

std::vector<MaterialData> G4VDNAPTBModel::fModelMaterialData
protected

Definition at line 117 of file G4VDNAPTBModel.hh.

const G4String G4VDNAPTBModel::fStringOfMaterials
protected

Definition at line 93 of file G4VDNAPTBModel.hh.

TableMapData G4VDNAPTBModel::fTableData
protected

Definition at line 95 of file G4VDNAPTBModel.hh.

std::map<G4int, G4DNACrossSectionDataSet*> G4VDNAPTBModel::fTableDataRuntime
protected

Definition at line 129 of file G4VDNAPTBModel.hh.


The documentation for this class was generated from the following files: