Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4WentzelVIRelModel.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//
26//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class file
30//
31//
32// File name: G4WentzelVIRelModel
33//
34// Author: V.Ivanchenko
35//
36// Creation date: 08.06.2012 from G4WentzelVIRelModel
37//
38// Modifications:
39//
40// Class Description:
41//
42// Implementation of the model of multiple scattering based on
43// G.Wentzel, Z. Phys. 40 (1927) 590.
44// H.W.Lewis, Phys Rev 78 (1950) 526.
45// J.M. Fernandez-Varea et al., NIM B73 (1993) 447.
46// L.Urban, CERN-OPEN-2006-077.
47
48// -------------------------------------------------------------------
49//
50
51//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
52//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
53
58#include "G4SystemOfUnits.hh"
59#include "G4Material.hh"
60#include "G4ElementVector.hh"
62#include "G4NistManager.hh"
63#include "G4EmParameters.hh"
64#include "G4AutoLock.hh"
65
66//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
67
68std::vector<G4double> G4WentzelVIRelModel::effMass;
69
70namespace
71{
72 G4Mutex theWVIRelMutex = G4MUTEX_INITIALIZER;
73}
74
76 G4WentzelVIModel(true, "WentzelVIRel")
77{
78 fNistManager = G4NistManager::Instance();
79 auto ptr = new G4WentzelVIRelXSection();
81}
82
83//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
84
86
87//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
88
90 const G4DataVector& cuts)
91{
92 // Access to materials
93 const G4ProductionCutsTable* theCoupleTable =
95 if(theCoupleTable->GetTableSize() != effMass.size()) {
96 ComputeEffectiveMass();
97 }
98
100}
101
102//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
103
105{
106 if(cup != currentCouple) {
107 currentCouple = cup;
108 SetCurrentCouple(cup);
112 }
113}
114
115//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
116
118 const G4ParticleDefinition* p,
119 G4double kinEnergy,
121 G4double cutEnergy, G4double)
122{
123 G4double cross = 0.0;
124 if(p != particle) { SetupParticle(p); }
125 if(kinEnergy < lowEnergyLimit) { return cross; }
126 if(nullptr == CurrentCouple()) {
127 G4Exception("G4WentzelVIRelModel::ComputeCrossSectionPerAtom", "em0011",
128 FatalException, " G4MaterialCutsCouple is not defined");
129 return cross;
130 }
132 G4int iz = G4lrint(Z);
133 G4double tmass = (1 == iz) ? CLHEP::proton_mass_c2
134 : fNistManager->GetAtomicMassAmu(iz)*amu_c2;
135 wokvi->SetTargetMass(tmass);
137 if(cosTetMaxNuc < 1.0) {
138 G4double cost = wokvi->SetupTarget(iz, cutEnergy);
140 /*
141 //if(p->GetParticleName() == "e-")
142 G4cout << "G4WentzelVIRelModel::CS: Z= " << G4int(Z)
143 << " e(MeV)= " << kinEnergy
144 << " 1-cosN= " << 1 - cosTetMaxNuc << " cross(bn)= " << cross/barn
145 << " " << particle->GetParticleName() << G4endl;
146 */
147 }
148 return cross;
149}
150
151//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
152
153void G4WentzelVIRelModel::ComputeEffectiveMass()
154{
155 const G4ProductionCutsTable* theCoupleTable =
157 std::size_t ncouples = (G4int)theCoupleTable->GetTableSize();
158 if(ncouples == effMass.size()) { return; }
159
160 G4AutoLock l(&theWVIRelMutex);
161 if(ncouples != effMass.size()) {
162 effMass.resize(ncouples, 0.0);
163 for(std::size_t i=0; i<ncouples; ++i) {
164 const G4Material* mat =
165 theCoupleTable->GetMaterialCutsCouple((G4int)i)->GetMaterial();
166 const G4ElementVector* elmVector = mat->GetElementVector();
167 std::size_t nelm = mat->GetNumberOfElements();
168 G4double sum = 0.0;
169 G4double norm= 0.0;
170 for(std::size_t j=0; j<nelm; ++j) {
171 G4int Z = (*elmVector)[j]->GetZasInt();
172 G4double mass = fNistManager->GetAtomicMassAmu(Z)*CLHEP::amu_c2;
173 G4int Z2 = Z*Z;
174 sum += mass*Z2;
175 norm += Z2;
176 }
177 effMass[i] = sum/norm;
178 }
179 }
180 l.unlock();
181}
182
183//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
std::vector< const G4Element * > G4ElementVector
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:85
std::mutex G4Mutex
Definition: G4Threading.hh:81
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4Material * GetMaterial() const
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:185
size_t GetNumberOfElements() const
Definition: G4Material.hh:181
static G4NistManager * Instance()
G4double GetAtomicMassAmu(const G4String &symb) const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:468
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:486
G4double SetupTarget(G4int Z, G4double cut)
G4double ComputeTransportCrossSectionPerAtom(G4double CosThetaMax)
virtual G4double SetupKinematic(G4double kinEnergy, const G4Material *mat)
const G4MaterialCutsCouple * currentCouple
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
void SetWVICrossSection(G4WentzelOKandVIxSection *)
G4WentzelOKandVIxSection * GetWVICrossSection()
const G4ParticleDefinition * particle
const G4Material * currentMaterial
void SetupParticle(const G4ParticleDefinition *)
G4WentzelOKandVIxSection * wokvi
~G4WentzelVIRelModel() override
void DefineMaterial(const G4MaterialCutsCouple *cup)
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double KineticEnergy, G4double AtomicNumber, G4double AtomicWeight=0., G4double cut=DBL_MAX, G4double emax=DBL_MAX) override
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
int G4lrint(double ad)
Definition: templates.hh:134