Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4WentzelVIRelModel.cc
이 파일의 문서화 페이지로 가기
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 // $Id: G4WentzelVIRelModel.cc 104487 2017-06-01 13:41:55Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4WentzelVIRelModel
34 //
35 // Author: V.Ivanchenko
36 //
37 // Creation date: 08.06.2012 from G4WentzelVIRelModel
38 //
39 // Modifications:
40 //
41 // Class Description:
42 //
43 // Implementation of the model of multiple scattering based on
44 // G.Wentzel, Z. Phys. 40 (1927) 590.
45 // H.W.Lewis, Phys Rev 78 (1950) 526.
46 // J.M. Fernandez-Varea et al., NIM B73 (1993) 447.
47 // L.Urban, CERN-OPEN-2006-077.
48 
49 // -------------------------------------------------------------------
50 //
51 
52 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
53 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
54 
55 #include "G4WentzelVIRelModel.hh"
58 #include "G4PhysicalConstants.hh"
59 #include "G4SystemOfUnits.hh"
60 #include "G4Material.hh"
61 #include "G4ElementVector.hh"
62 #include "G4ProductionCutsTable.hh"
63 #include "G4NistManager.hh"
64 #include "G4EmParameters.hh"
65 
66 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
67 
68 std::vector<G4double> G4WentzelVIRelModel::effMass;
69 
70 #ifdef G4MULTITHREADED
71 G4Mutex G4WentzelVIRelModel::WentzelVIRelModelMutex;
72 #endif
73 
75  G4WentzelVIModel(true, "WentzelVIRel")
76 {
79  SetWVICrossSection(static_cast<G4WentzelOKandVIxSection*>(ptr));
80 }
81 
82 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
83 
85 {
86 }
87 
88 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
89 
91  const G4DataVector& cuts)
92 {
93  // Access to materials
94  const G4ProductionCutsTable* theCoupleTable =
96  if(theCoupleTable->GetTableSize() != effMass.size()) {
98  }
99 
101 }
102 
103 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
104 
106 {
107  if(cup != currentCouple) {
108  currentCouple = cup;
109  SetCurrentCouple(cup);
110  currentMaterial = cup->GetMaterial();
113  }
114 }
115 
116 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
117 
119  const G4ParticleDefinition* p,
120  G4double kinEnergy,
122  G4double cutEnergy, G4double)
123 {
124  G4double cross = 0.0;
125  if(p != particle) { SetupParticle(p); }
126  if(kinEnergy < lowEnergyLimit) { return cross; }
127  if(!CurrentCouple()) {
128  G4Exception("G4WentzelVIRelModel::ComputeCrossSectionPerAtom", "em0011",
129  FatalException, " G4MaterialCutsCouple is not defined");
130  return cross;
131  }
133  G4int iz = G4lrint(Z);
134  G4double tmass = (1 == iz) ? CLHEP::proton_mass_c2
136  wokvi->SetTargetMass(tmass);
138  if(cosTetMaxNuc < 1.0) {
139  G4double cost = wokvi->SetupTarget(iz, cutEnergy);
141  /*
142  //if(p->GetParticleName() == "e-")
143  G4cout << "G4WentzelVIRelModel::CS: Z= " << G4int(Z)
144  << " e(MeV)= " << kinEnergy
145  << " 1-cosN= " << 1 - cosTetMaxNuc << " cross(bn)= " << cross/barn
146  << " " << particle->GetParticleName() << G4endl;
147  */
148  }
149  return cross;
150 }
151 
152 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
153 
155 {
156 #ifdef G4MULTITHREADED
157  G4MUTEXLOCK(&G4WentzelVIRelModel::WentzelVIRelModelMutex);
158 #endif
159  const G4ProductionCutsTable* theCoupleTable =
161  size_t ncouples = theCoupleTable->GetTableSize();
162  if(ncouples != effMass.size()) {
163  effMass.resize(ncouples, 0.0);
164  for(size_t i=0; i<ncouples; ++i) {
165  const G4Material* mat =
166  theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
167  const G4ElementVector* elmVector = mat->GetElementVector();
168  G4int nelm = mat->GetNumberOfElements();
169  G4double sum = 0.0;
170  G4double norm= 0.0;
171  for(G4int j=0; j<nelm; ++j) {
172  G4int Z = (*elmVector)[j]->GetZasInt();
174  G4int Z2 = Z*Z;
175  sum += mass*Z2;
176  norm += Z2;
177  }
178  effMass[i] = sum/norm;
179  }
180  }
181 #ifdef G4MULTITHREADED
182  G4MUTEXUNLOCK(&G4WentzelVIRelModel::WentzelVIRelModelMutex);
183 #endif
184 }
185 
186 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4WentzelOKandVIxSection * GetWVICrossSection()
const char * p
Definition: xmltok.h:285
virtual void DefineMaterial(const G4MaterialCutsCouple *cup)
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:191
static std::vector< G4double > effMass
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:455
G4double GetAtomicMassAmu(const G4String &symb) const
G4NistManager * fNistManager
void SetWVICrossSection(G4WentzelOKandVIxSection *)
static constexpr double proton_mass_c2
Float_t Z
double G4double
Definition: G4Types.hh:76
static constexpr double amu_c2
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4WentzelOKandVIxSection * wokvi
Float_t mat
const G4MaterialCutsCouple * currentCouple
#define G4MUTEXUNLOCK(mutex)
Definition: G4Threading.hh:234
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:462
G4double SetupTarget(G4int Z, G4double cut)
Double_t Z2
static G4ProductionCutsTable * GetProductionCutsTable()
const G4ParticleDefinition * particle
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.hh:65
int G4lrint(double ad)
Definition: templates.hh:151
std::vector< G4Element * > G4ElementVector
int G4int
Definition: G4Types.hh:78
Float_t norm
#define G4MUTEXLOCK(mutex)
Definition: G4Threading.hh:233
void SetupParticle(const G4ParticleDefinition *)
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double KineticEnergy, G4double AtomicNumber, G4double AtomicWeight=0., G4double cut=DBL_MAX, G4double emax=DBL_MAX) override
const G4Material * GetMaterial() const
G4double ComputeTransportCrossSectionPerAtom(G4double CosThetaMax)
const G4Material * currentMaterial
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
size_t GetNumberOfElements() const
Definition: G4Material.hh:187
static G4NistManager * Instance()
G4double SetupKinematic(G4double kinEnergy, const G4Material *mat)
std::mutex G4Mutex
Definition: G4Threading.hh:84