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

#include <G4WHadronElasticProcess.hh>

+ Inheritance diagram for G4WHadronElasticProcess:

Public Member Functions

 G4WHadronElasticProcess (const G4String &procName="hadElastic")
 
virtual ~G4WHadronElasticProcess ()
 
virtual G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep)
 
void SetLowestEnergy (G4double)
 
void SetLowestEnergyNeutron (G4double)
 
virtual void Description () const
 
- Public Member Functions inherited from G4HadronicProcess
 G4HadronicProcess (const G4String &processName="Hadronic", G4ProcessType procType=fHadronic)
 
 G4HadronicProcess (const G4String &processName, G4HadronicProcessType subType)
 
virtual ~G4HadronicProcess ()
 
void RegisterMe (G4HadronicInteraction *a)
 
G4double GetElementCrossSection (const G4DynamicParticle *part, const G4Element *elm, const G4Material *mat=0)
 
G4double GetMicroscopicCrossSection (const G4DynamicParticle *part, const G4Element *elm, const G4Material *mat=0)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
void DumpPhysicsTable (const G4ParticleDefinition &p)
 
void AddDataSet (G4VCrossSectionDataSet *aDataSet)
 
G4EnergyRangeManagerGetManagerPointer ()
 
G4double GetMeanFreePath (const G4Track &aTrack, G4double, G4ForceCondition *)
 
const G4NucleusGetTargetNucleus () const
 
const G4IsotopeGetTargetIsotope ()
 
virtual void ProcessDescription (std::ostream &outFile) const
 
void BiasCrossSectionByFactor (G4double aScale)
 
void SetEpReportLevel (G4int level)
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
G4CrossSectionDataStoreGetCrossSectionDataStore ()
 
void MultiplyCrossSectionBy (G4double factor)
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
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 &)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4int operator== (const G4VProcess &right) const
 
G4int operator!= (const G4VProcess &right) const
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
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 void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (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 void StartTracking (G4Track *)
 
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
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Member Functions inherited from G4HadronicProcess
G4HadronicInteractionChooseHadronicInteraction (G4double kineticEnergy, G4Material *aMaterial, G4Element *anElement)
 
G4NucleusGetTargetNucleusPointer ()
 
void DumpState (const G4Track &, const G4String &, G4ExceptionDescription &)
 
const G4EnergyRangeManagerGetEnergyRangeManager () const
 
void SetEnergyRangeManager (const G4EnergyRangeManager &value)
 
G4HadronicInteractionGetHadronicInteraction () const
 
G4double GetLastCrossSection ()
 
void FillResult (G4HadFinalState *aR, const G4Track &aT)
 
G4HadFinalStateCheckResult (const G4HadProjectile &thePro, const G4Nucleus &targetNucleus, G4HadFinalState *result) const
 
void CheckEnergyMomentumConservation (const G4Track &, const G4Nucleus &)
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- Protected Attributes inherited from G4HadronicProcess
G4HadProjectile thePro
 
G4ParticleChangetheTotalResult
 
G4int epReportLevel
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChangepParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 

Detailed Description

Definition at line 52 of file G4WHadronElasticProcess.hh.

Constructor & Destructor Documentation

◆ G4WHadronElasticProcess()

G4WHadronElasticProcess::G4WHadronElasticProcess ( const G4String procName = "hadElastic")

Definition at line 60 of file G4WHadronElasticProcess.cc.

63 theNeutron = G4Neutron::Neutron();
64 lowestEnergy = 1.*keV;
65 lowestEnergyNeutron = 1.e-6*eV;
66 G4HadronicDeprecate("G4WHadronElasticProcess");
67}
#define G4HadronicDeprecate(name)
@ fHadronElastic
void AddDataSet(G4VCrossSectionDataSet *aDataSet)
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104

◆ ~G4WHadronElasticProcess()

G4WHadronElasticProcess::~G4WHadronElasticProcess ( )
virtual

Definition at line 69 of file G4WHadronElasticProcess.cc.

70{}

Member Function Documentation

◆ Description()

void G4WHadronElasticProcess::Description ( ) const
virtual

Definition at line 72 of file G4WHadronElasticProcess.cc.

