Geant4-11
Public Member Functions | Private Member Functions | Private Attributes | Static Private Attributes
G4ICRU90StoppingData Class Reference

#include <G4ICRU90StoppingData.hh>

Public Member Functions

 G4ICRU90StoppingData ()
 
 G4ICRU90StoppingData (const G4ICRU90StoppingData &)=delete
 
G4double GetElectronicDEDXforAlpha (const G4Material *, G4double scaledKinEnergy) const
 
G4double GetElectronicDEDXforAlpha (G4int idx, G4double scaledKinEnergy) const
 
G4double GetElectronicDEDXforProton (const G4Material *, G4double kinEnergy) const
 
G4double GetElectronicDEDXforProton (G4int idx, G4double kinEnergy) const
 
G4int GetIndex (const G4Material *) const
 
G4int GetIndex (const G4String &) const
 
void Initialise ()
 
G4bool IsApplicable (const G4Material *) const
 
G4ICRU90StoppingDataoperator= (const G4ICRU90StoppingData &right)=delete
 
 ~G4ICRU90StoppingData ()
 

Private Member Functions

G4PhysicsFreeVectorAddData (G4int n, const G4double *e, const G4float *dedx)
 
void FillData ()
 
G4double GetDEDX (G4PhysicsFreeVector *, G4double e) const
 

Private Attributes

G4bool isInitialized
 
const G4Materialmaterials [nvectors]
 
G4PhysicsFreeVectorsdata_alpha [nvectors]
 
G4PhysicsFreeVectorsdata_proton [nvectors]
 

Static Private Attributes

static constexpr G4int nvectors = 3
 

Detailed Description

Definition at line 56 of file G4ICRU90StoppingData.hh.

Constructor & Destructor Documentation

◆ G4ICRU90StoppingData() [1/2]

G4ICRU90StoppingData::G4ICRU90StoppingData ( )
explicit

Definition at line 53 of file G4ICRU90StoppingData.cc.

53 : isInitialized(false)
54{
55 // 1st initialisation
56 for(size_t i=0; i<nvectors; ++i) {
57 materials[i] = nullptr;
58 sdata_proton[i] = nullptr;
59 sdata_alpha[i] = nullptr;
60 }
61 FillData();
62
63 Initialise();
64}
G4PhysicsFreeVector * sdata_proton[nvectors]
static constexpr G4int nvectors
const G4Material * materials[nvectors]
G4PhysicsFreeVector * sdata_alpha[nvectors]

References FillData(), Initialise(), materials, nvectors, sdata_alpha, and sdata_proton.

◆ ~G4ICRU90StoppingData()

G4ICRU90StoppingData::~G4ICRU90StoppingData ( )

Definition at line 68 of file G4ICRU90StoppingData.cc.

69{
70 for(size_t i=0; i<nvectors; ++i) {
71 delete sdata_proton[i];
72 delete sdata_alpha[i];
73 }
74}

References nvectors, sdata_alpha, and sdata_proton.

◆ G4ICRU90StoppingData() [2/2]

G4ICRU90StoppingData::G4ICRU90StoppingData ( const G4ICRU90StoppingData )
delete

Member Function Documentation

◆ AddData()

G4PhysicsFreeVector * G4ICRU90StoppingData::AddData ( G4int  n,
const G4double e,
const G4float dedx 
)
private

Definition at line 163 of file G4ICRU90StoppingData.cc.

165{
167
168 G4PhysicsFreeVector* data = new G4PhysicsFreeVector(n, e[0], e[n-1], true);
169 for(G4int i=0; i<n; ++i) {
170 data->PutValues(i, e[i]*CLHEP::MeV, ((G4double)dedx[i])*fac);
171 }
172 data->FillSecondDerivatives();
173 return data;
174}
static const G4double fac
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
void PutValues(const std::size_t index, const G4double energy, const G4double value)
void FillSecondDerivatives(const G4SplineType=G4SplineType::Base, const G4double dir1=0.0, const G4double dir2=0.0)
static constexpr double MeV
static constexpr double g
static constexpr double cm2

References CLHEP::cm2, fac, G4PhysicsVector::FillSecondDerivatives(), CLHEP::g, CLHEP::MeV, CLHEP::detail::n, and G4PhysicsFreeVector::PutValues().

Referenced by FillData().

◆ FillData()

void G4ICRU90StoppingData::FillData ( )
private

Definition at line 134 of file G4ICRU90StoppingData.cc.

