Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4VDNAPTBModel.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 #include "G4VDNAPTBModel.hh"
27 #include "G4SystemOfUnits.hh"
28 
29 G4VDNAPTBModel::G4VDNAPTBModel(const G4String& nam, const G4String& applyToMaterial)
30  : G4VEmModel(nam),
31  fStringOfMaterials(applyToMaterial)
32 {
33 
34 }
35 
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 }
52 
53 void G4VDNAPTBModel::AddCrossSectionData(const G4String& materialName, const G4String& particleName, const G4String& fileCS, const G4String& fileDiffCS, G4double scaleFactor)
54 {
55  fModelMaterialData.push_back(MaterialData(materialName, particleName, fileCS, fileDiffCS, scaleFactor) );
56 }
57 
58 void G4VDNAPTBModel::AddCrossSectionData(const G4String& materialName, const G4String& particleName, const G4String& fileCS, G4double scaleFactor)
59 {
60  fModelMaterialData.push_back(MaterialData(materialName, particleName, fileCS, "", scaleFactor) );
61 }
62 
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 }
135 
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 }
143 
144 void G4VDNAPTBModel::EnableForMaterialAndParticle(const G4String& materialName, const G4String& particleName)
145 {
146  fTableData[particleName][materialName] = 0;
147 }
148 
149 std::vector<G4String> G4VDNAPTBModel::BuildApplyToMatVect(const G4String& materials)
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 }
183 
185  const G4String& particleName,
186  const G4String& file, G4double scaleFactor)
187 {
188  fTableData[particleName][materialName] = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, eV, scaleFactor);
189  fTableData[particleName][materialName]->LoadData(file);
190 }
191 
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 }
244 
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 }
264 
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 }
278 
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 }
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
const G4String & GetName() const
Definition: G4Material.hh:178
std::map< G4String, std::map< G4String, G4double > > fLowEnergyLimits
G4bool IsParticleExistingInModel(const G4String &particleName)
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:587
std::vector< G4Material * > G4MaterialTable
std::map< G4int, G4double > fLowEnergyLimitsRuntime
void LoadCrossSectionData(const G4String &particleName)
G4bool IsMaterialDefine(const G4String &materialName)
int G4int
Definition: G4Types.hh:78
G4bool IsMaterialExistingInModelForParticle(const G4String &particleName, const G4String &materialName)
#define G4UniformRand()
Definition: Randomize.hh:97
void AddCrossSectionData(const G4String &materialName, const G4String &particleName, const G4String &fileCS, const G4String &fileDiffCS, G4double scaleFactor)
const XML_Char int const XML_Char * value
Definition: expat.h:331
std::vector< G4String > BuildApplyToMatVect(const G4String &materials)
bool G4bool
Definition: G4Types.hh:79
static constexpr double eV
Definition: G4SIunits.hh:215
void EnableForMaterialAndParticle(const G4String &materialName, const G4String &particleName)
virtual size_t NumberOfComponents(void) const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4int RandomSelectShell(G4double k, const G4Material *material)
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
const G4String & GetName() const
Definition: G4VEmModel.hh:794
double G4double
Definition: G4Types.hh:76
std::vector< MaterialData > fModelMaterialData
virtual ~G4VDNAPTBModel()
G4VDNAPTBModel(const G4String &nam, const G4String &applyToMaterial)
static const G4double pos