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

#include <G4PenelopeCrossSection.hh>

Public Member Functions

 G4PenelopeCrossSection (size_t nOfEnergyPoints, size_t nOfShells=0)
 
 ~G4PenelopeCrossSection ()
 
G4double GetTotalCrossSection (G4double energy) const
 Returns total cross section at the given energy. More...
 
G4double GetHardCrossSection (G4double energy) const
 Returns hard cross section at the given energy. More...
 
G4double GetSoftStoppingPower (G4double energy) const
 Returns the total stopping power due to soft collisions. More...
 
G4double GetShellCrossSection (size_t shellID, G4double energy) const
 Returns the hard cross section for the given shell (per molecule) More...
 
G4double GetNormalizedShellCrossSection (size_t shellID, G4double energy) const
 Returns the hard cross section for the given shell (normalized to 1) More...
 
size_t GetNumberOfShells () const
 
void AddCrossSectionPoint (size_t binNumber, G4double energy, G4double XH0, G4double XH1, G4double XH2, G4double XS0, G4double XS1, G4double XS2)
 
void AddShellCrossSectionPoint (size_t binNumber, size_t shellID, G4double energy, G4double xs)
 
void NormalizeShellCrossSections ()
 

Detailed Description

Definition at line 72 of file G4PenelopeCrossSection.hh.

Constructor & Destructor Documentation

G4PenelopeCrossSection::G4PenelopeCrossSection ( size_t  nOfEnergyPoints,
size_t  nOfShells = 0 
)

Definition at line 44 of file G4PenelopeCrossSection.cc.

References FatalException, G4endl, G4Exception(), and G4PhysicsTable::push_back().