135{
136 G4double T0_proton[57] = { 0.0010, 0.001500, 0.0020, 0.0030, 0.0040, 0.0050, 0.0060, 0.0080, 0.010, 0.0150, 0.020, 0.030, 0.040, 0.050, 0.060, 0.080, 0.10, 0.150, 0.20, 0.30, 0.40, 0.50, 0.60, 0.80, 1.00, 1.50, 2.00, 3.00, 4.00, 5.00, 6.00, 8.00, 10.00, 15.00, 20.00, 30.00, 40.00, 50.00, 60.00, 80.00, 100.00, 150.00, 200.00, 300.00, 400.00, 500.00, 600.00, 800.00, 1000.00, 1500.00, 2000.00, 3000.00, 4000.00, 5000.00, 6000.00, 8000.00, 10000.00 };
137
138 G4double T0_alpha[49] = { 0.0010, 0.001500, 0.0020, 0.0030, 0.0040, 0.0050, 0.0060, 0.0080, 0.010, 0.0150, 0.020, 0.030, 0.040, 0.050, 0.060, 0.080, 0.10, 0.150, 0.20, 0.30, 0.40, 0.50, 0.60, 0.80, 1.00, 1.50, 2.00, 3.00, 4.00, 5.00, 6.00, 8.00, 10.00, 15.00, 20.00, 30.00, 40.00, 50.00, 60.00, 80.00, 100.00, 150.00, 200.00, 300.00, 400.00, 500.00, 600.00, 800.00, 1000.00};
139
140 static const G4float e0_proton[57] = { 119.70f, 146.70f, 169.30f, 207.40f, 239.50f, 267.80f, 293.30f, 338.70f, 378.70f, 450.40f, 506.70f, 590.50f, 648.30f, 687.70f, 713.20f, 734.10f, 729.00f, 667.20f, 592.20f, 476.50f, 401.20f, 349.80f, 312.10f, 258.70f, 222.70f, 168.20f, 137.00f, 101.70f, 81.920f, 69.050f, 59.940f, 47.810f, 40.040f, 28.920f, 22.930f, 16.520f, 13.120f, 10.980f, 9.5140f, 7.6180f, 6.4410f, 4.8150f, 3.9750f, 3.1170f, 2.6860f, 2.4310f, 2.2650f, 2.0690f, 1.9620f, 1.850f, 1.820f, 1.8280f, 1.8610f, 1.8980f, 1.9340f, 1.9980f, 2.0520f };
141
142 static const G4float e0_alpha[49] = { 87.50f, 108.60f, 126.70f, 157.30f, 183.50f, 206.70f, 227.90f, 265.90f, 299.60f, 372.30f, 434.30f, 539.50f, 629.00f, 708.10f, 779.80f, 906.80f, 1018.00f, 1247.00f, 1429.00f, 1693.00f, 1861.00f, 1961.00f, 2008.00f, 2002.00f, 1922.00f, 1626.00f, 1381.00f, 1071.00f, 885.80f, 760.60f, 669.60f, 545.30f, 463.40f, 342.30f, 274.70f, 200.10f, 159.30f, 133.20f, 115.00f, 91.190f, 76.140f, 54.930f, 43.690f, 31.840f, 25.630f, 21.790f, 19.170f, 15.830f, 13.790f };
143
144 static const G4float e1_proton[57] = { 133.70f, 163.80f, 189.10f, 231.60f, 267.50f, 299.00f, 327.60f, 378.20f, 422.90f, 503.60f, 567.30f, 662.80f, 729.00f, 774.00f, 802.60f, 824.10f, 814.50f, 736.00f, 658.50f, 543.50f, 464.30f, 406.50f, 362.40f, 299.70f, 257.40f, 193.40f, 156.90f, 116.00f, 93.190f, 78.420f, 68.010f, 54.170f, 45.320f, 32.690f, 25.890f, 18.640f, 14.790f, 12.380f, 10.720f, 8.5780f, 7.250f, 5.4170f, 4.470f, 3.5040f, 3.0180f, 2.7310f, 2.5440f, 2.3230f, 2.2030f, 2.0650f, 2.0170f, 1.9980f, 2.010f, 2.0290f, 2.050f, 2.090f, 2.1240f };
145
146 static const G4float e1_alpha[49] = { 98.910f, 122.80f, 143.10f, 177.60f, 206.90f, 233.00f, 256.80f, 299.30f, 337.00f, 418.10f, 487.10f, 603.90f, 703.00f, 790.50f, 869.60f, 1009.00f, 1131.00f, 1383.00f, 1582.00f, 1873.00f, 2062.00f, 2178.00f, 2240.00f, 2256.00f, 2190.00f, 1877.00f, 1599.00f, 1239.00f, 1021.00f, 874.70f, 768.70f, 623.90f, 529.00f, 389.40f, 311.80f, 226.70f, 180.20f, 150.60f, 130.00f, 103.00f, 85.930f, 61.940f, 49.230f, 35.860f, 28.850f, 24.520f, 21.560f, 17.80f, 15.50f };
147
148 static const G4float e2_proton[57] = { 118.50f, 145.10f, 167.60f, 205.30f, 237.00f, 265.00f, 290.30f, 335.20f, 374.80f, 435.30f, 481.10f, 550.30f, 611.00f, 663.00f, 697.60f, 726.70f, 726.60f, 671.20f, 598.60f, 483.30f, 407.10f, 354.60f, 315.90f, 262.40f, 226.10f, 171.10f, 139.40f, 103.50f, 83.260f, 70.10f, 60.810f, 48.450f, 40.550f, 29.250f, 23.170f, 16.680f, 13.230f, 11.070f, 9.5910f, 7.6750f, 6.4860f, 4.8430f, 3.9940f, 3.1250f, 2.6870f, 2.4260f, 2.2550f, 2.050f, 1.9360f, 1.8060f, 1.760f, 1.7440f, 1.7560f, 1.7750f, 1.7960f, 1.8340f, 1.8680f };
149
150 static const G4float e2_alpha[49] = { 192.30f, 228.90f, 259.00f, 308.30f, 348.90f, 384.00f, 415.30f, 469.90f, 517.00f, 615.10f, 695.50f, 826.20f, 932.70f, 1024.00f, 1104.00f, 1240.00f, 1354.00f, 1574.00f, 1731.00f, 1929.00f, 2027.00f, 2063.00f, 2060.00f, 1993.00f, 1891.00f, 1615.00f, 1387.00f, 1085.00f, 898.20f, 772.00f, 680.40f, 554.40f, 471.20f, 347.80f, 278.80f, 202.80f, 161.20f, 134.80f, 116.30f, 92.140f, 76.890f, 55.430f, 44.050f, 32.090f, 25.810f, 21.930f, 19.280f, 15.910f, 13.840f };
151
152 sdata_proton[0] = AddData(57, T0_proton, e0_proton);
153 sdata_proton[1] = AddData(57, T0_proton, e1_proton);
154 sdata_proton[2] = AddData(57, T0_proton, e2_proton);
155
156 sdata_alpha[0] = AddData(49, T0_alpha, e0_alpha);
157 sdata_alpha[1] = AddData(49, T0_alpha, e1_alpha);
158 sdata_alpha[2] = AddData(49, T0_alpha, e2_alpha);
159}
float G4float
Definition: G4Types.hh:84
G4PhysicsFreeVector * AddData(G4int n, const G4double *e, const G4float *dedx)

