Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
XrayFluoHPGeDetectorType.cc
이 파일의 문서화 페이지로 가기
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 // $Id: XrayFluoHPGeDetectorType.cc
28 // GEANT4 tag $Name:
29 //
30 // Author: Alfonso Mantero (Alfonso.Mantero@ge.infn.it)
31 //
32 // History:
33 // -----------
34 // 16 Jul 2003 Alfonso Mantero Created
35 //
36 // -------------------------------------------------------------------
37 
38 #include <fstream>
39 #include <sstream>
40 
42 #include "XrayFluoDataSet.hh"
43 #include "G4SystemOfUnits.hh"
44 #include "G4UnitsTable.hh"
45 #include "G4DataVector.hh"
46 #include "G4LogLogInterpolation.hh"
47 #include "G4ios.hh"
48 #include "Randomize.hh"
49 
50 
52  detectorMaterial("HPGe"),efficiencySet(0)
53 {
54  LoadResponseData("response");
55 
56  LoadEfficiencyData("efficienza");
57 
58 
59 }
61 {
62  std::map<G4int,G4DataVector*,std::less<G4int> >::iterator pos;
63 
64  for (pos = energyMap.begin(); pos != energyMap.end(); pos++)
65  {
66  G4DataVector* dataSet = (*pos).second;
67  delete dataSet;
68  dataSet = 0;
69  }
70  for (pos = dataMap.begin(); pos != dataMap.end(); pos++)
71  {
72  G4DataVector* dataSet = (*pos).second;
73  delete dataSet;
74  dataSet = 0;
75  }
76 
77  delete interpolation4;
78 
79 }
81 {
82  return detectorMaterial;
83 }
84 
86 
88 
89 {
90  if (instance == 0)
91  {
93 
94  }
95  return instance;
96 }
97 
99 {
100 
101 
102  G4double eMin = 1* keV;
103  G4double eMax = 10*keV;
104  G4double value = 0.;
105  G4double efficiency = 1.;
106 
107  const XrayFluoDataSet* dataSet = efficiencySet;
108  G4int id = 0;
109 
110  G4double random = G4UniformRand();
111 
112  if (energy>=eMin && energy <=eMax)
113  {
114  G4double infEnergy = (G4int)(energy/keV)* keV;
115 
116  G4double supEnergy = ((G4int)(energy/keV) + 1)*keV;
117 
118  G4double infData = GetInfData(energy, random, 0);// 0 is not used
119 
120  G4double supData = GetSupData(energy,random, 0); // 0 is not used
121 
122  value = (std::log10(infData)*std::log10(supEnergy/energy) +
123  std::log10(supData)*std::log10(energy/infEnergy)) /
124  std::log10(supEnergy/infEnergy);
125  value = std::pow(10,value);
126  }
127  else if (energy<eMin)
128  {
129  G4double infEnergy = eMin;
130  G4double supEnergy = eMin/keV +1*keV;
131 
132  G4double infData = GetInfData(eMin, random, 0);
133  G4double supData = GetSupData(eMin,random, 0);
134  value = (std::log10(infData)*std::log10(supEnergy/eMin) +
135  std::log10(supData)*std::log10(eMin/infEnergy)) /
136  std::log10(supEnergy/infEnergy);
137  value = std::pow(10,value);
138  value = value-eMin+ energy;
139 
140 
141  }
142  else if (energy>eMax)
143  {
144  G4double infEnergy = eMax/keV - 1. *keV;
145  G4double supEnergy = eMax;
146 
147  G4double infData = GetInfData(eMax, random, 0);
148  G4double supData = GetSupData(eMax,random, 0);
149  value = (std::log10(infData)*std::log10(supEnergy/eMax) +
150  std::log10(supData)*std::log10(eMax/infEnergy)) /
151  std::log10(supEnergy/infEnergy);
152  value = std::pow(10,value);
153  value = value+energy- eMax;
154  }
155  G4double RandomNum = G4UniformRand();
156 
157  efficiency = dataSet->FindValue(value,id);
158  if ( RandomNum>efficiency )
159  {
160  value = 0.;
161  }
162 
163  return value;
164 }
166 {
167  G4double value = 0.;
168  G4int zMin = 1;
169  G4int zMax = 10;
170 
171  G4int Z = ((G4int)(energy/keV));
172 
173  if (Z<zMin) {Z=zMin;}
174  if (Z>zMax) {Z=zMax;}
175 
176  if (Z >= zMin && Z <= zMax)
177  {
178  std::map<G4int,G4DataVector*,std::less<G4int> >::const_iterator pos;
179  pos = energyMap.find(Z);
180  std::map<G4int,G4DataVector*,std::less<G4int> >::const_iterator posData;
181  posData = dataMap.find(Z);
182  if (pos!= energyMap.end())
183  {
184  G4DataVector energySet = *((*pos).second);
185  G4DataVector dataSet = *((*posData).second);
186  G4int nData = energySet.size();
187 
188  G4double partSum = 0;
189  G4int index = 0;
190 
191  while (random> partSum)
192  {
193  partSum += dataSet[index];
194  index++;
195  }
196 
197 
198  if (index >= 0 && index < nData)
199  {
200  value = energySet[index];
201 
202  }
203 
204  }
205  }
206  return value;
207 
208 }
210 {
211  G4double value = 0.;
212  G4int zMin = 1;
213  G4int zMax = 10;
214  G4int Z = ((G4int)(energy/keV)+1);
215 
216  if (Z<zMin) {Z=zMin;}
217  if (Z>zMax) {Z=zMax;}
218  if (Z >= zMin && Z <= zMax)
219  {
220  std::map<G4int,G4DataVector*,std::less<G4int> >::const_iterator pos;
221  pos = energyMap.find(Z);
222  std::map<G4int,G4DataVector*,std::less<G4int> >::const_iterator posData;
223  posData = dataMap.find(Z);
224  if (pos!= energyMap.end())
225  {
226  G4DataVector energySet = *((*pos).second);
227  G4DataVector dataSet = *((*posData).second);
228  G4int nData = energySet.size();
229  G4double partSum = 0;
230  G4int index = 0;
231 
232  while (random> partSum)
233  {
234  partSum += dataSet[index];
235  index++;
236  }
237 
238 
239  if (index >= 0 && index < nData)
240  {
241  value = energySet[index];
242  }
243  }
244  }
245  return value;
246 
247 }
249 {
250  std::ostringstream ost;
251 
252  ost << fileName<<".dat";
253 
254  G4String name = ost.str();
255 
256  char* path = getenv("XRAYDATA");
257 
258  G4String pathString(path);
259  G4String dirFile = pathString + "/" + name;
260  std::ifstream file(dirFile);
261  std::filebuf* lsdp = file.rdbuf();
262 
263  if (! (lsdp->is_open()) )
264  {
266  execp << "XrayFluoHPGeDetectorType - data file: " + dirFile + " not found";
267  G4Exception("XrayFluoHPGeDetectorType::LoadResponseData()","example-xray_fluorescence02",
268  FatalException, execp);
269  }
270  G4double a = 0;
271  G4int k = 1;
272  G4int q = 0;
273 
274  G4int Z = 1;
275  G4DataVector* energies = new G4DataVector;
277 
278  do
279  {
280  file >> a;
281  G4int nColumns = 2;
282  if (a == -1)
283  {
284  if (q == 0)
285  {
286  // End of a data set
287  energyMap[Z] = energies;
288  dataMap[Z] = data;
289  // Start of new shell data set
290  energies = new G4DataVector;
291  data = new G4DataVector;
292  Z++;
293  }
294  q++;
295  if (q == nColumns)
296  {
297  q = 0;
298  }
299  }
300  else if (a == -2)
301  {
302  // End of file; delete the empty vectors
303  //created when encountering the last -1 -1 row
304  delete energies;
305  delete data;
306 
307  }
308  else
309  {
310  // 1st column is energy
311  if(k%nColumns != 0)
312  {
313  G4double e = a * keV;
314  energies->push_back(e);
315  k++;
316  }
317  else if (k%nColumns == 0)
318  {
319  // 2nd column is data
320 
321  data->push_back(a);
322  k = 1;
323  }
324  }
325  } while (a != -2); // end of file
326  file.close();
327 }
328 
330 {
331  /*
332  char* path = getenv("XRAYDATA");
333  G4String dirFile;
334  if (path) {
335  G4String pathString(path);
336  dirFile = pathString + "/" + fileName;
337  }
338  else{
339  path = getenv("PWD");
340  G4String pathString(path);
341  dirFile = pathString + "/" + fileName;
342  }
343  */
345  efficiencySet = new XrayFluoDataSet(1,fileName,interpolation4,keV,1);
346 }
347 
const XML_Char * name
Definition: expat.h:151
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:45
std::vector< ExP01TrackerHit * > a
Definition: ExP01Classes.hh:33
static const G4double pos
static constexpr double keV
Definition: G4SIunits.hh:216
std::map< G4int, G4DataVector *, std::less< G4int > > energyMap
G4double GetSupData(G4double, G4double, G4int)
std::map< G4int, G4DataVector *, std::less< G4int > > dataMap
const XML_Char const XML_Char * data
Definition: expat.h:268
static XrayFluoHPGeDetectorType * GetInstance()
Float_t Z
static XrayFluoHPGeDetectorType * instance
double G4double
Definition: G4Types.hh:76
const XML_Char int const XML_Char * value
Definition: expat.h:331
const XrayFluoDataSet * efficiencySet
#define G4UniformRand()
Definition: Randomize.hh:53
double energy
Definition: plottest35.C:25
G4double GetInfData(G4double, G4double, G4int)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.hh:65
int G4int
Definition: G4Types.hh:78
G4double FindValue(G4double e, G4int) const
TFile * file
G4VDataSetAlgorithm * interpolation4