Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4eCoulombScatteringModel.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: G4eCoulombScatteringModel.cc 104802 2017-06-19 07:11:40Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4eCoulombScatteringModel
34 //
35 // Author: Vladimir Ivanchenko
36 //
37 // Creation date: 22.08.2005
38 //
39 // Modifications:
40 //
41 // 01.08.06 V.Ivanchenko extend upper limit of table to TeV and review the
42 // logic of building - only elements from G4ElementTable
43 // 08.08.06 V.Ivanchenko build internal table in ekin scale, introduce faclim
44 // 19.08.06 V.Ivanchenko add inline function ScreeningParameter
45 // 09.10.07 V.Ivanchenko reorganized methods, add cut dependence in scattering off e-
46 // 09.06.08 V.Ivanchenko add SelectIsotope and sampling of the recoil ion
47 // 16.06.09 C.Consolandi fixed computation of effective mass
48 // 27.05.10 V.Ivanchenko added G4WentzelOKandVIxSection class to
49 // compute cross sections and sample scattering angle
50 //
51 //
52 // Class Description:
53 //
54 // -------------------------------------------------------------------
55 //
56 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
57 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
58 
60 #include "G4PhysicalConstants.hh"
61 #include "G4SystemOfUnits.hh"
62 #include "Randomize.hh"
63 #include "G4DataVector.hh"
64 #include "G4ElementTable.hh"
66 #include "G4Proton.hh"
67 #include "G4ParticleTable.hh"
68 #include "G4IonTable.hh"
69 #include "G4ProductionCutsTable.hh"
70 #include "G4NucleiProperties.hh"
71 #include "G4Pow.hh"
72 #include "G4LossTableManager.hh"
73 #include "G4LossTableBuilder.hh"
74 #include "G4NistManager.hh"
75 
76 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
77 
78 using namespace std;
79 
81  : G4VEmModel("eCoulombScattering"),
82  cosThetaMin(1.0),
83  cosThetaMax(-1.0),
84  isCombined(combined)
85 {
86  fParticleChange = nullptr;
90  currentMaterial = nullptr;
91  fixedCut = -1.0;
92 
93  pCuts = nullptr;
94 
95  recoilThreshold = 0.0; // by default does not work
96 
97  particle = nullptr;
98  currentCouple = nullptr;
99  wokvi = nullptr;
100 
103  elecRatio = 0.0;
104 }
105 
106 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
107 
109 {
110  delete wokvi;
111 }
112 
113 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
114 
116  const G4DataVector& cuts)
117 {
118  if(!wokvi) { wokvi = new G4WentzelOKandVIxSection(); }
119 
120  SetupParticle(part);
121  currentCouple = nullptr;
122 
123  // defined theta limit between single and multiple scattering
124  if(isCombined) {
125  cosThetaMin = 1.0;
127  if(tet >= pi) { cosThetaMin = -1.0; }
128  else if(tet > 0.0) { cosThetaMin = cos(tet); }
129  }
130 
131  wokvi->Initialise(part, cosThetaMin);
132  /*
133  G4cout << "G4eCoulombScatteringModel: " << particle->GetParticleName()
134  << " 1-cos(ThetaLimit)= " << 1 - cosThetaMin
135  << " cos(thetaMax)= " << cosThetaMax
136  << G4endl;
137  */
138  pCuts = &cuts;
139  //G4ProductionCutsTable::GetProductionCutsTable()->GetEnergyCutsVector(3);
140  /*
141  G4cout << "!!! G4eCoulombScatteringModel::Initialise for "
142  << part->GetParticleName() << " cos(TetMin)= " << cosThetaMin
143  << " cos(TetMax)= " << cosThetaMax <<G4endl;
144  G4cout << "cut= " << (*pCuts)[0] << " cut1= " << (*pCuts)[1] << G4endl;
145  */
146  if(!fParticleChange) {
148  }
149  if(IsMaster() && mass < GeV && part->GetParticleName() != "GenericIon") {
150  InitialiseElementSelectors(part, cuts);
151  }
152 }
153 
154 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
155 
157  G4VEmModel* masterModel)
158 {
160 }
161 
162 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
163 
164 G4double
166  const G4ParticleDefinition* part,
167  G4double)
168 {
169  SetupParticle(part);
170 
171  // define cut using cuts for proton
172  G4double cut =
173  std::max(recoilThreshold, (*pCuts)[CurrentCouple()->GetIndex()]);
174 
175  // find out lightest element
176  const G4ElementVector* theElementVector = material->GetElementVector();
177  G4int nelm = material->GetNumberOfElements();
178 
179  // select lightest element
180  G4int Z = 300;
181  for (G4int j=0; j<nelm; ++j) {
182  Z = std::min(Z,(*theElementVector)[j]->GetZasInt());
183  }
185  G4double targetMass = G4NucleiProperties::GetNuclearMass(A, Z);
186  G4double t = std::max(cut, 0.5*(cut + sqrt(2*cut*targetMass)));
187 
188  return t;
189 }
190 
191 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
192 
194  const G4ParticleDefinition* p,
195  G4double kinEnergy,
197  G4double cutEnergy, G4double)
198 {
199  //G4cout << "### G4eCoulombScatteringModel::ComputeCrossSectionPerAtom for "
200  //<< p->GetParticleName()<<" Z= "<<Z<<" e(MeV)= "<< kinEnergy/MeV << G4endl;
201  G4double cross = 0.0;
202  elecRatio = 0.0;
203  if(p != particle) { SetupParticle(p); }
204 
205  // cross section is set to zero to avoid problems in sample secondary
206  if(kinEnergy <= 0.0) { return cross; }
208  G4double costmin = wokvi->SetupKinematic(kinEnergy, currentMaterial);
209 
210  if(cosThetaMax < costmin) {
211  G4int iz = G4lrint(Z);
212  G4double cut = (0.0 < fixedCut) ? fixedCut : cutEnergy;
213  costmin = wokvi->SetupTarget(iz, cut);
214  G4double costmax = (1 == iz && particle == theProton && cosThetaMax < 0.0)
215  ? 0.0 : cosThetaMax;
216  if(costmin > costmax) {
217  cross = wokvi->ComputeNuclearCrossSection(costmin, costmax)
218  + wokvi->ComputeElectronCrossSection(costmin, costmax);
219  }
220  /*
221  if(p->GetParticleName() == "e-")
222  G4cout << "Z= " << Z << " e(MeV)= " << kinEnergy/MeV
223  << " cross(b)= " << cross/barn << " 1-costmin= " << 1-costmin
224  << " 1-costmax= " << 1-costmax
225  << " 1-cosThetaMax= " << 1-cosThetaMax
226  << " " << currentMaterial->GetName()
227  << G4endl;
228  */
229  }
230  return cross;
231 }
232 
233 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
234 
236  std::vector<G4DynamicParticle*>* fvect,
237  const G4MaterialCutsCouple* couple,
238  const G4DynamicParticle* dp,
239  G4double cutEnergy,
240  G4double)
241 {
242  G4double kinEnergy = dp->GetKineticEnergy();
244  DefineMaterial(couple);
245  /*
246  G4cout << "G4eCoulombScatteringModel::SampleSecondaries e(MeV)= "
247  << kinEnergy << " " << particle->GetParticleName()
248  << " cut= " << cutEnergy<< G4endl;
249  */
250  // Choose nucleus
251  G4double cut = (0.0 < fixedCut) ? fixedCut : cutEnergy;
252 
253  wokvi->SetupKinematic(kinEnergy, currentMaterial);
254 
255  const G4Element* currentElement =
256  SelectRandomAtom(couple,particle,kinEnergy,cut,kinEnergy);
257 
258  G4int iz = currentElement->GetZasInt();
259 
260  G4double costmin = wokvi->SetupTarget(iz, cut);
261  G4double costmax = (1 == iz && particle == theProton && cosThetaMax < 0.0)
262  ? 0.0 : cosThetaMax;
263  if(costmin <= costmax) { return; }
264 
265  G4double cross = wokvi->ComputeNuclearCrossSection(costmin, costmax);
266  G4double ecross = wokvi->ComputeElectronCrossSection(costmin, costmax);
267  G4double ratio = ecross/(cross + ecross);
268 
269  G4int ia = SelectIsotopeNumber(currentElement);
270  G4double targetMass = G4NucleiProperties::GetNuclearMass(ia, iz);
271  wokvi->SetTargetMass(targetMass);
272 
273  G4ThreeVector newDirection =
274  wokvi->SampleSingleScattering(costmin, costmax, ratio);
275  G4double cost = newDirection.z();
276  /*
277  G4cout << "SampleSec: e(MeV)= " << kinEnergy/MeV
278  << " 1-costmin= " << 1-costmin
279  << " 1-costmax= " << 1-costmax
280  << " 1-cost= " << 1-cost
281  << " ratio= " << ratio
282  << G4endl;
283  */
284  G4ThreeVector direction = dp->GetMomentumDirection();
285  newDirection.rotateUz(direction);
286 
288 
289  // recoil sampling assuming a small recoil
290  // and first order correction to primary 4-momentum
291  G4double mom2 = wokvi->GetMomentumSquare();
292  G4double trec = mom2*(1.0 - cost)
293  /(targetMass + (mass + kinEnergy)*(1.0 - cost));
294 
295  // the check likely not needed
296  if(trec > kinEnergy) { trec = kinEnergy; }
297  G4double finalT = kinEnergy - trec;
298  G4double edep = 0.0;
299  /*
300  G4cout<<"G4eCoulombScatteringModel: finalT= "<<finalT<<" Trec= "
301  <<trec << " Z= " << iz << " A= " << ia
302  << " tcut(keV)= " << (*pCuts)[currentMaterialIndex]/keV << G4endl;
303  */
304  G4double tcut = recoilThreshold;
305  if(pCuts) { tcut= std::max(tcut,(*pCuts)[currentMaterialIndex]); }
306 
307  if(trec > tcut) {
309  G4ThreeVector dir = (direction*sqrt(mom2) -
310  newDirection*sqrt(finalT*(2*mass + finalT))).unit();
311  G4DynamicParticle* newdp = new G4DynamicParticle(ion, dir, trec);
312  fvect->push_back(newdp);
313  } else {
314  edep = trec;
316  }
317 
318  // finelize primary energy and energy balance
319  // this threshold may be applied only because for low-enegry
320  // e+e- msc model is applied
321  if(finalT < 0.0) {
322  edep += finalT;
323  finalT = 0.0;
324  }
325  edep = std::max(edep, 0.0);
328 }
329 
330 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:792
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4ThreeVector & SampleSingleScattering(G4double CosThetaMin, G4double CosThetaMax, G4double elecRatio)
const std::vector< G4double > * pCuts
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
static G4ParticleTable * GetParticleTable()
const G4ParticleDefinition * particle
G4int SelectIsotopeNumber(const G4Element *)
Definition: G4VEmModel.cc:335
const G4ThreeVector & GetMomentumDirection() const
const char * p
Definition: xmltok.h:285
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:191
G4IonTable * GetIonTable() const
const G4ParticleDefinition * theProton
void SetProposedKineticEnergy(G4double proposedKinEnergy)
double z() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:38
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:552
static G4Proton * Proton()
Definition: G4Proton.cc:93
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:784
G4double GetAtomicMassAmu(const G4String &symb) const
void SetupParticle(const G4ParticleDefinition *)
void ProposeNonIonizingEnergyDeposit(G4double anEnergyPart)
static constexpr double proton_mass_c2
Float_t Z
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax) override
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
TString part[npart]
Definition: Style.C:32
G4ParticleDefinition * GetDefinition() const
G4double PolarAngleLimit() const
Definition: G4VEmModel.hh:637
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:146
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4WentzelOKandVIxSection * wokvi
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
G4double ComputeNuclearCrossSection(G4double CosThetaMin, G4double CosThetaMax)
double A(double temperature)
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4eCoulombScatteringModel(G4bool combined=true)
Double_t edep
G4int GetZasInt() const
Definition: G4Element.hh:132
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:462
G4double SetupTarget(G4int Z, G4double cut)
void Initialise(const G4ParticleDefinition *, G4double CosThetaLim)
G4ParticleChangeForGamma * fParticleChange
int G4lrint(double ad)
Definition: templates.hh:151
std::vector< G4Element * > G4ElementVector
int G4int
Definition: G4Types.hh:78
void DefineMaterial(const G4MaterialCutsCouple *)
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double) final
TDirectory * dir
G4double GetKineticEnergy() const
G4bool IsMaster() const
Definition: G4VEmModel.hh:700
static constexpr double pi
Definition: G4SIunits.hh:75
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:504
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:131
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
size_t GetNumberOfElements() const
Definition: G4Material.hh:187
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4double ComputeElectronCrossSection(G4double CosThetaMin, G4double CosThetaMax)
const G4MaterialCutsCouple * currentCouple
static G4double tet[DIM]
static G4NistManager * Instance()
G4double SetupKinematic(G4double kinEnergy, const G4Material *mat)