Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4HadronElasticProcess Class Reference

#include <G4HadronElasticProcess.hh>

+ Inheritance diagram for G4HadronElasticProcess:

Public Member Functions

 G4HadronElasticProcess (const G4String &procName="hadElastic")
 
 ~G4HadronElasticProcess () override
 
G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep) override
 
virtual void SetLowestEnergy (G4double)
 
virtual void SetLowestEnergyNeutron (G4double)
 
void ProcessDescription (std::ostream &outFile) const override
 
void SetDiffraction (G4HadronicInteraction *, G4VCrossSectionRatio *)
 
- Public Member Functions inherited from G4HadronicProcess
 G4HadronicProcess (const G4String &processName="Hadronic", G4ProcessType procType=fHadronic)
 
 G4HadronicProcess (const G4String &processName, G4HadronicProcessType subType)
 
 ~G4HadronicProcess () override
 
void RegisterMe (G4HadronicInteraction *a)
 
G4double GetElementCrossSection (const G4DynamicParticle *part, const G4Element *elm, const G4Material *mat=nullptr)
 
G4double GetMicroscopicCrossSection (const G4DynamicParticle *part, const G4Element *elm, const G4Material *mat=nullptr)
 
void StartTracking (G4Track *track) override
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double, G4ForceCondition *) override
 
G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep) override
 
void PreparePhysicsTable (const G4ParticleDefinition &) override
 
void BuildPhysicsTable (const G4ParticleDefinition &) override
 
void DumpPhysicsTable (const G4ParticleDefinition &p)
 
void AddDataSet (G4VCrossSectionDataSet *aDataSet)
 
std::vector< G4HadronicInteraction * > & GetHadronicInteractionList ()
 
G4HadronicInteractionGetHadronicModel (const G4String &)
 
G4HadronicInteractionGetHadronicInteraction () const
 
G4double GetMeanFreePath (const G4Track &aTrack, G4double, G4ForceCondition *) override
 
const G4NucleusGetTargetNucleus () const
 
G4NucleusGetTargetNucleusPointer ()
 
const G4IsotopeGetTargetIsotope ()
 
G4double ComputeCrossSection (const G4ParticleDefinition *, const G4Material *, const G4double kinEnergy)
 
G4HadXSType CrossSectionType () const
 
void SetCrossSectionType (G4HadXSType val)
 
void ProcessDescription (std::ostream &outFile) const override
 
void BiasCrossSectionByFactor (G4double aScale)
 
void MultiplyCrossSectionBy (G4double factor)
 
G4double CrossSectionFactor () const
 
void SetIntegral (G4bool val)
 
void SetEpReportLevel (G4int level)
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
G4CrossSectionDataStoreGetCrossSectionDataStore ()
 
std::vector< G4TwoPeaksHadXS * > * TwoPeaksXS () const
 
std::vector< G4double > * EnergyOfCrossSectionMax () const
 
G4HadronicProcessoperator= (const G4HadronicProcess &right)=delete
 
 G4HadronicProcess (const G4HadronicProcess &)=delete
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &aName, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
G4VDiscreteProcessoperator= (const G4VDiscreteProcess &)=delete
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double GetCrossSection (const G4double, const G4MaterialCutsCouple *)
 
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4VProcessoperator= (const G4VProcess &)=delete
 
G4bool operator== (const G4VProcess &right) const
 
G4bool operator!= (const G4VProcess &right) const
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual const G4VProcessGetCreatorProcess () const
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Member Functions inherited from G4HadronicProcess
G4HadronicInteractionChooseHadronicInteraction (const G4HadProjectile &aHadProjectile, G4Nucleus &aTargetNucleus, const G4Material *aMaterial, const G4Element *anElement)
 
G4double GetLastCrossSection ()
 
void FillResult (G4HadFinalState *aR, const G4Track &aT)
 
void DumpState (const G4Track &, const G4String &, G4ExceptionDescription &)
 
