Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4hCoulombScatteringModel.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: G4hCoulombScatteringModel.cc 104802 2017-06-19 07:11:40Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4hCoulombScatteringModel
34 //
35 // Author: Vladimir Ivanchenko
36 //
37 // Creation date: 08.06.2012 from G4eCoulombScatteringModel
38 //
39 // Modifications:
40 //
41 //
42 // Class Description:
43 //
44 // -------------------------------------------------------------------
45 //
46 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
47 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
48 
50 #include "G4PhysicalConstants.hh"
51 #include "G4SystemOfUnits.hh"
52 #include "Randomize.hh"
53 #include "G4DataVector.hh"
54 #include "G4ElementTable.hh"
56 #include "G4Proton.hh"
57 #include "G4ParticleTable.hh"
58 #include "G4IonTable.hh"
59 #include "G4ProductionCutsTable.hh"
60 #include "G4NucleiProperties.hh"
61 #include "G4Pow.hh"
62 #include "G4LossTableManager.hh"
63 #include "G4LossTableBuilder.hh"
64 #include "G4NistManager.hh"
65 
66 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
67 
68 using namespace std;
69 
71  : G4VEmModel("hCoulombScattering"),
72  cosThetaMin(1.0),
73  cosThetaMax(-1.0),
74  isCombined(combined)
75 {
76  fParticleChange = nullptr;
80  currentMaterial = nullptr;
81  fixedCut = -1.0;
82 
83  pCuts = nullptr;
84 
85  recoilThreshold = 0.0; // by default does not work
86 
87  particle = nullptr;
88  currentCouple = nullptr;
90 
93  elecRatio = 0.0;
94 }
95 
96 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
97 
99 {
100  delete wokvi;
101 }
102 
103 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
104 
106  const G4DataVector& cuts)
107 {
108  SetupParticle(part);
109  currentCouple = nullptr;
110 
111  // defined theta limit between single and multiple scattering
112  isCombined = true;
114 
115  if(tet <= 0.0) {
116  cosThetaMin = 1.0;
117  isCombined = false;
118  } else if(tet >= CLHEP::pi) {
119  cosThetaMin = -1.0;
120  } else {
121  cosThetaMin = cos(tet);
122  }
123 
124  wokvi->Initialise(part, cosThetaMin);
125  /*
126  G4cout << "G4hCoulombScatteringModel: " << particle->GetParticleName()
127  << " 1-cos(ThetaLimit)= " << 1 - cosThetaMin
128  << " cos(thetaMax)= " << cosThetaMax
129  << G4endl;
130  */
131  pCuts = &cuts;
132  //G4ProductionCutsTable::GetProductionCutsTable()->GetEnergyCutsVector(3);
133  /*
134  G4cout << "!!! G4hCoulombScatteringModel::Initialise for "
135  << part->GetParticleName() << " cos(TetMin)= " << cosThetaMin
136  << " cos(TetMax)= " << cosThetaMax <<G4endl;
137  G4cout << "cut= " << (*pCuts)[0] << " cut1= " << (*pCuts)[1] << G4endl;
138  */
139  if(!fParticleChange) {
141  }
142  if(IsMaster() && mass < GeV && part->GetParticleName() != "GenericIon") {
143  InitialiseElementSelectors(part, cuts);
144  }
145 }
146 
147 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
148 
150  G4VEmModel* masterModel)
151 {
153 }
154 
155 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
156 
157 G4double
159  const G4ParticleDefinition* part,
160  G4double)
161 {
162  SetupParticle(part);
163 
164  // define cut using cuts for proton
165  G4double cut =
166  std::max(recoilThreshold, (*pCuts)[CurrentCouple()->GetIndex()]);
167 
168  // find out lightest element
169  const G4ElementVector* theElementVector = material->GetElementVector();
170  G4int nelm = material->GetNumberOfElements();
171 
172  // select lightest element
173  G4int Z = 300;
174  for (G4int j=0; j<nelm; ++j) {
175  Z = std::min(Z,(*theElementVector)[j]->GetZasInt());
176  }
178  G4double targetMass = G4NucleiProperties::GetNuclearMass(A, Z);
179  G4double t = std::max(cut, 0.5*(cut + sqrt(2*cut*targetMass)));
180 
181  return t;
182 }
183 
184 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
185 
187  const G4ParticleDefinition* p,
188  G4double kinEnergy,
190  G4double cutEnergy, G4double)
191 {
192  //G4cout << "### G4hCoulombScatteringModel::ComputeCrossSectionPerAtom for "
193  //<< p->GetParticleName()<<" Z= "<<Z<<" e(MeV)= "<< kinEnergy/MeV << G4endl;
194  G4double cross = 0.0;
195  elecRatio = 0.0;
196  if(p != particle) { SetupParticle(p); }
197 
198  // cross section is set to zero to avoid problems in sample secondary
199  if(kinEnergy <= 0.0) { return cross; }
201 
202  G4int iz = G4lrint(Z);
203  G4double tmass = (1 == iz) ? proton_mass_c2 :
205  wokvi->SetTargetMass(tmass);
206 
207  G4double costmin =
208  wokvi->SetupKinematic(kinEnergy, currentMaterial);
209 
210  if(cosThetaMax < costmin) {
211  G4double cut = (0.0 < fixedCut) ? fixedCut : cutEnergy;
212  costmin = wokvi->SetupTarget(iz, cut);
213  G4double costmax =
214  (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() == "mu+")
222  G4cout << "e(MeV)= " << kinEnergy/MeV << " cross(b)= " << cross/barn
223  << " 1-costmin= " << 1-costmin
224  << " 1-costmax= " << 1-costmax
225  << " 1-cosThetaMax= " << 1-cosThetaMax
226  << G4endl;
227  */
228  }
229  return cross;
230 }
231 
232 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
233 
235  std::vector<G4DynamicParticle*>* fvect,
236  const G4MaterialCutsCouple* couple,
237  const G4DynamicParticle* dp,
238  G4double cutEnergy,
239  G4double)
240 {
241  G4double kinEnergy = dp->GetKineticEnergy();
243  DefineMaterial(couple);
244 
245  // Choose nucleus
246  G4double cut = (0.0 < fixedCut) ? fixedCut : cutEnergy;
247 
248  const G4Element* elm = SelectRandomAtom(couple,particle,
249  kinEnergy,cut,kinEnergy);
250 
251  G4int iz = elm->GetZasInt();
252  G4int ia = SelectIsotopeNumber(elm);
254 
255  wokvi->SetTargetMass(mass2);
256  wokvi->SetupKinematic(kinEnergy, currentMaterial);
257  G4double costmin = wokvi->SetupTarget(iz, cut);
258  G4double costmax = (1 == iz && particle == theProton && cosThetaMax < 0.0)
259  ? 0.0 : cosThetaMax;
260  if(costmin <= costmax) { return; }
261 
262  G4double cross = wokvi->ComputeNuclearCrossSection(costmin, costmax);
263  G4double ecross = wokvi->ComputeElectronCrossSection(costmin, costmax);
264  G4double ratio = ecross/(cross + ecross);
265 
266  G4ThreeVector newDirection =
267  wokvi->SampleSingleScattering(costmin, costmax, ratio);
268 
269  // kinematics in the Lab system
270  G4double ptot = dp->GetTotalMomentum();
271  G4double e1 = dp->GetTotalEnergy();
272 
273  // Lab. system kinematics along projectile direction
274  G4LorentzVector v0 = G4LorentzVector(0, 0, ptot, e1 + mass2);
275  G4double bet = ptot/v0.e();
276  G4double gam = 1.0/sqrt((1.0 - bet)*(1.0 + bet));
277 
278  // CM projectile
279  G4double momCM = gam*(ptot - bet*e1);
280  G4double eCM = gam*(e1 - bet*ptot);
281  // energy & momentum after scattering of incident particle
282  G4double pxCM = momCM*newDirection.x();
283  G4double pyCM = momCM*newDirection.y();
284  G4double pzCM = momCM*newDirection.z();
285 
286  // CM--->Lab
287  G4LorentzVector v1(pxCM , pyCM, gam*(pzCM + bet*eCM), gam*(eCM + bet*pzCM));
288 
290  newDirection = v1.vect().unit();
291  newDirection.rotateUz(dir);
292 
294 
295  // recoil
296  v0 -= v1;
297  G4double trec = v0.e() - mass2;
298  G4double edep = 0.0;
299 
300  G4double tcut = recoilThreshold;
301  if(pCuts) { tcut= std::max(tcut,(*pCuts)[currentMaterialIndex]); }
302 
303  if(trec > tcut) {
305  newDirection = v0.vect().unit();
306  newDirection.rotateUz(dir);
307  G4DynamicParticle* newdp = new G4DynamicParticle(ion, newDirection, trec);
308  fvect->push_back(newdp);
309  } else if(trec > 0.0) {
310  edep = trec;
312  }
313 
314  // finelize primary energy and energy balance
315  G4double finalT = v1.e() - mass;
316  if(finalT < 0.0) {
317  edep += finalT;
318  finalT = 0.0;
319  }
320  edep = std::max(edep, 0.0);
323 }
324 
325 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
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)
static G4ParticleTable * GetParticleTable()
const G4ParticleDefinition * particle
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double) final
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax) override
void SetupParticle(const G4ParticleDefinition *)
G4int SelectIsotopeNumber(const G4Element *)
Definition: G4VEmModel.cc:335
const G4ThreeVector & GetMomentumDirection() const
const char * p
Definition: xmltok.h:285
const std::vector< G4double > * pCuts
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:191
G4IonTable * GetIonTable() const
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
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
static G4Proton * Proton()
Definition: G4Proton.cc:93
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:784
G4double GetAtomicMassAmu(const G4String &symb) const
void ProposeNonIonizingEnergyDeposit(G4double anEnergyPart)
static constexpr double proton_mass_c2
Float_t Z
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
static constexpr double amu_c2
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)
G4double ComputeNuclearCrossSection(G4double CosThetaMin, G4double CosThetaMax)
double A(double temperature)
static G4double GetNuclearMass(const G4double A, const G4double Z)
virtual G4double SetupKinematic(G4double kinEnergy, const G4Material *mat)
G4hCoulombScatteringModel(G4bool combined=true)
Double_t edep
const G4ParticleDefinition * theProton
G4int GetZasInt() const
Definition: G4Element.hh:132
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:462
Hep3Vector unit() const
G4double SetupTarget(G4int Z, G4double cut)
void Initialise(const G4ParticleDefinition *, G4double CosThetaLim)
int G4lrint(double ad)
Definition: templates.hh:151
std::vector< G4Element * > G4ElementVector
int G4int
Definition: G4Types.hh:78
TDirectory * dir
void DefineMaterial(const G4MaterialCutsCouple *)
G4double GetKineticEnergy() const
double x() const
G4bool IsMaster() const
Definition: G4VEmModel.hh:700
const G4MaterialCutsCouple * currentCouple
Hep3Vector vect() const
G4double GetTotalMomentum() const
CLHEP::HepLorentzVector G4LorentzVector
G4double GetTotalEnergy() const
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:504
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:131
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
double y() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4ParticleChangeForGamma * fParticleChange
static constexpr double pi
Definition: SystemOfUnits.h:54
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)
static G4double tet[DIM]
static G4NistManager * Instance()