Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Static Public Member Functions
XrayFluoSiLiDetectorType Class Reference

#include <XrayFluoSiLiDetectorType.hh>

Inheritance diagram for XrayFluoSiLiDetectorType:
XrayFluoVDetectorType

Public Member Functions

 ~XrayFluoSiLiDetectorType ()
 
G4String GetDetectorMaterial ()
 
G4double ResponseFunction (G4double)
 
G4double GetInfData (G4double, G4double, G4int)
 
G4double GetSupData (G4double, G4double, G4int)
 
void LoadResponseData (G4String)
 
void LoadEfficiencyData (G4String)
 
- Public Member Functions inherited from XrayFluoVDetectorType
virtual ~XrayFluoVDetectorType ()
 

Static Public Member Functions

static XrayFluoSiLiDetectorTypeGetInstance ()
 

Additional Inherited Members

- Protected Member Functions inherited from XrayFluoVDetectorType
 XrayFluoVDetectorType ()
 

Detailed Description

Definition at line 50 of file XrayFluoSiLiDetectorType.hh.

Constructor & Destructor Documentation

XrayFluoSiLiDetectorType::~XrayFluoSiLiDetectorType ( )

Definition at line 59 of file XrayFluoSiLiDetectorType.cc.

60 {
61  std::map<G4int,G4DataVector*,std::less<G4int> >::iterator pos;
62 
63  for (pos = energyMap.begin(); pos != energyMap.end(); pos++)
64  {
65  G4DataVector* dataSet = (*pos).second;
66  delete dataSet;
67  dataSet = 0;
68  }
69  for (pos = dataMap.begin(); pos != dataMap.end(); pos++)
70  {
71  G4DataVector* dataSet = (*pos).second;
72  delete dataSet;
73  dataSet = 0;
74  }
75 
76  delete interpolation4;
77 
78 }

Member Function Documentation

G4String XrayFluoSiLiDetectorType::GetDetectorMaterial ( )
virtual

Implements XrayFluoVDetectorType.

Definition at line 79 of file XrayFluoSiLiDetectorType.cc.

80 {
81  return detectorMaterial;
82 }
G4double XrayFluoSiLiDetectorType::GetInfData ( G4double  ,
G4double  random,
G4int  posIndex 
)
virtual

Implements XrayFluoVDetectorType.

Definition at line 202 of file XrayFluoSiLiDetectorType.cc.

Referenced by ResponseFunction().

203 {
204  G4double value = 0.;
205  G4int zMin = 1;
206  G4int zMax = 6;
207 
208  G4int Z = posIndex;
209 
210  if (Z<zMin) {Z=zMin;}
211  if (Z>zMax) {Z=zMax;}
212 
213  if (Z >= zMin && Z <= zMax)
214  {
215  std::map<G4int,G4DataVector*,std::less<G4int> >::const_iterator pos;
216  pos = energyMap.find(Z);
217  std::map<G4int,G4DataVector*,std::less<G4int> >::const_iterator posData;
218  posData = dataMap.find(Z);
219 
220 
221  if (pos!= energyMap.end())
222  {
223  G4DataVector energySet = *((*pos).second);
224  G4DataVector dataSet = *((*posData).second);
225  G4int nData = energySet.size();
226 
227 
228  G4double dataSum = 0;
229  G4double partSum = 0;
230  G4int index = 0;
231 
232  // if data is not perfectly normalized (it may happen)
233  // rnadom number is renormalized, in case it is higer
234  //than the sum of all energies => segmentation fault.
235 
236  for (G4int i = 0; i<nData; i++){
237  dataSum += dataSet[i];
238  }
239 
240  G4double normRandom = random*dataSum;
241 
242 
243  while (normRandom> partSum)
244  {
245 
246  partSum += dataSet[index];
247  index++;
248  }
249 
250 
251  if (index >= 0 && index < nData)
252  {
253  value = energySet[index];
254 
255  }
256 
257  }
258  }
259  return value;
260 }
int G4int
Definition: G4Types.hh:78
const XML_Char int const XML_Char * value
double G4double
Definition: G4Types.hh:76
XrayFluoSiLiDetectorType * XrayFluoSiLiDetectorType::GetInstance ( void  )
static
G4double XrayFluoSiLiDetectorType::GetSupData ( G4double  ,
G4double  random,
G4int  posIndex 
)
virtual

Implements XrayFluoVDetectorType.