G4HadFinalStateCheckResult (const G4HadProjectile &thePro, const G4Nucleus &targetNucleus, G4HadFinalState *result)
 
void CheckEnergyMomentumConservation (const G4Track &, const G4Nucleus &)
 
- Protected Member Functions inherited from G4VDiscreteProcess
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- Protected Attributes inherited from G4HadronicProcess
G4HadProjectile thePro
 
G4ParticleChangetheTotalResult
 
G4CrossSectionDataStoretheCrossSectionDataStore
 
G4double fWeight = 1.0
 
G4double aScaleFactor = 1.0
 
G4double theLastCrossSection = 0.0
 
G4double mfpKinEnergy = DBL_MAX
 
G4int epReportLevel = 0
 
G4HadXSType fXSType = fHadNoIntegral
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft = -1.0
 
G4double currentInteractionLength = -1.0
 
G4double theInitialNumberOfInteractionLength = -1.0
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType = fNotDefined
 
G4int theProcessSubType = -1
 
G4double thePILfactor = 1.0
 
G4int verboseLevel = 0
 
G4bool enableAtRestDoIt = true
 
G4bool enableAlongStepDoIt = true
 
G4bool enablePostStepDoIt = true
 

Detailed Description

Definition at line 49 of file G4HadronElasticProcess.hh.

Constructor & Destructor Documentation

◆ G4HadronElasticProcess()

G4HadronElasticProcess::G4HadronElasticProcess ( const G4String & procName = "hadElastic")
explicit

Definition at line 46 of file G4HadronElasticProcess.cc.

48 fDiffraction(nullptr), fDiffractionRatio(nullptr)
49{}
G4HadronicProcess(const G4String &processName="Hadronic", G4ProcessType procType=fHadronic)

◆ ~G4HadronElasticProcess()

G4HadronElasticProcess::~G4HadronElasticProcess ( )
override

Definition at line 51 of file G4HadronElasticProcess.cc.

52{}

Member Function Documentation

◆ PostStepDoIt()

G4VParticleChange * G4HadronElasticProcess::PostStepDoIt ( const G4Track & aTrack,
const G4Step & aStep )
overridevirtual

Reimplemented from G4VDiscreteProcess.

Definition at line 62 of file G4HadronElasticProcess.cc.

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);
233 t->SetTouchableHandle(track.GetTouchableHandle());
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}
@ 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)
G4DynamicParticle * GetParticle()
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 &)
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
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()
void SetWeight(G4double aValue)
void SetTouchableHandle(const G4TouchableHandle &apValue)
void SetCreatorModelID(const G4int id)
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

◆ ProcessDescription()

void G4HadronElasticProcess::ProcessDescription ( std::ostream & outFile) const
overridevirtual

Reimplemented from G4VProcess.

Definition at line 54 of file G4HadronElasticProcess.cc.

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}

◆ SetDiffraction()

void G4HadronElasticProcess::SetDiffraction ( G4HadronicInteraction * hi,
G4VCrossSectionRatio * xsr )

Definition at line 261 of file G4HadronElasticProcess.cc.

263{
264 if(hi && xsr) {
265 fDiffraction = hi;
266 fDiffractionRatio = xsr;
267 }
268}

Referenced by G4HadronHElasticPhysics::ConstructProcess().

◆ SetLowestEnergy()

void G4HadronElasticProcess::SetLowestEnergy ( G4double )
virtual

Definition at line 250 of file G4HadronElasticProcess.cc.

251{
252 PrintWarning("G4HadronElasticProcess::SetLowestEnergy(..) ");
253}

◆ SetLowestEnergyNeutron()

void G4HadronElasticProcess::SetLowestEnergyNeutron ( G4double )
virtual

Definition at line 256 of file G4HadronElasticProcess.cc.

257{
258 PrintWarning("G4HadronElasticProcess::SetLowestEnergyNeutron(..) ");
259}

The documentation for this class was generated from the following files: