Geant4 11.2.2
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
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 // ClearNumberOfInteractionLengthLeft();
73
74 const G4DynamicParticle* dynParticle = track.GetDynamicParticle();
75 G4double kineticEnergy = dynParticle->GetKineticEnergy();
76 G4TrackStatus status = track.GetTrackStatus();
77 if(kineticEnergy == 0.0 || track.GetTrackStatus() != fAlive) {
78 return theTotalResult;
79 }
80
81 const G4Material* material = track.GetMaterial();
82
83 // check only for charged particles
84 if(fXSType != fHadNoIntegral) {
87 theCrossSectionDataStore->ComputeCrossSection(dynParticle, material);
89 // No interaction
90 return theTotalResult;
91 }
92 }
93
94 const G4ParticleDefinition* part = dynParticle->GetDefinition();
95 G4Nucleus* targNucleus = GetTargetNucleusPointer();
96
97 // Select element
98 const G4Element* elm =
99 theCrossSectionDataStore->SampleZandA(dynParticle, material, *targNucleus);
100
101 // Initialize the hadronic projectile from the track
102 G4HadProjectile theProj(track);
103 G4HadronicInteraction* hadi = nullptr;
104 G4HadFinalState* result = nullptr;
105
106 if(nullptr != fDiffraction) {
107 G4double ratio =
108 fDiffractionRatio->ComputeRatio(part, kineticEnergy,
109 targNucleus->GetZ_asInt(),
110 targNucleus->GetA_asInt());
111 // diffraction is chosen
112 if(ratio > 0.0 && G4UniformRand() < ratio)
113 {
114 try
115 {
116 result = fDiffraction->ApplyYourself(theProj, *targNucleus);
117 }
118 catch(G4HadronicException & aR)
119 {
121 aR.Report(ed);
122 ed << "Call for " << fDiffraction->GetModelName() << G4endl;
123 ed << part->GetParticleName()
124 << " off target element " << elm->GetName() << " Z= "
125 << targNucleus->GetZ_asInt()
126 << " A= " << targNucleus->GetA_asInt() << G4endl;
127 DumpState(track,"ApplyYourself",ed);
128 ed << " ApplyYourself failed" << G4endl;
129 G4Exception("G4HadronElasticProcess::PostStepDoIt", "had006",
130 FatalException, ed);
131 }
132 // Check the result for catastrophic energy non-conservation
133 result = CheckResult(theProj, *targNucleus, result);
134 result->SetTrafoToLab(theProj.GetTrafoToLab());
135
136 // The following method of the base class takes care also of setting
137 // the creator model ID for the secondaries that are created
138 FillResult(result, track);
139
140 if (epReportLevel != 0) {
141 CheckEnergyMomentumConservation(track, *targNucleus);
142 }
143 return theTotalResult;
144 }
145 }
146
147 // ordinary elastic scattering
148 hadi = ChooseHadronicInteraction( theProj, *targNucleus, material, elm );
149 if(nullptr == hadi) {
151 ed << part->GetParticleName()
152 << " off target element " << elm->GetName() << " Z= "
153 << targNucleus->GetZ_asInt() << " A= "
154 << targNucleus->GetA_asInt() << G4endl;
155 DumpState(track,"ChooseHadronicInteraction",ed);
156 ed << " No HadronicInteraction found out" << G4endl;
157 G4Exception("G4HadronElasticProcess::PostStepDoIt", "had005",
158 FatalException, ed);
159 return theTotalResult;
160 }
161
162 size_t idx = track.GetMaterialCutsCouple()->GetIndex();
164 ->GetEnergyCutsVector(3)))[idx];
165 hadi->SetRecoilEnergyThreshold(tcut);
166 /*
167 if(verboseLevel>1) {
168 G4cout << "G4HadronElasticProcess::PostStepDoIt for "
169 << part->GetParticleName()
170 << " in " << material->GetName()
171 << " Target Z= " << targNucleus->GetZ_asInt()
172 << " A= " << targNucleus->GetA_asInt()
173 << " Tcut(MeV)= " << tcut << G4endl;
174 }
175 */
176 result = hadi->ApplyYourself( theProj, *targNucleus);
177
178 // Check the result for catastrophic energy non-conservation
179 // cannot be applied because is not guranteed that recoil
180 // nucleus is created
181 // result = CheckResult(theProj, targNucleus, result);
182
183 // directions
184 G4ThreeVector indir = track.GetMomentumDirection();
185 G4ThreeVector outdir = result->GetMomentumChange();
186 /*
187 if(verboseLevel>1) {
188 G4cout << "Efin= " << result->GetEnergyChange()
189 << " de= " << result->GetLocalEnergyDeposit()
190 << " nsec= " << result->GetNumberOfSecondaries()
191 << " dir= " << outdir
192 << G4endl;
193 }
194 */
195 // energies
196 G4double edep = std::max(result->GetLocalEnergyDeposit(), 0.0);
197 G4double efinal = std::max(result->GetEnergyChange(), 0.0);
198
199 // primary change
201
202 if(efinal > 0.0) {
203 outdir.rotateUz(indir);
205 } else {
206 if(part->GetProcessManager()->GetAtRestProcessVector()->size() > 0)
207 { status = fStopButAlive; }
208 else { status = fStopAndKill; }
210 }
211 /*
212 G4cout << "Efinal= " << efinal << " TrackStatus= " << status
213 << " time(ns)=" << track.GetGlobalTime()/ns << G4endl;
214 */
216
217 // recoil
218 if(result->GetNumberOfSecondaries() > 0) {
219 G4DynamicParticle* p = result->GetSecondary(0)->GetParticle();
220
221 if(p->GetKineticEnergy() > tcut) {
224 // G4cout << "recoil " << pdir << G4endl;
225 pdir.rotateUz(indir);
226 // G4cout << "recoil rotated " << pdir << G4endl;
227 p->SetMomentumDirection(pdir);
228
229 // in elastic scattering time and weight are not changed
230 G4Track* t = new G4Track(p, track.GetGlobalTime(),
231 track.GetPosition());
232 t->SetWeight(weight);
234 G4int secID = G4PhysicsModelCatalog::GetModelID( "model_" + hadi->GetModelName() );
235 if ( secID > 0 ) t->SetCreatorModelID(secID);
237
238 } else {
239 edep += p->GetKineticEnergy();
240 delete p;
241 }
242 }
245 result->Clear();
246
247 return theTotalResult;
248}
249
251{
252 PrintWarning("G4HadronElasticProcess::SetLowestEnergy(..) ");
253}
254
255void
257{
258 PrintWarning("G4HadronElasticProcess::SetLowestEnergyNeutron(..) ");
259}
260
263{
264 if(hi && xsr) {
265 fDiffraction = hi;
266 fDiffractionRatio = xsr;
267 }
268}
269
270void G4HadronElasticProcess::PrintWarning(const G4String& tit) const
271{
272 G4Exception(tit, "had003", JustWarning,
273 " method is obsolete and will be removed in the next release");
274}
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
@ fHadNoIntegral
G4TrackStatus
@ fAlive
@ fStopAndKill
@ fStopButAlive
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
#define G4UniformRand()
Definition Randomize.hh:52
Hep3Vector & rotateUz(const Hep3Vector &)
G4double ComputeCrossSection(const G4DynamicParticle *, const G4Material *)
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:115
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)
G4CrossSectionDataStore * theCrossSectionDataStore
void CheckEnergyMomentumConservation(const G4Track &, const G4Nucleus &)
void DumpState(const G4Track &, const G4String &, G4ExceptionDescription &)
G4int GetA_asInt() const
Definition G4Nucleus.hh:99
G4int GetZ_asInt() const
Definition G4Nucleus.hh:105
void AddSecondary(G4Track *aSecondary)
void Initialize(const G4Track &) override
void ProposeEnergy(G4double finalEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4ProcessManager * GetProcessManager() const
const G4String & GetParticleName() const
static G4int GetModelID(const G4int modelIndex)
G4ProcessVector * GetAtRestProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
std::size_t size() const
const std::vector< G4double > * GetEnergyCutsVector(std::size_t pcIdx) const
static G4ProductionCutsTable * GetProductionCutsTable()
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
void SetCreatorModelID(const G4int id)
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)
G4double theNumberOfInteractionLengthLeft
#define DBL_MAX
Definition templates.hh:62