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