source: trunk/examples/advanced/xray_fluorescence/src/XrayFluoSiLiDetectorType.cc @ 1187

Last change on this file since 1187 was 807, checked in by garnier, 16 years ago

update

File size: 9.8 KB
Line 
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: XrayFluoVdetectorType.cc
28// GEANT4 tag $Name:
29//
30// Author: Alfonso Mantero (Alfonso.Mantero@ge.infn.it)
31//
32// History:
33// -----------
34//  19 Jun 2003  Alfonso Mantero Created
35//
36// -------------------------------------------------------------------
37
38#include "XrayFluoSiLiDetectorType.hh"
39#include "XrayFluoDataSet.hh"
40#include "G4DataVector.hh"
41#include "G4LogLogInterpolation.hh"
42#include "G4ios.hh"
43#include <fstream>
44#include <sstream>
45#include "G4UnitsTable.hh"
46#include "Randomize.hh"
47
48
49XrayFluoSiLiDetectorType::XrayFluoSiLiDetectorType():
50  detectorMaterial("Silicon"),efficiencySet(0)
51{
52  LoadResponseData("SILIresponse");
53
54  LoadEfficiencyData("SILIefficiency");
55
56 
57}
58XrayFluoSiLiDetectorType::~XrayFluoSiLiDetectorType()
59{
60  std::map<G4int,G4DataVector*,std::less<G4int> >::iterator pos;
61
62  for (pos = energyMap.begin(); pos != energyMap.end(); pos++)
63    {
64      G4DataVector* dataSet = (*pos).second;
65      delete dataSet;
66      dataSet = 0;
67    }
68  for (pos = dataMap.begin(); pos != dataMap.end(); pos++)
69    {
70      G4DataVector* dataSet = (*pos).second;
71      delete dataSet;
72      dataSet = 0;
73    }
74
75  delete interpolation4; 
76
77}
78G4String XrayFluoSiLiDetectorType::GetDetectorMaterial()
79{
80  return detectorMaterial;
81}
82
83XrayFluoSiLiDetectorType* XrayFluoSiLiDetectorType::instance = 0;
84
85XrayFluoSiLiDetectorType* XrayFluoSiLiDetectorType::GetInstance()
86
87{
88  if (instance == 0)
89    {
90      instance = new XrayFluoSiLiDetectorType;
91     
92    }
93  return instance;
94}
95
96
97G4double XrayFluoSiLiDetectorType::ResponseFunction(G4double energy)
98{
99 
100  G4double eMin = 1.500 *keV;
101  G4double eMax = 6.403 *keV; 
102  G4double value = 0.;
103  G4double efficiency = 1.;
104 
105  const XrayFluoDataSet* dataSet = efficiencySet;
106  G4int id = 0;
107  G4DataVector energyVector;
108  energyVector.push_back(1.486* keV);
109  energyVector.push_back(1.740* keV);
110  energyVector.push_back(3.688* keV);
111  energyVector.push_back(4.510* keV);
112  energyVector.push_back(5.414* keV);
113  energyVector.push_back(6.404* keV);
114
115  G4double infEnergy = 0 *keV;
116  G4double supEnergy = 10* keV;
117  G4int energyNumber = 0;
118  G4double random = G4UniformRand();
119 
120  if (energy>=eMin && energy <=eMax)
121    {
122
123      for (G4int i=0; i<(G4int)energyVector.size(); i++){
124        if (energyVector[i]/keV < energy/keV){
125
126          infEnergy = energyVector[i];
127          supEnergy = energyVector[i+1];
128
129          energyNumber = i+1;
130
131        }
132      }
133     
134
135      G4double infData = GetInfData(energy, random, energyNumber);
136
137      G4double supData = GetSupData(energy,random, energyNumber);
138
139      value = (std::log10(infData)*std::log10(supEnergy/energy) +
140               std::log10(supData)*std::log10(energy/infEnergy)) / 
141        std::log10(supEnergy/infEnergy);
142      value = std::pow(10,value);
143
144
145    }
146//    else if (energy<eMin || energy>eMax)
147//      {
148//        G4double infEnergy = eMin;
149//        G4double supEnergy = eMin/keV +1*keV;
150 
151//        G4double infData = GetInfData(eMin, random);
152//        G4double supData = GetSupData(eMin,random);
153//        value = (std::log10(infData)*std::log10(supEnergy/eMin) +
154//             std::log10(supData)*std::log10(eMin/infEnergy)) /
155//      std::log10(supEnergy/infEnergy);
156//        value = std::pow(10,value);
157//        value = value-eMin+ energy;
158
159
160//      }
161
162
163  else if (energy > eMax)
164    {
165
166      energyNumber = 5;
167
168      value = (GetSupData(energy, random, energyNumber))+(energy - 6.404* keV);
169
170    }
171
172
173
174  else 
175    {
176 
177      energyNumber = 1;
178
179      value = (GetInfData(energy, random, energyNumber))+(energy - 1.486* keV);
180
181
182      //      G4double mean = -14.03 * eV + 1.0047*energy/eV;
183      //      G4double stdDev = 35.38 * eV + 0.004385*energy/eV;
184     
185     
186      //      value = (G4RandGauss::shoot(mean,stdDev))*eV;
187
188    }
189  G4double  RandomNum = G4UniformRand(); 
190 
191  efficiency = dataSet->FindValue(value,id);
192  if ( RandomNum>efficiency )
193    {
194      value = 0.;
195    }
196
197  //  G4cout << value << G4endl;
198  return value;
199
200}
201G4double XrayFluoSiLiDetectorType::GetInfData(G4double, G4double random, G4int posIndex)
202{
203  G4double value = 0.;
204  G4int zMin = 1;
205  G4int zMax = 6; 
206 
207  G4int Z = posIndex;
208 
209  if (Z<zMin) {Z=zMin;}
210  if (Z>zMax) {Z=zMax;}
211 
212  if (Z >= zMin && Z <= zMax)
213    {
214      std::map<G4int,G4DataVector*,std::less<G4int> >::const_iterator pos;
215      pos = energyMap.find(Z);
216      std::map<G4int,G4DataVector*,std::less<G4int> >::const_iterator posData;
217      posData = dataMap.find(Z);
218
219
220      if (pos!= energyMap.end())
221        {
222          G4DataVector energySet = *((*pos).second);
223          G4DataVector dataSet = *((*posData).second);
224          G4int nData = energySet.size();
225
226
227          G4double dataSum = 0;
228          G4double partSum = 0;
229          G4int index = 0;
230
231          // if data is not perfectly normalized (it may happen)
232          // rnadom number is renormalized, in case it is higer
233          //than the sum of all energies => segmentation fault.
234         
235          for  (G4int i = 0; i<nData; i++){
236            dataSum += dataSet[i];
237          }
238
239          G4double normRandom = random*dataSum;
240
241
242          while (normRandom> partSum)
243            {
244             
245              partSum += dataSet[index];
246              index++;
247            }
248         
249         
250          if (index >= 0 && index < nData)
251            {
252              value = energySet[index];
253             
254            }
255         
256        }
257    }
258  return value;
259}
260G4double XrayFluoSiLiDetectorType::GetSupData(G4double, G4double random, G4int posIndex)
261{
262  G4double value = 0.;
263  G4int zMin = 1;
264  G4int zMax = 6;
265  G4int Z = (posIndex+1);
266 
267  if (Z<zMin) {Z=zMin;}
268  if (Z>zMax) {Z=zMax;}
269  if (Z >= zMin && Z <= zMax)
270    {
271      std::map<G4int,G4DataVector*,std::less<G4int> >::const_iterator pos;
272      pos = energyMap.find(Z);
273      std::map<G4int,G4DataVector*,std::less<G4int> >::const_iterator posData;
274      posData = dataMap.find(Z);
275      if (pos!= energyMap.end())
276        {
277          G4DataVector energySet = *((*pos).second);
278          G4DataVector dataSet = *((*posData).second);
279          G4int nData = energySet.size();
280          G4double dataSum = 0;
281          G4double partSum = 0;
282          G4int index = 0;
283
284          // if data is not perfectly normalized (it may happen)
285          // rnadom number is renormalized, in case it is higer
286          //than the sum of all energies => segmentation fault.
287
288          for  (G4int i = 0; i<nData; i++){
289            dataSum += dataSet[i];
290          }       
291
292          G4double normRandom = random*dataSum;
293
294
295          while (normRandom> partSum)
296            {
297              partSum += dataSet[index];
298              index++;
299            }
300         
301         
302          if (index >= 0 && index < nData)
303            {
304              value = energySet[index];
305            }
306        }
307    }
308  return value;
309}
310void XrayFluoSiLiDetectorType::LoadResponseData(G4String fileName)
311{
312  std::ostringstream ost;
313   
314  ost << fileName<<".dat";
315 
316  G4String name = ost.str();
317 
318  char* path = getenv("XRAYDATA");
319 
320  G4String pathString(path);
321  G4String dirFile = pathString + "/" + name;
322  std::ifstream file(dirFile);
323  std::filebuf* lsdp = file.rdbuf();
324 
325  if (! (lsdp->is_open()) )
326    {
327      G4String excep = "XrayFluoSiLiDetectorType - data file: " + dirFile + " not found";
328      G4Exception(excep);
329        }
330  G4double a = 0;
331  G4int k = 1;
332  G4int s = 0;
333 
334  G4int Z = 1;
335  G4DataVector* energies = new G4DataVector;
336  G4DataVector* data = new G4DataVector;
337
338  do
339    {
340      file >> a;
341      G4int nColumns = 2;
342      if (a == -1)
343        {
344          if (s == 0)
345            {
346              // End of a  data set
347              energyMap[Z] = energies;
348              dataMap[Z] = data;
349              // Start of new shell data set
350              energies = new G4DataVector;
351              data = new G4DataVector;
352              Z++;         
353            }     
354          s++;
355          if (s == nColumns)
356            {
357              s = 0;
358            }
359        }
360      else if (a == -2)
361        {
362          // End of file; delete the empty vectors
363          //created when encountering the last -1 -1 row
364          delete energies;
365          delete data;
366         
367        }
368      else
369        {
370          // 1st column is energy
371          if(k%nColumns != 0)
372            {   
373              G4double e = a * keV;
374              energies->push_back(e);
375              k++;
376            }
377          else if (k%nColumns == 0)
378            {
379              // 2nd column is data
380             
381              data->push_back(a);
382              k = 1;
383            }
384        }
385    } while (a != -2); // end of file
386  file.close();   
387}
388
389void XrayFluoSiLiDetectorType::LoadEfficiencyData(G4String fileName)
390{
391  char* path = getenv("XRAYDATA");
392  G4String dirFile;
393  if (path) {
394    G4String pathString(path);
395    dirFile = pathString + "/" + fileName;
396  }
397  else{ 
398    path = getenv("PWD");
399    G4String pathString(path);
400    dirFile = pathString + "/" + fileName;
401  }
402
403  interpolation4 = new G4LogLogInterpolation();
404  efficiencySet = new XrayFluoDataSet(1,fileName,interpolation4,keV,1);
405}
406
Note: See TracBrowser for help on using the repository browser.