Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4HadronElasticProcess.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// Geant4 Hadron Elastic Scattering Process
27//
28// Created 26 July 2012 V.Ivanchenko from G4WHadronElasticProcess
29//
30// Modified:
31// 14-Sep-12 M.Kelsey -- Pass subType code to base ctor
32
33#include <iostream>
34
36#include "G4SystemOfUnits.hh"
37#include "G4Nucleus.hh"
38#include "G4ProcessManager.hh"
45
48 fDiffraction(nullptr), fDiffractionRatio(nullptr)
49{}
50
52{}
53
54void G4HadronElasticProcess::ProcessDescription(std::ostream& outFile) const
55{
56 outFile << "G4HadronElasticProcess handles the elastic scattering of \n"
57 << "hadrons by invoking the following hadronic model(s) and \n"
58 << "hadronic cross section(s).\n";
59}
60
63 const G4Step&)
64{
67 G4double weight = track.GetWeight();
69
70 // For elastic scattering, _any_ result is considered an interaction
72
73 G4double kineticEnergy = track.GetKineticEnergy();
74 G4TrackStatus status = track.GetTrackStatus();
75 if(kineticEnergy == 0.0 || track.GetTrackStatus() != fAlive) {
76 return theTotalResult;
77 }
78
79 const G4DynamicParticle* dynParticle = track.GetDynamicParticle();
80 const G4ParticleDefinition* part = dynParticle->GetDefinition();
81 const G4Material* material = track.GetMaterial();
82 G4Nucleus* targNucleus = GetTargetNucleusPointer();
83
84 // Select element
85 const G4Element* elm =
86 GetCrossSectionDataStore()->SampleZandA(dynParticle, material, *targNucleus);
87
88 // Initialize the hadronic projectile from the track
89 G4HadProjectile theProj(track);
90 G4HadronicInteraction* hadi = nullptr;
91 G4HadFinalState* result = nullptr;
92
93 if(fDiffraction)
94 {
95 G4double ratio =
96 fDiffractionRatio->ComputeRatio(part, kineticEnergy,
97 targNucleus->GetZ_asInt(),
98 targNucleus->GetA_asInt());
99 // diffraction is chosen
100 if(ratio > 0.0 && G4UniformRand() < ratio)
101 {
102 try
103 {
104 result = fDiffraction->ApplyYourself(theProj, *targNucleus);
105 }
106 catch(G4HadronicException & aR)
107 {
109 aR.Report(ed);
110 ed << "Call for " << fDiffraction->GetModelName() << G4endl;
111 ed << "Target element "<< elm->GetName()<<" Z= "
112 << targNucleus->GetZ_asInt()
113 << " A= " << targNucleus->GetA_asInt() << G4endl;
114 DumpState(track,"ApplyYourself",ed);
115 ed << " ApplyYourself failed" << G4endl;
116 G4Exception("G4HadronElasticProcess::PostStepDoIt", "had006",
117 FatalException, ed);
118 }
119 // Check the result for catastrophic energy non-conservation
120 result = CheckResult(theProj, *targNucleus, result);
121
122 result->SetTrafoToLab(theProj.GetTrafoToLab());
124
125 FillResult(result, track);
126
127 if (epReportLevel != 0) {
128 CheckEnergyMomentumConservation(track, *targNucleus);
129 }
130 return theTotalResult;
131 }
132 }
133
134 // ordinary elastic scattering
135 try
136 {
137 hadi = ChooseHadronicInteraction( theProj, *targNucleus, material, elm );
138 }
139 catch(G4HadronicException & aE)
140 {
142 aE.Report(ed);
143 ed << "Target element "<< elm->GetName()<<" Z= "
144 << targNucleus->GetZ_asInt() << " A= "
145 << targNucleus->GetA_asInt() << G4endl;
146 DumpState(track,"ChooseHadronicInteraction",ed);
147 ed << " No HadronicInteraction found out" << G4endl;
148 G4Exception("G4HadronElasticProcess::PostStepDoIt", "had005",
149 FatalException, ed);
150 }
151
152 size_t idx = track.GetMaterialCutsCouple()->GetIndex();
154 ->GetEnergyCutsVector(3)))[idx];
155 hadi->SetRecoilEnergyThreshold(tcut);
156
157 if(verboseLevel>1) {
158 G4cout << "G4HadronElasticProcess::PostStepDoIt for "
159 << part->GetParticleName()
160 << " in " << material->GetName()
161 << " Target Z= " << targNucleus->GetZ_asInt()
162 << " A= " << targNucleus->GetA_asInt()
163 << " Tcut(MeV)= " << tcut << G4endl;
164 }
165
166 try
167 {
168 result = hadi->ApplyYourself( theProj, *targNucleus);
169 }
170 catch(G4HadronicException & aR)
171 {
173 aR.Report(ed);
174 ed << "Call for " << hadi->GetModelName() << G4endl;
175 ed << "Target element "<< elm->GetName()<<" Z= "
176 << targNucleus->GetZ_asInt()
177 << " A= " << targNucleus->GetA_asInt() << G4endl;
178 DumpState(track,"ApplyYourself",ed);
179 ed << " ApplyYourself failed" << G4endl;
180 G4Exception("G4HadronElasticProcess::PostStepDoIt", "had006",
181 FatalException, ed);
182 }
183
184 // Check the result for catastrophic energy non-conservation
185 // cannot be applied because is not guranteed that recoil
186 // nucleus is created
187 // result = CheckResult(theProj, targNucleus, result);
188
189 // directions
190 G4ThreeVector indir = track.GetMomentumDirection();
191 G4ThreeVector outdir = result->GetMomentumChange();
192
193 if(verboseLevel>1) {
194 G4cout << "Efin= " << result->GetEnergyChange()
195 << " de= " << result->GetLocalEnergyDeposit()
196 << " nsec= " << result->GetNumberOfSecondaries()
197 << " dir= " << outdir
198 << G4endl;
199 }
200
201 // energies
202 G4double edep = std::max(result->GetLocalEnergyDeposit(), 0.0);
203 G4double efinal = std::max(result->GetEnergyChange(), 0.0);
204
205 // primary change
207
208 if(efinal > 0.0) {
209 outdir.rotateUz(indir);
211 } else {
212 if(part->GetProcessManager()->GetAtRestProcessVector()->size() > 0)
213 { status = fStopButAlive; }
214 else { status = fStopAndKill; }
216 }
217
218 //G4cout << "Efinal= " << efinal << " TrackStatus= " << status << G4endl;
220
221 // recoil
222 if(result->GetNumberOfSecondaries() > 0) {
223 G4DynamicParticle* p = result->GetSecondary(0)->GetParticle();
224
225 if(p->GetKineticEnergy() > tcut) {
228 // G4cout << "recoil " << pdir << G4endl;
229 pdir.rotateUz(indir);
230 // G4cout << "recoil rotated " << pdir << G4endl;
231 p->SetMomentumDirection(pdir);
232
233 // in elastic scattering time and weight are not changed
234 G4Track* t = new G4Track(p, track.GetGlobalTime(),
235 track.GetPosition());
236 t->SetWeight(weight);
239
240 } else {
241 edep += p->GetKineticEnergy();
242 delete p;
243 }
244 }
247 result->Clear();
248
249 return theTotalResult;
250}
251
253{
254 PrintWarning("G4HadronElasticProcess::SetLowestEnergy(..) ");
255}
256
257void
259{
260 PrintWarning("G4HadronElasticProcess::SetLowestEnergyNeutron(..) ");
261}
262
265{
266 if(hi && xsr) {
267 fDiffraction = hi;
268 fDiffractionRatio = xsr;
269 }
270}
271
272void G4HadronElasticProcess::PrintWarning(const G4String& tit) const
273{
274 G4Exception(tit, "had003", JustWarning,
275 " method is obsolete and will be removed in the next release");
276}
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
@ fHadronElastic
G4TrackStatus
@ fAlive
@ fStopAndKill
@ fStopButAlive
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
const G4Element * SampleZandA(const G4DynamicParticle *, const G4Material *, G4Nucleus &target)
void SetMomentumDirection(const G4ThreeVector &aDirection)
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4String & GetName() const
Definition: G4Element.hh:126
G4double GetEnergyChange() const
void SetTrafoToLab(const G4LorentzRotation &aT)
G4double GetLocalEnergyDeposit() const
const G4ThreeVector & GetMomentumChange() const
std::size_t GetNumberOfSecondaries() const
G4HadSecondary * GetSecondary(size_t i)
G4LorentzRotation & GetTrafoToLab()
G4DynamicParticle * GetParticle()
void SetDiffraction(G4HadronicInteraction *, G4VCrossSectionRatio *)
G4HadronElasticProcess(const G4String &procName="hadElastic")
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
void ProcessDescription(std::ostream &outFile) const override
virtual void SetLowestEnergy(G4double)
virtual void SetLowestEnergyNeutron(G4double)
void Report(std::ostream &aS) const
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
const G4String & GetModelName() const
void SetRecoilEnergyThreshold(G4double val)
void FillResult(G4HadFinalState *aR, const G4Track &aT)
G4Nucleus * GetTargetNucleusPointer()
G4HadFinalState * CheckResult(const G4HadProjectile &thePro, const G4Nucleus &targetNucleus, G4HadFinalState *result)
G4ParticleChange * theTotalResult
G4HadronicInteraction * ChooseHadronicInteraction(const G4HadProjectile &aHadProjectile, G4Nucleus &aTargetNucleus, const G4Material *aMaterial, const G4Element *anElement)
void CheckEnergyMomentumConservation(const G4Track &, const G4Nucleus &)
G4CrossSectionDataStore * GetCrossSectionDataStore()
void DumpState(const G4Track &, const G4String &, G4ExceptionDescription &)
const G4String & GetName() const
Definition: G4Material.hh:175
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
void AddSecondary(G4Track *aSecondary)
void ProposeEnergy(G4double finalEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
virtual void Initialize(const G4Track &)
G4ProcessManager * GetProcessManager() const
const G4String & GetParticleName() const
G4ProcessVector * GetAtRestProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
std::size_t size() const
const std::vector< G4double > * GetEnergyCutsVector(std::size_t pcIdx) const
static G4ProductionCutsTable * GetProductionCutsTable()
Definition: G4Step.hh:62
G4TrackStatus GetTrackStatus() const
G4double GetWeight() const
void SetWeight(G4double aValue)
const G4ThreeVector & GetPosition() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4double GetGlobalTime() const
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
const G4TouchableHandle & GetTouchableHandle() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
virtual G4double ComputeRatio(const G4ParticleDefinition *, G4double kinEnergy, G4int, G4int)=0
void ProposeTrackStatus(G4TrackStatus status)
void ProposeNonIonizingEnergyDeposit(G4double anEnergyPart)
void ProposeWeight(G4double finalWeight)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetNumberOfSecondaries(G4int totSecondaries)
void ClearNumberOfInteractionLengthLeft()
Definition: G4VProcess.hh:424
G4int verboseLevel
Definition: G4VProcess.hh:356