73{
74 char* dirName = getenv("G4PhysListDocDir");
75 if (dirName) {
76 std::ofstream outFile;
77 G4String outFileName = GetProcessName() + ".html";
78 G4String pathName = G4String(dirName) + "/" + outFileName;
79 outFile.open(pathName);
80 outFile << "<html>\n";
81 outFile << "<head>\n";
82
83 outFile << "<title>Description of G4WHadronElasticProcess</title>\n";
84 outFile << "</head>\n";
85 outFile << "<body>\n";
86
87 outFile << "G4WHadronElasticProcess handles the elastic scattering of\n"
88 << "hadrons by invoking one or more hadronic models and one or\n"
89 << "more hadronic cross sections.\n";
90
91 outFile << "</body>\n";
92 outFile << "</html>\n";
93 outFile.close();
94 }
95}
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379

◆ PostStepDoIt()

G4VParticleChange * G4WHadronElasticProcess::PostStepDoIt ( const G4Track aTrack,
const G4Step aStep 
)
virtual

! is not needed for models inheriting G4HadronElastic

Reimplemented from G4HadronicProcess.

Definition at line 99 of file G4WHadronElasticProcess.cc.

101{
104 G4double weight = track.GetWeight();
106
107 // For elastic scattering, _any_ result is considered an interaction
109
110 G4double kineticEnergy = track.GetKineticEnergy();
111 const G4DynamicParticle* dynParticle = track.GetDynamicParticle();
112 const G4ParticleDefinition* part = dynParticle->GetDefinition();
113
114 // NOTE: Very low energy scatters were causing numerical (FPE) errors
115 // in earlier releases; these limits have not been changed since.
116 if (part == theNeutron) {
117 if(kineticEnergy <= lowestEnergyNeutron) return theTotalResult;
118 } else if(kineticEnergy <= lowestEnergy) return theTotalResult;
119
120 G4Material* material = track.GetMaterial();
121 G4Nucleus* targNucleus = GetTargetNucleusPointer();
122
123 // Select element
124 G4Element* elm = 0;
125 try
126 {
127 elm = GetCrossSectionDataStore()->SampleZandA(dynParticle, material,
128 *targNucleus);
129 }
130 catch(G4HadronicException & aR)
131 {
133 DumpState(track,"SampleZandA",ed);
134 ed << " PostStepDoIt failed on element selection" << G4endl;
135 G4Exception("G4WHadronElasticProcess::PostStepDoIt", "had003",
136 FatalException, ed);
137 }
138 G4HadronicInteraction* hadi = 0;
139 try
140 {
141 hadi = ChooseHadronicInteraction( kineticEnergy, material, elm );
142 }
143 catch(G4HadronicException & aE)
144 {
146 ed << "Target element "<< elm->GetName()<<" Z= "
147 << targNucleus->GetZ_asInt() << " A= "
148 << targNucleus->GetA_asInt() << G4endl;
149 DumpState(track,"ChooseHadronicInteraction",ed);
150 ed << " No HadronicInteraction found out" << G4endl;
151 G4Exception("G4WHadronElasticProcess::PostStepDoIt", "had005",
152 FatalException, ed);
153 }
154
155 size_t idx = track.GetMaterialCutsCouple()->GetIndex();
157 ->GetEnergyCutsVector(3)))[idx];
158 hadi->SetRecoilEnergyThreshold(tcut);
159
160 // Initialize the hadronic projectile from the track
161 // G4cout << "track " << track.GetDynamicParticle()->Get4Momentum()<<G4endl;
162 G4HadProjectile theProj(track);
163 if(verboseLevel>1) {
164 G4cout << "G4WHadronElasticProcess::PostStepDoIt for "
165 << part->GetParticleName()
166 << " in " << material->GetName()
167 << " Target Z= " << targNucleus->GetZ_asInt()
168 << " A= " << targNucleus->GetA_asInt() << G4endl;
169 }
170
171 G4HadFinalState* result = 0;
172 try
173 {
174 result = hadi->ApplyYourself( theProj, *targNucleus);
175 }
176 catch(G4HadronicException aR)
177 {
179 ed << "Call for " << hadi->GetModelName() << G4endl;
180 ed << "Target element "<< elm->GetName()<<" Z= "
181 << targNucleus->GetZ_asInt()
182 << " A= " << targNucleus->GetA_asInt() << G4endl;
183 DumpState(track,"ApplyYourself",ed);
184 ed << " ApplyYourself failed" << G4endl;
185 G4Exception("G4WHadronElasticProcess::PostStepDoIt", "had006",
186 FatalException, ed);
187 }
188
189 // Check the result for catastrophic energy non-conservation
190 // cannot be applied because is not guranteed that recoil
191 // nucleus is created
192 // result = CheckResult(theProj, targNucleus, result);
193
194 // directions
195 G4ThreeVector indir = track.GetMomentumDirection();
196 G4double phi = CLHEP::twopi*G4UniformRand();
197 G4ThreeVector it(0., 0., 1.);
198 G4ThreeVector outdir = result->GetMomentumChange();
199
200 if(verboseLevel>1) {
201 G4cout << "Efin= " << result->GetEnergyChange()
202 << " de= " << result->GetLocalEnergyDeposit()
203 << " nsec= " << result->GetNumberOfSecondaries()
204 << " dir= " << outdir
205 << G4endl;
206 }
207
208 // energies
209 G4double edep = result->GetLocalEnergyDeposit();
210 G4double efinal = result->GetEnergyChange();
211 if(efinal < 0.0) { efinal = 0.0; }
212 if(edep < 0.0) { edep = 0.0; }
213
214 // NOTE: Very low energy scatters were causing numerical (FPE) errors
215 // in earlier releases; these limits have not been changed since.
216 if((part == theNeutron && efinal <= lowestEnergyNeutron) ||
217 (part != theNeutron && efinal <= lowestEnergy)) {
218 edep += efinal;
219 efinal = 0.0;
220 }
221
222 // primary change
224
225 G4TrackStatus status = track.GetTrackStatus();
226 if(efinal > 0.0) {
227 outdir.rotate(phi, it);
228 outdir.rotateUz(indir);
230 } else {
231 if(part->GetProcessManager()->GetAtRestProcessVector()->size() > 0)
232 { status = fStopButAlive; }
233 else { status = fStopAndKill; }
235 }
236
237 //G4cout << "Efinal= " << efinal << " TrackStatus= " << status << G4endl;
238
240
241 // recoil
242 if(result->GetNumberOfSecondaries() > 0) {
243 G4DynamicParticle* p = result->GetSecondary(0)->GetParticle();
244
245 if(p->GetKineticEnergy() > lowestEnergy) {
248 // G4cout << "recoil " << pdir << G4endl;
249 //!! is not needed for models inheriting G4HadronElastic
250 pdir.rotate(phi, it);
251 pdir.rotateUz(indir);
252 // G4cout << "recoil rotated " << pdir << G4endl;
253 p->SetMomentumDirection(pdir);
254
255 // in elastic scattering time and weight are not changed
256 G4Track* t = new G4Track(p, track.GetGlobalTime(),
257 track.GetPosition());
258 t->SetWeight(weight);
259 t->SetTouchableHandle(track.GetTouchableHandle());
261
262 } else {
263 edep += p->GetKineticEnergy();
264 delete p;
265 }
266 }
269 result->Clear();
270
271 return theTotalResult;
272}
@ FatalException
G4TrackStatus
@ fStopAndKill
@ fStopButAlive
double G4double
Definition: G4Types.hh:64
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
Hep3Vector & rotate(double, const Hep3Vector &)
Definition: ThreeVectorR.cc:28
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
G4double GetLocalEnergyDeposit() const
G4int GetNumberOfSecondaries() const
const G4ThreeVector & GetMomentumChange() const
G4HadSecondary * GetSecondary(size_t i)
G4DynamicParticle * GetParticle()
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)=0
const G4String & GetModelName() const
void SetRecoilEnergyThreshold(G4double val)
G4Nucleus * GetTargetNucleusPointer()
G4ParticleChange * theTotalResult
G4CrossSectionDataStore * GetCrossSectionDataStore()
G4HadronicInteraction * ChooseHadronicInteraction(G4double kineticEnergy, G4Material *aMaterial, G4Element *anElement)
void DumpState(const G4Track &, const G4String &, G4ExceptionDescription &)
const G4String & GetName() const
Definition: G4Material.hh:177
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:576
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
G4int size() const
const std::vector< G4double > * GetEnergyCutsVector(size_t pcIdx) const
static G4ProductionCutsTable * GetProductionCutsTable()
void SetWeight(G4double aValue)
void SetTouchableHandle(const G4TouchableHandle &apValue)
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:418
G4int verboseLevel
Definition: G4VProcess.hh:368
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76

◆ SetLowestEnergy()

void G4WHadronElasticProcess::SetLowestEnergy ( G4double  val)
inline

Definition at line 82 of file G4WHadronElasticProcess.hh.

83{
84 lowestEnergy = val;
85}

◆ SetLowestEnergyNeutron()

void G4WHadronElasticProcess::SetLowestEnergyNeutron ( G4double  val)
inline

Definition at line 88 of file G4WHadronElasticProcess.hh.

89{
90 lowestEnergyNeutron = val;
91}

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