Definition at line 261 of file XrayFluoSiLiDetectorType.cc.

Referenced by ResponseFunction().

262 {
263  G4double value = 0.;
264  G4int zMin = 1;
265  G4int zMax = 6;
266  G4int Z = (posIndex+1);
267 
268  if (Z<zMin) {Z=zMin;}
269  if (Z>zMax) {Z=zMax;}
270  if (Z >= zMin && Z <= zMax)
271  {
272  std::map<G4int,G4DataVector*,std::less<G4int> >::const_iterator pos;
273  pos = energyMap.find(Z);
274  std::map<G4int,G4DataVector*,std::less<G4int> >::const_iterator posData;
275  posData = dataMap.find(Z);
276  if (pos!= energyMap.end())
277  {
278  G4DataVector energySet = *((*pos).second);
279  G4DataVector dataSet = *((*posData).second);
280  G4int nData = energySet.size();
281  G4double dataSum = 0;
282  G4double partSum = 0;
283  G4int index = 0;
284 
285  // if data is not perfectly normalized (it may happen)
286  // rnadom number is renormalized, in case it is higer
287  //than the sum of all energies => segmentation fault.
288 
289  for (G4int i = 0; i<nData; i++){
290  dataSum += dataSet[i];
291  }
292 
293  G4double normRandom = random*dataSum;
294 
295 
296  while (normRandom> partSum)
297  {
298  partSum += dataSet[index];
299  index++;
300  }
301 
302 
303  if (index >= 0 && index < nData)
304  {
305  value = energySet[index];
306  }
307  }
308  }
309  return value;
310 }
int G4int
Definition: G4Types.hh:78
const XML_Char int const XML_Char * value
double G4double
Definition: G4Types.hh:76
void XrayFluoSiLiDetectorType::LoadEfficiencyData ( G4String  fileName)
virtual

Implements XrayFluoVDetectorType.

Definition at line 392 of file XrayFluoSiLiDetectorType.cc.

References python.hepunit::keV.

393 {
394  char* path = getenv("XRAYDATA");
395  G4String dirFile;
396  if (path) {
397  G4String pathString(path);
398  dirFile = pathString + "/" + fileName;
399  }
400  else{
401  path = getenv("PWD");
402  G4String pathString(path);
403  dirFile = pathString + "/" + fileName;
404  }
405 
406  interpolation4 = new G4LogLogInterpolation();
407  efficiencySet = new XrayFluoDataSet(1,fileName,interpolation4,keV,1);
408 }
void XrayFluoSiLiDetectorType::LoadResponseData ( G4String  fileName)
virtual

Implements XrayFluoVDetectorType.

Definition at line 311 of file XrayFluoSiLiDetectorType.cc.

References test::a, FatalException, G4Exception(), and python.hepunit::keV.

