Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4BremsstrahlungParameters Class Reference

#include <G4BremsstrahlungParameters.hh>

Public Member Functions

 G4BremsstrahlungParameters (const G4String &name, size_t num, G4int minZ=1, G4int maxZ=99)
 
 ~G4BremsstrahlungParameters ()
 
G4double Parameter (G4int parameterIndex, G4int Z, G4double energy) const
 
G4double ParameterC (G4int index) const
 
void PrintData () const
 

Detailed Description

Definition at line 61 of file G4BremsstrahlungParameters.hh.

Constructor & Destructor Documentation

◆ G4BremsstrahlungParameters()

G4BremsstrahlungParameters::G4BremsstrahlungParameters ( const G4String name,
size_t  num,
G4int  minZ = 1,
G4int  maxZ = 99 
)

Definition at line 53 of file G4BremsstrahlungParameters.cc.

55 : zMin(minZ),
56 zMax(maxZ),
57 length(num)
58{
59 LoadData(name);
60}

◆ ~G4BremsstrahlungParameters()

G4BremsstrahlungParameters::~G4BremsstrahlungParameters ( )

Definition at line 63 of file G4BremsstrahlungParameters.cc.

64{
65 // Reset the map of data sets: remove the data sets from the map
66 std::map<G4int,G4VEMDataSet*,std::less<G4int> >::iterator pos;
67
68 for (pos = param.begin(); pos != param.end(); ++pos)
69 {
70 G4VEMDataSet* dataSet = (*pos).second;
71 delete dataSet;
72 }
73
74 activeZ.clear();
75 paramC.clear();
76}

Member Function Documentation

◆ Parameter()

G4double G4BremsstrahlungParameters::Parameter ( G4int  parameterIndex,
G4int  Z,
G4double  energy 
) const

Definition at line 79 of file G4BremsstrahlungParameters.cc.

82{
83 G4double value = 0.;
84 G4int id = Z*length + parameterIndex;
85 std::map<G4int,G4VEMDataSet*,std::less<G4int> >::const_iterator pos;
86
87 pos = param.find(id);
88 if (pos!= param.end()) {
89
90 G4VEMDataSet* dataSet = (*pos).second;
91 const G4DataVector ener = dataSet->GetEnergies(0);
92 G4double ee = std::max(ener.front(),std::min(ener.back(),energy));
93 value = dataSet->FindValue(ee);
94
95 } else {
96 G4cout << "WARNING: G4BremsstrahlungParameters::FindValue "
97 << "did not find ID = "
98 << id << G4endl;
99 }
100
101 return value;
102}
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
virtual const G4DataVector & GetEnergies(G4int componentId) const =0
virtual G4double FindValue(G4double x, G4int componentId=0) const =0

Referenced by G4eBremsstrahlungSpectrum::AverageEnergy(), G4eBremsstrahlungSpectrum::Probability(), and G4eBremsstrahlungSpectrum::SampleEnergy().

◆ ParameterC()

G4double G4BremsstrahlungParameters::ParameterC ( G4int  index) const

Definition at line 248 of file G4BremsstrahlungParameters.cc.

249{
250 G4int n = paramC.size();
251 if (id < 0 || id >= n)
252 {
253 G4String stringConversion2(id);
254 G4String ex = "Wrong id " + stringConversion2;
255 G4Exception("G4BremsstrahlungParameters::ParameterC",
256 "em1002",FatalException,ex);
257
258 }
259
260 return paramC[id];
261}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35

Referenced by G4eBremsstrahlungSpectrum::AverageEnergy().

◆ PrintData()

void G4BremsstrahlungParameters::PrintData ( ) const

Definition at line 264 of file G4BremsstrahlungParameters.cc.

265{
266
267 G4cout << G4endl;
268 G4cout << "===== G4BremsstrahlungParameters =====" << G4endl;
269 G4cout << G4endl;
270 G4cout << "===== Parameters =====" << G4endl;
271 G4cout << G4endl;
272
273 size_t nZ = activeZ.size();
274 std::map<G4int,G4VEMDataSet*,std::less<G4int> >::const_iterator pos;
275
276 for (size_t j=0; j<nZ; j++) {
277 G4int Z = (G4int)activeZ[j];
278
279 for (size_t i=0; i<length; i++) {
280
281 pos = param.find(Z*length + i);
282 if (pos!= param.end()) {
283
284 G4cout << "===== Z= " << Z
285 << " parameter[" << i << "] ====="
286 << G4endl;
287 G4VEMDataSet* dataSet = (*pos).second;
288 dataSet->PrintData();
289 }
290 }
291 }
292
293 G4cout << "==========================================" << G4endl;
294}
virtual void PrintData(void) const =0

Referenced by G4eBremsstrahlungSpectrum::PrintData().


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