References AddData(), sdata_alpha, and sdata_proton.

Referenced by G4ICRU90StoppingData().

◆ GetDEDX()

G4double G4ICRU90StoppingData::GetDEDX ( G4PhysicsFreeVector data,
G4double  e 
) const
inlineprivate

Definition at line 134 of file G4ICRU90StoppingData.hh.

135{
136 G4double emin = data->Energy(0);
137 return (e <= emin) ? (*data)[0]*std::sqrt(e/emin) : data->Value(e);
138}
G4double Energy(const std::size_t index) const

References G4PhysicsVector::Energy(), and G4PhysicsVector::Value().

Referenced by GetElectronicDEDXforAlpha(), and GetElectronicDEDXforProton().

◆ GetElectronicDEDXforAlpha() [1/2]

G4double G4ICRU90StoppingData::GetElectronicDEDXforAlpha ( const G4Material mat,
G4double  scaledKinEnergy 
) const

Definition at line 125 of file G4ICRU90StoppingData.cc.

127{
128 G4int idx = GetIndex(mat);
129 return (idx < 0) ? 0.0 : GetDEDX(sdata_alpha[idx], scaledKinEnergy);
130}
G4double GetDEDX(G4PhysicsFreeVector *, G4double e) const
G4int GetIndex(const G4Material *) const

References GetDEDX(), GetIndex(), and sdata_alpha.

Referenced by G4BetheBlochModel::ComputeDEDXPerVolume(), and G4BraggIonModel::DEDX().

◆ GetElectronicDEDXforAlpha() [2/2]

G4double G4ICRU90StoppingData::GetElectronicDEDXforAlpha ( G4int  idx,
G4double  scaledKinEnergy 
) const
inline

Definition at line 151 of file G4ICRU90StoppingData.hh.

153{
154 return (idx < 0 || idx >= nvectors) ? 0.0
155 : GetDEDX(sdata_alpha[idx], scaledKinEnergy);
156}

References GetDEDX(), nvectors, and sdata_alpha.

◆ GetElectronicDEDXforProton() [1/2]

G4double G4ICRU90StoppingData::GetElectronicDEDXforProton ( const G4Material mat,
G4double  kinEnergy 
) const

Definition at line 116 of file G4ICRU90StoppingData.cc.

118{
119 G4int idx = GetIndex(mat);
120 return (idx < 0) ? 0.0 : GetDEDX(sdata_proton[idx], kinEnergy);
121}