312 {
313  std::ostringstream ost;
314 
315  ost << fileName<<".dat";
316 
317  G4String name = ost.str();
318 
319  char* path = getenv("XRAYDATA");
320 
321  G4String pathString(path);
322  G4String dirFile = pathString + "/" + name;
323  std::ifstream file(dirFile);
324  std::filebuf* lsdp = file.rdbuf();
325 
326  if (! (lsdp->is_open()) )
327  {
329  execp << "XrayFluoSiLiDetectorType - data file: " + dirFile + " not found";
330  G4Exception("XrayFluoSiLiDetectorType::LoadResponseData()","example-xray_fluorescence07",
331  FatalException, execp);
332  }
333  G4double a = 0;
334  G4int k = 1;
335  G4int q = 0;
336 
337  G4int Z = 1;
338  G4DataVector* energies = new G4DataVector;
340 
341  do
342  {
343  file >> a;
344  G4int nColumns = 2;
345  if (a == -1)
346  {
347  if (q == 0)
348  {
349  // End of a data set
350  energyMap[Z] = energies;
351  dataMap[Z] = data;
352  // Start of new shell data set
353  energies = new G4DataVector;
354  data = new G4DataVector;
355  Z++;
356  }
357  q++;
358  if (q == nColumns)
359  {
360  q = 0;
361  }
362  }
363  else if (a == -2)
364  {
365  // End of file; delete the empty vectors
366  //created when encountering the last -1 -1 row
367  delete energies;
368  delete data;
369 
370  }
371  else
372  {
373  // 1st column is energy
374  if(k%nColumns != 0)
375  {
376  G4double e = a * keV;
377  energies->push_back(e);
378  k++;
379  }
380  else if (k%nColumns == 0)
381  {
382  // 2nd column is data
383 
384  data->push_back(a);
385  k = 1;
386  }
387  }
388  } while (a != -2); // end of file
389  file.close();
390 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
const XML_Char * name
int G4int
Definition: G4Types.hh:78
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
double G4double
Definition: G4Types.hh:76
const XML_Char const XML_Char * data
G4double XrayFluoSiLiDetectorType::ResponseFunction ( G4double  energy)
virtual

Implements XrayFluoVDetectorType.

Definition at line 98 of file XrayFluoSiLiDetectorType.cc.

References XrayFluoDataSet::FindValue(), G4UniformRand, GetInfData(), GetSupData(), and python.hepunit::keV.

99 {
100 
101  G4double eMin = 1.500 *keV;
102  G4double eMax = 6.403 *keV;
103  G4double value = 0.;
104  G4double efficiency = 1.;
105 
106  const XrayFluoDataSet* dataSet = efficiencySet;
107  G4int id = 0;
108  G4DataVector energyVector;
109  energyVector.push_back(1.486* keV);
110  energyVector.push_back(1.740* keV);
111  energyVector.push_back(3.688* keV);
112  energyVector.push_back(4.510* keV);
113  energyVector.push_back(5.414* keV);
114  energyVector.push_back(6.404* keV);
115 
116  G4double infEnergy = 0 *keV;
117  G4double supEnergy = 10* keV;
118  G4int energyNumber = 0;
119  G4double random = G4UniformRand();
120 
121  if (energy>=eMin && energy <=eMax)
122  {
123 
124  for (G4int i=0; i<(G4int)energyVector.size(); i++){
125  if (energyVector[i]/keV < energy/keV){
126 
127  infEnergy = energyVector[i];
128  supEnergy = energyVector[i+1];
129 
130  energyNumber = i+1;
131 
132  }
133  }
134 
135 
136  G4double infData = GetInfData(energy, random, energyNumber);
137 
138  G4double supData = GetSupData(energy,random, energyNumber);
139 
140  value = (std::log10(infData)*std::log10(supEnergy/energy) +
141  std::log10(supData)*std::log10(energy/infEnergy)) /
142  std::log10(supEnergy/infEnergy);
143  value = std::pow(10,value);
144 
145 
146  }
147 // else if (energy<eMin || energy>eMax)
148 // {
149 // G4double infEnergy = eMin;
150 // G4double supEnergy = eMin/keV +1*keV;
151 
152 // G4double infData = GetInfData(eMin, random);
153 // G4double supData = GetSupData(eMin,random);
154 // value = (std::log10(infData)*std::log10(supEnergy/eMin) +
155 // std::log10(supData)*std::log10(eMin/infEnergy)) /
156 // std::log10(supEnergy/infEnergy);
157 // value = std::pow(10,value);
158 // value = value-eMin+ energy;
159 
160 
161 // }
162 
163 
164  else if (energy > eMax)
165  {
166 
167  energyNumber = 5;
168 
169  value = (GetSupData(energy, random, energyNumber))+(energy - 6.404* keV);
170 
171  }
172 
173 
174 
175  else
176  {
177 
178  energyNumber = 1;
179 
180  value = (GetInfData(energy, random, energyNumber))+(energy - 1.486* keV);
181 
182 
183  // G4double mean = -14.03 * eV + 1.0047*energy/eV;
184  // G4double stdDev = 35.38 * eV + 0.004385*energy/eV;
185 
186 
187  // value = (G4RandGauss::shoot(mean,stdDev))*eV;
188 
189  }
190  G4double RandomNum = G4UniformRand();
191 
192  efficiency = dataSet->FindValue(value,id);
193  if ( RandomNum>efficiency )
194  {
195  value = 0.;
196  }
197 
198  // G4cout << value << G4endl;
199  return value;
200 
201 }
int G4int
Definition: G4Types.hh:78
double precision function energy(A, Z)
Definition: dpm25nuc6.f:4106
#define G4UniformRand()
Definition: Randomize.hh:87
G4double GetInfData(G4double, G4double, G4int)
G4double GetSupData(G4double, G4double, G4int)
G4double FindValue(G4double e, G4int) const
const XML_Char int const XML_Char * value
double G4double
Definition: G4Types.hh:76

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