Geant4  10.02
G4DiscreteScatteringModel.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 //
27 #include "G4ParticleDefinition.hh"
28 #include "G4AtomicShells.hh"
29 #include "G4Track.hh"
30 #include "G4Material.hh"
31 #include "G4UnitsTable.hh"
32 #include "G4ProductionCutsTable.hh"
34 #include "G4ElementData.hh"
35 #include "G4Exp.hh"
36 #include "G4Log.hh"
37 
38 #include <iostream>
39 #include <string>
40 
41 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
42 
43 using namespace std;
44 using namespace CLHEP;
45 
48 
49 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
50 
52  : G4VEmModel("DiscrScat"), fParticleChange(nullptr),fAnalogModel("pwe"),
53  fNumAngles(iNumAngles), fLowEnergyLimit(2*keV)
54 {
55  SetHighEnergyLimit(100.*MeV);
57  if(IsMaster() && fCdf == nullptr) {
58  fCdf = new G4ElementData();
59  fTcs = new G4ElementData();
60  }
61 }
62 
63 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
64 
66 {
67  if(IsMaster()) {
68  delete fCdf;
69  delete fTcs;
70  fCdf = nullptr;
71  fTcs = nullptr;
72  }
73 }
74 
75 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
76 
78  const G4DataVector&)
79 {
80  if(nullptr == fParticleChange) {
82  }
83  if(!IsMaster()) { return; }
84 
85  G4cout << "G4DiscreteScatteringModel::Initialise start"<<G4endl;
86 
87  const G4int maxZ = 100;
88 
89  char *path = getenv("G4GBFPDATA");
90  if (!path)
91  {
92  G4Exception("G4DiscreteScatteringModel::Initialise","em0006",
93  FatalException,"G4GBFPDATA environment variable not set.");
94  return;
95  }
96 
97  std::ostringstream eFullFileName;
98  eFullFileName << path;
99 
100  G4ProductionCutsTable* theCoupleTable =
102 
103  G4int numOfCouples = theCoupleTable->GetTableSize();
104 
105  for(G4int i=0; i<numOfCouples; ++i)
106  {
107  const G4Material* material =
108  theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
109  const G4ElementVector* theElementVector = material->GetElementVector();
110  G4int nelm = material->GetNumberOfElements();
111 
112  for (G4int j=0; j<nelm; ++j)
113  {
114  G4int Z = G4lrint((*theElementVector)[j]->GetZ());
115  if(Z < 1) { Z = 1; }
116  else if(Z > maxZ) { Z = maxZ; }
117  if(!fTcs->GetElementData(Z)) { ReadData(Z, path); }
118  }
119  }
120  G4cout << "G4DiscreteScatteringModel::Initialise completed"<<G4endl;
121 }
122 
123 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
124 
126 {
127 
128  // Set strings with paths to tcs data file. Currently not a variable but this
129  // this will change
130  G4String fullFileName(argFileName);
131 
132  stringstream ss1;//create a stringstream
133  ss1 << Z;//add number to the stream
134  stringstream ss2;//create a stringstream
135  ss2 << fNumAngles;
136  G4String tcsPath = fullFileName+"/log_gbfp_"+fAnalogModel
137  +"_tcs_"+ss1.str()+"_"+ss2.str()+".dat";
138 
139  // Do some error checking and exiting if data does not exist
140  std::ifstream in(tcsPath, std::ifstream::binary|std::ifstream::in);
141  if (!in.is_open())
142  {
143  G4String message("Data file \"");
144  message+=tcsPath;
145  message+="\" not found";
146  G4Exception("G4DiscreteScatteringModel::LoadData","em0003",
147  FatalException,message);
148  return;
149  }
150  //G4cout<<"Reading in data file "<< tcsPath<<G4endl;
151 
152  // Create a temporary G4PhysicsVector object pointer
153  G4PhysicsVector* tempData = new G4PhysicsVector(false);
154 
155  // Use retrieve to read in the total cross section (tcs) data
156  tempData->Retrieve(in,true);
157 
158  // Convert tcs from cm^2 to mm^2
159  //tempData->ScaleVector(1.0, 100.0);
160 
161  // store pass this data to the tcs object and initialise for current element
162  fTcs->InitialiseForElement(Z,tempData);
163  in.close();
164 
165  // Set strings with paths to cdf data files. Currently not a variable
166  // but this this will change
167  G4String cdfPath = fullFileName+"/gbfp_"+fAnalogModel+"_cdf_"
168  +ss1.str()+"_"+ss2.str()+".dat";
169 
170  // Do some error checking and exiting if data does not exist
171  std::ifstream in2(cdfPath, std::ifstream::binary|std::ifstream::in);
172  if (!in2.is_open())
173  {
174  G4String message("Data file \"");
175  message+=cdfPath;
176  message+="\" not found";
177  G4Exception("G4DiscreteScatteringModel::LoadData","em0003",
178  FatalException,message);
179  return;
180  }
181 
182  //G4cout<<"Reading in data file "<< cdfPath<<G4endl;
183 
184  // The cumulative distribution functions (CDF) for energy E_j and X_i
185  // on (0,1) is C(E_j,X_i) and read-in/stored at this time.
186  // For the purposes of this model, each energy grid point where the CDF
187  // is evaluated is considered a component. The number of energy grid
188  // points is consistent with the fTcs data, so the following int is set
189  // by calling G4PhysicsVector::GetVectorLength().
190  G4int numEnergies = fTcs->GetElementData(Z)->GetVectorLength();
191 
192  // The ElementData object pointer is then initialized by
193  fCdf->InitialiseForComponent(Z, numEnergies);
194 
195  // Now the data files are read in for all energies. At each energy,
196  // there are fNumAngles angles and fNumAngles CDF values.
197 
198  std::vector<G4PhysicsVector*> tempDataCDF;
199  // Loop through each energy
200  for (int j=0; j<numEnergies; j++)
201  {
202  // Push back a new G4PhysicsVector for the jth energy
203  tempDataCDF.push_back(new G4PhysicsVector(false));
204 
205  // For use with David's PhysicsVector class
206  //tempDataCDF.push_back(new G4PhysicsVector(false,false,true));
207 
208  //tempDataCDF.push_back(new G4PhysicsVector(false,false,false));
209  //tempDataCDF.push_back(new G4PhysicsVector());
210  }
211 
212  // Open a temporary stream. The data for the jth energy group is copied
213  // to the temp file and then the temp file is sent to
214  // G4PhysicsVector::Retrieve(). Once the data is stored in the
215  // G4PhysicsVector, tempDataCDF, it is then passed to the ElementData, cdf,
216  // which is the container for all of the data for each element and energy.
217 
218  //static G4Mutex m = G4MUTEX_INITIALIZER;
219  //G4AutoLock l(&m);
220 
221  // Open the stream and call file "tempDataFile.dat"
222  std::ofstream file("tempDataFile.dat",
223  std::fstream::out | std::fstream::trunc);
224 
225  // Write to file the lower/upper bounds and the number of grid points
226  // for each column of data. The data in tempfile is 2xfNumAngles,
227  // hence fNumAngles fNumAngles.
228  // file<<lowerBound<<" "<<upperBound<<" "<<fNumAngles<<" "
229  // <<fNumAngles<<G4endl;
230  file<<"-1. 1. "<<fNumAngles<<" "<<fNumAngles<<G4endl;
231 
232  // Start while loop over the entire data file opened above
233  // e.g. pwe_cdf_79.dat. This data contains all of the data for each energy
234  G4int cntr = 0;
235  G4int j = 0;
236  G4double temp1, temp2;
237  while(in2>>temp1>>temp2)
238  {
239  // Write data to temporary file
240  file<<setprecision(16)<<temp1<<" "<<temp2<<G4endl;
241  cntr++;
242 
243  // When the first fNumAngles data points are copied to tempDataFile.dat,
244  // store data in tempDataCDF[j]. Then increment the energy index, j,
245  // close and clear the streams, and then reopen the streams such that
246  // the next fNumAngles data points are copied to tempDataFile.dat.
247  if (cntr==fNumAngles)
248  {
249  cntr=0;
250  std::ifstream inTemp("tempDataFile.dat",
251  std::ifstream::binary|std::ifstream::in);
252  tempDataCDF[j]->Retrieve(inTemp,true);
253  fCdf->AddComponent(Z,j,tempDataCDF[j]);
254  j++;
255  inTemp.close(), inTemp.clear(), file.close(), file.clear();
256  file.open("tempDataFile.dat",std::fstream::out | std::fstream::trunc);
257  file<<"-1. 1. "<<fNumAngles<<" "<<fNumAngles<<G4endl;
258  }
259  }
260 
261  in2.close();
262 
263  return;
264 }
265 
266 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
267 
269  std::vector<G4DynamicParticle*>*,
270  const G4MaterialCutsCouple* couple,
271  const G4DynamicParticle* p,
272  G4double cutEnergy,
273  G4double)
274 {
275  G4double E = p->GetKineticEnergy();
276 
277  if(E < fLowEnergyLimit) return;
278 
279  // Select random atom
281  E,cutEnergy,E)->GetZ());
282 
283  // Determine the energy bin
284  G4double logE = G4Log(E);
285  G4ThreeVector dir = p->GetMomentumDirection(); //old direction
286  G4int j = fTcs->GetElementData(Z)->FindBin(logE,0);
287 
288  //-------------------------------------------------------------------------
289  // it would be nice to have the following block of code in G4PhysicsVector.
290  // it could be a simple function.
291 
292  // v.GetVectorLength() - number of points
293  // v.Energy(size_t idx) - value x(i)
294  // v[i] - value y(i)
295 
296  // This is a monte carlo interpolation scheme
298  G4double e2 = fTcs->GetElementData(Z)->Energy(j+1);
299 
300  // This is a monte carlo interpolation scheme
301  G4double pie1 = (logE-e1)/(e2-e1);
302  G4double r = G4UniformRand();
303  if (r<pie1){ ++j; }
304 
305  //-------------------------------------------------------------------------
306  // it would be nice to have the following block of code in G4PhysicsVector.
307  //---------------
308  // Given the energy grid value associated with the DCS,
309  // sample a deflection cosine
310  r = G4UniformRand();
311  G4int k = -1;
312  // First test if the angle is the most probable angle
313  if (r>(*fCdf->GetComponentDataByIndex(Z,j)).Energy(fNumAngles-2))
314  { k = fNumAngles-1; }
315  // or the least probable... not sure why I do this (maybe because
316  // it is a simple check)
317  else if (r<=(*fCdf->GetComponentDataByIndex(Z,j)).Energy(0)) { k = 0; }
318  // if neither then loop through remaining angles, break when locating angle
319  else {
320  for (G4int i=fNumAngles-2; i>0; --i) {
321  if ( (r>(*fCdf->GetComponentDataByIndex(Z,j)).Energy(i-1))
322  && (r<= (*fCdf->GetComponentDataByIndex(Z,j)).Energy(i)) )
323  { k=i; break;}
324  }
325  }
326  // Throw an error if an angle was not sampled, data is probably no good
327  if(k<0)
328  {
329  G4cout << "G4DiscreteScatteringModel::SampleSecondaries():"
330  << " CDF was not inverted properly "<<k<<G4endl;
331  for (G4int i=0;i<fNumAngles;++i) {
332  G4cout<<i<<" "<<(*fCdf->GetComponentDataByIndex(Z,j)).Energy(i)<<" "<<r<<G4endl;
333  }
334  }
335  //-------------------------------------------------------------------------
336 
337  // Otherwise, go get the angle and pass it too local method GetNewDirection.
338  // Then do transformation and update fParticleChange.
339 
340  //G4cout<<"------------------"<<G4endl;
341  //G4cout<<"Sample Secondaries"<<G4endl;
342  //G4cout<<G4Exp(e1)<<" "<<E<<" "<<G4Exp(e2)<<" "<<pie1<<" "
343  // <<(*fCdf->GetComponentDataByIndex(Z,j)).Energy(k)<<G4endl;
344  //G4cout<<"------------------"<<G4endl;
345  //G4cout<<" "<<G4endl;
346 
347  G4ThreeVector newDirection =
349  newDirection.rotateUz(dir);
351 }
352 
353 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
354 
358 {
359  // super simple, very compact look-up for total cross section
360  if (E==0.) {return 1e30;}
361  //G4cout<<"------------------"<<G4endl;
362  //G4cout<<"ComputeCrossSectionPerAtom"<<G4endl;
363  //G4cout<< E<<" "<<G4Log(E)<<" "<< G4Exp(fTcs->GetValueForElement(Z,log(E)))
364  //<<" "<<fTcs->GetValueForElement(Z,log(E))<<G4endl;
365  //G4cout<<"------------------"<<G4endl;
366  //G4cout<<" "<<G4endl;
367 
368  //G4cout<< E<<" "<<G4Exp(fTcs->GetValueForElement(Z,log(E)))<<G4endl;
369 
370  return G4Exp(fTcs->GetValueForElement(Z,G4Log(E)));
371 }
372 
373 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
374 
376 {
377  G4ThreeVector dir(0.0,0.0,1.0);
378 
379  G4double sint = sin(acos(z1));
380  G4double cost = sqrt(1.0 - sint*sint);
381  G4double phi = twopi* G4UniformRand();
382  G4double dirx = sint*cos(phi);
383  G4double diry = sint*sin(phi);
384  G4double dirz = cost;
385 
386  dir.set(dirx,diry,dirz);
387  return dir;
388 }
389 
390 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
391 
392 
static const double MeV
Definition: G4SIunits.hh:211
G4PhysicsVector * GetComponentDataByIndex(G4int Z, size_t idx)
void ReadData(G4int, const G4String &argFileName)
void AddComponent(G4int Z, G4int id, G4PhysicsVector *v)
std::vector< G4Element * > G4ElementVector
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
G4ThreeVector GetNewDirection(G4double z1)
G4PhysicsVector * GetElementData(G4int Z)
static const G4double e2
G4bool IsMaster() const
Definition: G4VEmModel.hh:718
G4ParticleDefinition * GetDefinition() const
G4DiscreteScatteringModel(G4int iNumAngles=1)
size_t GetVectorLength() const
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:190
int G4int
Definition: G4Types.hh:78
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:725
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double E, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
void InitialiseForElement(G4int Z, G4PhysicsVector *v)
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetMomentumDirection() const
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double maxEnergy=DBL_MAX)
static const double twopi
Definition: G4SIunits.hh:75
size_t FindBin(G4double energy, size_t idx) const
G4double Energy(size_t index) const
static const G4double e1
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
static G4ProductionCutsTable * GetProductionCutsTable()
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
int G4lrint(double ad)
Definition: templates.hh:163
G4double GetValueForElement(G4int Z, G4double kinEnergy)
G4ParticleChangeForGamma * fParticleChange
virtual G4bool Retrieve(std::ifstream &fIn, G4bool ascii=false)
#define G4endl
Definition: G4ios.hh:61
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
static const double keV
Definition: G4SIunits.hh:213
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
void InitialiseForComponent(G4int Z, G4int nComponents=0)
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:544
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:134
const G4Material * GetMaterial() const