References GetDEDX(), GetIndex(), and sdata_proton.

Referenced by G4BetheBlochModel::ComputeDEDXPerVolume(), and G4BraggModel::DEDX().

◆ GetElectronicDEDXforProton() [2/2]

G4double G4ICRU90StoppingData::GetElectronicDEDXforProton ( G4int  idx,
G4double  kinEnergy 
) const
inline

Definition at line 142 of file G4ICRU90StoppingData.hh.

144{
145 return (idx < 0 || idx >= nvectors) ? 0.0
146 : GetDEDX(sdata_proton[idx], kinEnergy);
147}

References GetDEDX(), nvectors, and sdata_proton.

◆ GetIndex() [1/2]

G4int G4ICRU90StoppingData::GetIndex ( const G4Material mat) const
inline

Definition at line 111 of file G4ICRU90StoppingData.hh.

112{
113 G4int idx = -1;
114 if (mat == materials[0]) { idx = 0; }
115 else if(mat == materials[1]) { idx = 1; }
116 else if(mat == materials[2]) { idx = 2; }
117 return idx;
118}

References materials.

Referenced by G4BetheBlochModel::ComputeDEDXPerVolume(), G4BraggIonModel::DEDX(), G4BraggModel::DEDX(), GetElectronicDEDXforAlpha(), and GetElectronicDEDXforProton().

◆ GetIndex() [2/2]

G4int G4ICRU90StoppingData::GetIndex ( const G4String nam) const
inline

Definition at line 122 of file G4ICRU90StoppingData.hh.

123{
124 G4int idx = -1;
125 if (nam == materials[0]->GetName()) { idx = 0; }
126 else if(nam == materials[1]->GetName()) { idx = 1; }
127 else if(nam == materials[2]->GetName()) { idx = 2; }
128 return idx;
129}

References materials.

◆ Initialise()

void G4ICRU90StoppingData::Initialise ( )

Definition at line 78 of file G4ICRU90StoppingData.cc.

79{
80 if(isInitialized) { return; }
81 // this method may be called several times during initialisation
83 if(nmat == (G4int)nvectors) { return; }
84
85 static const G4String nameNIST_ICRU90[3] =
86 {"G4_AIR","G4_WATER","G4_GRAPHITE"};
87
88 // loop via material list to add extra data
89 for(G4int i=0; i<nmat; ++i) {
90 const G4Material* mat = (*(G4Material::GetMaterialTable()))[i];
91
92 G4bool isThere = false;
93 for(G4int j=0; j<nvectors; ++j) {
94 if(mat == materials[j]) {
95 isThere = true;
96 break;
97 }
98 }
99 if(!isThere) {
100 // check list of NIST materials
101 G4String mname = mat->GetName();
102 for(G4int j=0; j<nvectors; ++j) {
103 if(mname == nameNIST_ICRU90[j]) {
104 materials[j] = mat;
105 break;
106 }
107 }
108 }
109 isInitialized = (materials[0] && materials[1] && materials[2]);
110 if(isInitialized) { return; }
111 }
112}
bool G4bool
Definition: G4Types.hh:86
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:679
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:672
const G4String & GetName() const
Definition: G4Material.hh:173

References G4Material::GetMaterialTable(), G4Material::GetName(), G4Material::GetNumberOfMaterials(), isInitialized, materials, and nvectors.

Referenced by G4ICRU90StoppingData(), G4BetheBlochModel::Initialise(), G4BraggIonModel::Initialise(), and G4BraggModel::Initialise().

◆ IsApplicable()

G4bool G4ICRU90StoppingData::IsApplicable ( const G4Material mat) const
inline

Definition at line 106 of file G4ICRU90StoppingData.hh.

107{
108 return (mat == materials[0] || mat == materials[1] || mat == materials[2]);
109}

References materials.

◆ operator=()

G4ICRU90StoppingData & G4ICRU90StoppingData::operator= ( const G4ICRU90StoppingData right)
delete

Field Documentation

◆ isInitialized

G4bool G4ICRU90StoppingData::isInitialized
private

Definition at line 101 of file G4ICRU90StoppingData.hh.

Referenced by Initialise().

◆ materials

const G4Material* G4ICRU90StoppingData::materials[nvectors]
private

Definition at line 98 of file G4ICRU90StoppingData.hh.

Referenced by G4ICRU90StoppingData(), GetIndex(), Initialise(), and IsApplicable().

◆ nvectors

constexpr G4int G4ICRU90StoppingData::nvectors = 3
staticconstexprprivate

◆ sdata_alpha

G4PhysicsFreeVector* G4ICRU90StoppingData::sdata_alpha[nvectors]
private

◆ sdata_proton

G4PhysicsFreeVector* G4ICRU90StoppingData::sdata_proton[nvectors]
private

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