44  :
45  numberOfEnergyPoints(nPointsE),numberOfShells(nShells),softCrossSections(0),
46  hardCrossSections(0),shellCrossSections(0),shellNormalizedCrossSections(0)
47 {
48  //check the number of points is not zero
49  if (!numberOfEnergyPoints)
50  {
52  ed << "G4PenelopeCrossSection: invalid number of energy points " << G4endl;
53  G4Exception("G4PenelopeCrossSection::G4PenelopeCrossSection()",
54  "em2017",FatalException,ed);
55  }
56 
57  isNormalized = false;
58 
59  // 1) soft XS table
60  softCrossSections = new G4PhysicsTable();
61  //the table contains 3 G4PhysicsFreeVectors,
62  //(softCrossSections)[0] --> log XS0 vs. log E
63  //(softCrossSections)[1] --> log XS1 vs. log E
64  //(softCrossSections)[2] --> log XS2 vs. log E
65 
66  //everything is log-log
67  for (size_t i=0;i<3;i++)
68  softCrossSections->push_back(new G4PhysicsFreeVector(numberOfEnergyPoints));
69 
70  //2) hard XS table
71  hardCrossSections = new G4PhysicsTable();
72  //the table contains 3 G4PhysicsFreeVectors,
73  //(hardCrossSections)[0] --> log XH0 vs. log E
74  //(hardCrossSections)[1] --> log XH1 vs. log E
75  //(hardCrossSections)[2] --> log XH2 vs. log E
76 
77  //everything is log-log
78  for (size_t i=0;i<3;i++)
79  hardCrossSections->push_back(new G4PhysicsFreeVector(numberOfEnergyPoints));
80 
81  //3) shell XS table, if it is the case
82  if (numberOfShells)
83  {
84  shellCrossSections = new G4PhysicsTable();
85  shellNormalizedCrossSections = new G4PhysicsTable();
86  //the table has to contain numberofShells G4PhysicsFreeVectors,
87  //(theTable)[ishell] --> cross section for shell #ishell
88  for (size_t i=0;i<numberOfShells;i++)
89  {
90  shellCrossSections->push_back(new G4PhysicsFreeVector(numberOfEnergyPoints));
91  shellNormalizedCrossSections->push_back(new G4PhysicsFreeVector(numberOfEnergyPoints));
92  }
93  }
94 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
void push_back(G4PhysicsVector *)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61
G4PenelopeCrossSection::~G4PenelopeCrossSection ( )

Definition at line 97 of file G4PenelopeCrossSection.cc.

98 {
99  //clean up tables
100  if (shellCrossSections)
101  {
102  //shellCrossSections->clearAndDestroy();
103  delete shellCrossSections;
104  }
105  if (shellNormalizedCrossSections)
106  {
107  //shellNormalizedCrossSections->clearAndDestroy();
108  delete shellNormalizedCrossSections;
109  }
110  if (softCrossSections)
111  {
112  //softCrossSections->clearAndDestroy();
113  delete softCrossSections;
114  }
115  if (hardCrossSections)
116  {
117  //hardCrossSections->clearAndDestroy();
118  delete hardCrossSections;
119  }
120 }

Member Function Documentation

void G4PenelopeCrossSection::AddCrossSectionPoint ( size_t  binNumber,
G4double  energy,
G4double  XH0,
G4double  XH1,
G4double  XH2,
G4double  XS0,
G4double  XS1,
G4double  XS2 
)

Public interface for the master thread

Definition at line 123 of file G4PenelopeCrossSection.cc.

References python.hepunit::cm2, python.hepunit::eV, G4cout, G4endl, G4INCL::Math::max(), and G4PhysicsFreeVector::PutValue().

Referenced by G4PenelopeIonisationXSHandler::BuildXSTable().

128 {
129  if (!softCrossSections || !hardCrossSections)
130  {
131  G4cout << "Something wrong in G4PenelopeCrossSection::AddCrossSectionPoint" <<
132  G4endl;
133  G4cout << "Trying to fill un-initialized tables" << G4endl;
134  return;
135  }
136 
137  //fill vectors
138  G4PhysicsFreeVector* theVector = (G4PhysicsFreeVector*) (*softCrossSections)[0];
139 
140  if (binNumber >= numberOfEnergyPoints)
141  {
142  G4cout << "Something wrong in G4PenelopeCrossSection::AddCrossSectionPoint" <<
143  G4endl;
144  G4cout << "Trying to register more points than originally declared" << G4endl;
145  return;
146  }
147  G4double logEne = std::log(energy);
148 
149  //XS0
150  G4double val = std::log(std::max(XS0,1e-42*cm2)); //avoid log(0)
151  theVector->PutValue(binNumber,logEne,val);
152 
153  //XS1
154  theVector = (G4PhysicsFreeVector*) (*softCrossSections)[1];
155  val = std::log(std::max(XS1,1e-42*eV*cm2)); //avoid log(0)
156  theVector->PutValue(binNumber,logEne,val);
157 
158  //XS2
159  theVector = (G4PhysicsFreeVector*) (*softCrossSections)[2];
160  val = std::log(std::max(XS2,1e-42*eV*eV*cm2)); //avoid log(0)
161  theVector->PutValue(binNumber,logEne,val);
162 
163  //XH0
164  theVector = (G4PhysicsFreeVector*) (*hardCrossSections)[0];
165  val = std::log(std::max(XH0,1e-42*cm2)); //avoid log(0)
166  theVector->PutValue(binNumber,logEne,val);
167 
168  //XH1
169  theVector = (G4PhysicsFreeVector*) (*hardCrossSections)[1];
170  val = std::log(std::max(XH1,1e-42*eV*cm2)); //avoid log(0)
171  theVector->PutValue(binNumber,logEne,val);
172 
173  //XH2
174  theVector = (G4PhysicsFreeVector*) (*hardCrossSections)[2];
175  val = std::log(std::max(XH2,1e-42*eV*eV*cm2)); //avoid log(0)
176  theVector->PutValue(binNumber,logEne,val);
177 
178  return;
179 }
void PutValue(size_t binNumber, G4double binValue, G4double dataValue)
double precision function energy(A, Z)
Definition: dpm25nuc6.f:4106
G4GLOB_DLL std::ostream G4cout
T max(const T t1, const T t2)
brief Return the largest of the two arguments
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void G4PenelopeCrossSection::AddShellCrossSectionPoint ( size_t  binNumber,
size_t  shellID,
G4double  energy,
G4double  xs 
)

Definition at line 183 of file G4PenelopeCrossSection.cc.

References python.hepunit::cm2, G4cout, G4endl, G4INCL::Math::max(), and G4PhysicsFreeVector::PutValue().

Referenced by G4PenelopeIonisationXSHandler::BuildXSTable().

187 {
188  if (!shellCrossSections)
189  {
190  G4cout << "Something wrong in G4PenelopeCrossSection::AddShellCrossSectionPoint" <<
191  G4endl;
192  G4cout << "Trying to fill un-initialized table" << G4endl;
193  return;
194  }
195 
196  if (shellID >= numberOfShells)
197  {
198  G4cout << "Something wrong in G4PenelopeCrossSection::AddShellCrossSectionPoint" <<
199  G4endl;
200  G4cout << "Trying to fill shell #" << shellID << " while the maximum is "
201  << numberOfShells-1 << G4endl;
202  return;
203  }
204 
205  //fill vector
206  G4PhysicsFreeVector* theVector = (G4PhysicsFreeVector*) (*shellCrossSections)[shellID];
207 
208  if (binNumber >= numberOfEnergyPoints)
209  {
210  G4cout << "Something wrong in G4PenelopeCrossSection::AddShellCrossSectionPoint" <<
211  G4endl;
212  G4cout << "Trying to register more points than originally declared" << G4endl;
213  return;
214  }
215  G4double logEne = std::log(energy);
216  G4double val = std::log(std::max(xs,1e-42*cm2)); //avoid log(0)
217  theVector->PutValue(binNumber,logEne,val);
218 
219  return;
220 }
void PutValue(size_t binNumber, G4double binValue, G4double dataValue)
double precision function energy(A, Z)
Definition: dpm25nuc6.f:4106
G4GLOB_DLL std::ostream G4cout
T max(const T t1, const T t2)
brief Return the largest of the two arguments
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4double G4PenelopeCrossSection::GetHardCrossSection ( G4double  energy) const

Returns hard cross section at the given energy.

Definition at line 268 of file G4PenelopeCrossSection.cc.

References G4cout, G4endl, G4PhysicsVector::GetVectorLength(), and G4PhysicsVector::Value().

Referenced by G4PenelopeBremsstrahlungModel::CrossSectionPerVolume(), and G4PenelopeIonisationModel::CrossSectionPerVolume().

269 {
270  G4double result = 0;
271  //take here XH0
272  if (!hardCrossSections)
273  {
274  G4cout << "Something wrong in G4PenelopeCrossSection::GetHardCrossSection" <<
275  G4endl;
276  G4cout << "Trying to retrieve from un-initialized tables" << G4endl;
277  return result;
278  }
279 
280  G4PhysicsFreeVector* theVector = (G4PhysicsFreeVector*) (*hardCrossSections)[0];
281  if (theVector->GetVectorLength() < numberOfEnergyPoints)
282  {
283  G4cout << "Something wrong in G4PenelopeCrossSection::GetHardCrossSection" <<
284  G4endl;
285  G4cout << "Hard cross section table looks not filled" << G4endl;
286  return result;
287  }
288  G4double logene = std::log(energy);
289  G4double logXS = theVector->Value(logene);
290  result = std::exp(logXS);
291 
292  return result;
293 }
size_t GetVectorLength() const
double precision function energy(A, Z)
Definition: dpm25nuc6.f:4106
G4GLOB_DLL std::ostream G4cout
G4double Value(G4double theEnergy, size_t &lastidx) const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4double G4PenelopeCrossSection::GetNormalizedShellCrossSection ( size_t  shellID,
G4double  energy 
) const

Returns the hard cross section for the given shell (normalized to 1)

Definition at line 363 of file G4PenelopeCrossSection.cc.

References G4cout, G4endl, G4PhysicsVector::GetVectorLength(), and G4PhysicsVector::Value().

364 {
365  G4double result = 0;
366  if (!shellNormalizedCrossSections)
367  {
368  G4cout << "Something wrong in G4PenelopeCrossSection::GetShellCrossSection" <<
369  G4endl;
370  G4cout << "Trying to retrieve from un-initialized tables" << G4endl;
371  return result;
372  }
373 
374  if (!isNormalized)
375  {
376  G4cout << "Something wrong in G4PenelopeCrossSection::GetShellCrossSection" << G4endl;
377  G4cout << "The table of normalized cross section is not initialized" << G4endl;
378  }
379 
380 
381  if (shellID >= numberOfShells)
382  {
383  G4cout << "Something wrong in G4PenelopeCrossSection::GetShellCrossSection" <<
384  G4endl;
385  G4cout << "Trying to retrieve shell #" << shellID << " while the maximum is "
386  << numberOfShells-1 << G4endl;
387  return result;
388  }
389 
390  const G4PhysicsFreeVector* theVector =
391  (G4PhysicsFreeVector*) (*shellNormalizedCrossSections)[shellID];
392 
393  if (theVector->GetVectorLength() < numberOfEnergyPoints)
394  {
395  G4cout << "Something wrong in G4PenelopeCrossSection::GetShellCrossSection" <<
396  G4endl;
397  G4cout << "Shell cross section table looks not filled" << G4endl;
398  return result;
399  }
400  G4double logene = std::log(energy);
401  G4double logXS = theVector->Value(logene);
402  result = std::exp(logXS);
403 
404  return result;
405 }
size_t GetVectorLength() const
double precision function energy(A, Z)
Definition: dpm25nuc6.f:4106
G4GLOB_DLL std::ostream G4cout
G4double Value(G4double theEnergy, size_t &lastidx) const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
size_t G4PenelopeCrossSection::GetNumberOfShells ( ) const
inline

Definition at line 94 of file G4PenelopeCrossSection.hh.

94 {return numberOfShells;};
G4double G4PenelopeCrossSection::GetShellCrossSection ( size_t  shellID,
G4double  energy 
) const

Returns the hard cross section for the given shell (per molecule)

Definition at line 327 of file G4PenelopeCrossSection.cc.

References G4cout, G4endl, G4PhysicsVector::GetVectorLength(), and G4PhysicsVector::Value().

Referenced by G4PenelopeIonisationCrossSection::CrossSection().

328 {
329  G4double result = 0;
330  if (!shellCrossSections)
331  {
332  G4cout << "Something wrong in G4PenelopeCrossSection::GetShellCrossSection" <<
333  G4endl;
334  G4cout << "Trying to retrieve from un-initialized tables" << G4endl;
335  return result;
336  }
337  if (shellID >= numberOfShells)
338  {
339  G4cout << "Something wrong in G4PenelopeCrossSection::GetShellCrossSection" <<
340  G4endl;
341  G4cout << "Trying to retrieve shell #" << shellID << " while the maximum is "
342  << numberOfShells-1 << G4endl;
343  return result;
344  }
345 
346  G4PhysicsFreeVector* theVector = (G4PhysicsFreeVector*) (*shellCrossSections)[shellID];
347 
348  if (theVector->GetVectorLength() < numberOfEnergyPoints)
349  {
350  G4cout << "Something wrong in G4PenelopeCrossSection::GetShellCrossSection" <<
351  G4endl;
352  G4cout << "Shell cross section table looks not filled" << G4endl;
353  return result;
354  }
355  G4double logene = std::log(energy);
356  G4double logXS = theVector->Value(logene);
357  result = std::exp(logXS);
358 
359  return result;
360 }
size_t GetVectorLength() const
double precision function energy(A, Z)
Definition: dpm25nuc6.f:4106
G4GLOB_DLL std::ostream G4cout
G4double Value(G4double theEnergy, size_t &lastidx) const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4double G4PenelopeCrossSection::GetSoftStoppingPower ( G4double  energy) const

Returns the total stopping power due to soft collisions.

Definition at line 298 of file G4PenelopeCrossSection.cc.

References G4cout, G4endl, G4PhysicsVector::GetVectorLength(), and G4PhysicsVector::Value().

Referenced by G4PenelopeBremsstrahlungModel::ComputeDEDXPerVolume(), and G4PenelopeIonisationModel::ComputeDEDXPerVolume().

299 {
300  G4double result = 0;
301  //take here XH0
302  if (!softCrossSections)
303  {
304  G4cout << "Something wrong in G4PenelopeCrossSection::GetSoftStoppingPower" <<
305  G4endl;
306  G4cout << "Trying to retrieve from un-initialized tables" << G4endl;
307  return result;
308  }
309 
310  G4PhysicsFreeVector* theVector = (G4PhysicsFreeVector*) (*softCrossSections)[1];
311  if (theVector->GetVectorLength() < numberOfEnergyPoints)
312  {
313  G4cout << "Something wrong in G4PenelopeCrossSection::GetSoftStoppingPower" <<
314  G4endl;
315  G4cout << "Soft cross section table looks not filled" << G4endl;
316  return result;
317  }
318  G4double logene = std::log(energy);
319  G4double logXS = theVector->Value(logene);
320  result = std::exp(logXS);
321 
322  return result;
323 }
size_t GetVectorLength() const
double precision function energy(A, Z)
Definition: dpm25nuc6.f:4106
G4GLOB_DLL std::ostream G4cout
G4double Value(G4double theEnergy, size_t &lastidx) const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4double G4PenelopeCrossSection::GetTotalCrossSection ( G4double  energy) const

Returns total cross section at the given energy.

Definition at line 224 of file G4PenelopeCrossSection.cc.

References G4cout, G4endl, G4PhysicsVector::GetVectorLength(), and G4PhysicsVector::Value().

Referenced by G4PenelopeIonisationModel::CrossSectionPerVolume().

225 {
226  G4double result = 0;
227  //take here XS0 + XH0
228  if (!softCrossSections || !hardCrossSections)
229  {
230  G4cout << "Something wrong in G4PenelopeCrossSection::GetTotalCrossSection" <<
231  G4endl;
232  G4cout << "Trying to retrieve from un-initialized tables" << G4endl;
233  return result;
234  }
235 
236  // 1) soft part
237  G4PhysicsFreeVector* theVector = (G4PhysicsFreeVector*) (*softCrossSections)[0];
238  if (theVector->GetVectorLength() < numberOfEnergyPoints)
239  {
240  G4cout << "Something wrong in G4PenelopeCrossSection::GetTotalCrossSection" <<
241  G4endl;
242  G4cout << "Soft cross section table looks not filled" << G4endl;
243  return result;
244  }
245  G4double logene = std::log(energy);
246  G4double logXS = theVector->Value(logene);
247  G4double softXS = std::exp(logXS);
248 
249  // 2) hard part
250  theVector = (G4PhysicsFreeVector*) (*hardCrossSections)[0];
251  if (theVector->GetVectorLength() < numberOfEnergyPoints)
252  {
253  G4cout << "Something wrong in G4PenelopeCrossSection::GetTotalCrossSection" <<
254  G4endl;
255  G4cout << "Hard cross section table looks not filled" << G4endl;
256  return result;
257  }
258  logXS = theVector->Value(logene);
259  G4double hardXS = std::exp(logXS);
260 
261  result = hardXS + softXS;
262  return result;
263 
264 }
size_t GetVectorLength() const
double precision function energy(A, Z)
Definition: dpm25nuc6.f:4106
G4GLOB_DLL std::ostream G4cout
G4double Value(G4double theEnergy, size_t &lastidx) const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void G4PenelopeCrossSection::NormalizeShellCrossSections ( )

Definition at line 411 of file G4PenelopeCrossSection.cc.

References G4cout, G4endl, G4PhysicsVector::GetLowEdgeEnergy(), and G4PhysicsFreeVector::PutValue().

Referenced by G4PenelopeIonisationXSHandler::BuildXSTable().

412 {
413  if (isNormalized) //already done!
414  {
415  G4cout << "G4PenelopeCrossSection::NormalizeShellCrossSections()" << G4endl;
416  G4cout << "already invoked. Ignore it" << G4endl;
417  return;
418  }
419 
420  if (!shellNormalizedCrossSections)
421  {
422  G4cout << "Something wrong in G4PenelopeCrossSection::GetShellCrossSection" <<
423  G4endl;
424  G4cout << "Trying to retrieve from un-initialized tables" << G4endl;
425  return;
426  }
427 
428  for (size_t i=0;i<numberOfEnergyPoints;i++) //loop on energy
429  {
430  //energy grid is the same for all shells
431 
432  //Recalculate manually the XS factor, to avoid problems with
433  //underflows
434  G4double normFactor = 0.;
435  for (size_t shellID=0;shellID<numberOfShells;shellID++)
436  {
437  G4PhysicsFreeVector* theVec =
438  (G4PhysicsFreeVector*) (*shellCrossSections)[shellID];
439 
440  normFactor += std::exp((*theVec)[i]);
441  }
442  G4double logNormFactor = std::log(normFactor);
443  //Normalize
444  for (size_t shellID=0;shellID<numberOfShells;shellID++)
445  {
446  G4PhysicsFreeVector* theVec =
447  (G4PhysicsFreeVector*) (*shellNormalizedCrossSections)[shellID];
448  G4PhysicsFreeVector* theFullVec =
449  (G4PhysicsFreeVector*) (*shellCrossSections)[shellID];
450  G4double previousValue = (*theFullVec)[i]; //log(XS)
451  G4double logEnergy = theFullVec->GetLowEdgeEnergy(i);
452  //log(XS/normFactor) = log(XS) - log(normFactor)
453  theVec->PutValue(i,logEnergy,previousValue-logNormFactor);
454  }
455  }
456 
457  isNormalized = true;
458 
459 
460  /*
461  //TESTING
462  for (size_t shellID=0;shellID<numberOfShells;shellID++)
463  {
464  G4cout << "SHELL " << shellID << G4endl;
465  G4PhysicsFreeVector* theVec =
466  (G4PhysicsFreeVector*) (*shellCrossSections)[shellID];
467  for (size_t i=0;i<numberOfEnergyPoints;i++) //loop on energy
468  {
469  G4double logene = theVec->GetLowEdgeEnergy(i);
470  G4cout << std::exp(logene)/MeV << " " << std::exp((*theVec)[i]) << G4endl;
471  }
472  }
473  */
474 
475  return;
476 }
void PutValue(size_t binNumber, G4double binValue, G4double dataValue)
G4double GetLowEdgeEnergy(size_t binNumber) const
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76

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