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

#include <G4ErrorEnergyLoss.hh>

+ Inheritance diagram for G4ErrorEnergyLoss:

Public Member Functions

 G4ErrorEnergyLoss (const G4String &processName="G4ErrorEnergyLoss", G4ProcessType type=fElectromagnetic)
 
 ~G4ErrorEnergyLoss () override
 
G4bool IsApplicable (const G4ParticleDefinition &aParticleType) override
 
G4double GetContinuousStepLimit (const G4Track &aTrack, G4double, G4double currentMinimumStep, G4double &) override
 
G4VParticleChangeAlongStepDoIt (const G4Track &aTrack, const G4Step &aStep) override
 
G4double GetStepLimit () const
 
void SetStepLimit (G4double val)
 
 G4ErrorEnergyLoss (G4ErrorEnergyLoss &)=delete
 
G4ErrorEnergyLossoperator= (const G4ErrorEnergyLoss &right)=delete
 
- Public Member Functions inherited from G4VContinuousProcess
 G4VContinuousProcess (const G4String &aName, G4ProcessType aType=fNotDefined)
 
 G4VContinuousProcess (G4VContinuousProcess &)
 
virtual ~G4VContinuousProcess ()
 
G4VContinuousProcessoperator= (const G4VContinuousProcess &)=delete
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &, G4double, G4ForceCondition *)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
- 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
 
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 const G4VProcessGetCreatorProcess () const
 
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
 
virtual void ProcessDescription (std::ostream &outfile) 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 G4VContinuousProcess
virtual G4double GetContinuousStepLimit (const G4Track &aTrack, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)=0
 
void SetGPILSelection (G4GPILSelection selection)
 
G4GPILSelection GetGPILSelection () const
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- 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 47 of file G4ErrorEnergyLoss.hh.

Constructor & Destructor Documentation

◆ G4ErrorEnergyLoss() [1/2]

G4ErrorEnergyLoss::G4ErrorEnergyLoss ( const G4String processName = "G4ErrorEnergyLoss",
G4ProcessType  type = fElectromagnetic 
)
explicit

Definition at line 32 of file G4ErrorEnergyLoss.cc.

34 : G4VContinuousProcess(processName, type)
35{
36 if (verboseLevel>2) {
37 G4cout << GetProcessName() << " is created " << G4endl;
38 }
39
40 theELossForExtrapolator = new G4EnergyLossForExtrapolator;
41 theStepLimit = 1.*CLHEP::mm;
42}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4int verboseLevel
Definition: G4VProcess.hh:360
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386

◆ ~G4ErrorEnergyLoss()

G4ErrorEnergyLoss::~G4ErrorEnergyLoss ( )
override

Definition at line 49 of file G4ErrorEnergyLoss.cc.

50{
51 delete theELossForExtrapolator;
52}

◆ G4ErrorEnergyLoss() [2/2]

G4ErrorEnergyLoss::G4ErrorEnergyLoss ( G4ErrorEnergyLoss )
delete

Member Function Documentation

◆ AlongStepDoIt()

G4VParticleChange * G4ErrorEnergyLoss::AlongStepDoIt ( const G4Track aTrack,
const G4Step aStep 
)
overridevirtual

Reimplemented from G4VContinuousProcess.

Definition at line 63 of file G4ErrorEnergyLoss.cc.

64{
66
68
69 G4double kinEnergyStart = aTrack.GetKineticEnergy();
70 G4double step_length = aStep.GetStepLength();
71
72 const G4Material* aMaterial = aTrack.GetMaterial();
73 const G4ParticleDefinition* aParticleDef = aTrack.GetDynamicParticle()->GetDefinition();
74 G4double kinEnergyEnd = kinEnergyStart;
75
76 // backward - energy increased
77 if( g4edata->GetMode() == G4ErrorMode(G4ErrorMode_PropBackwards) ) {
78 kinEnergyEnd = theELossForExtrapolator->EnergyBeforeStep( kinEnergyStart,
79 step_length,
80 aMaterial,
81 aParticleDef );
82 G4double kinEnergyHalfStep = (kinEnergyStart + kinEnergyEnd) * 0.5;
83
84#ifdef G4VERBOSE
86 G4cout << " G4ErrorEnergyLoss FWD end " << kinEnergyEnd
87 << " halfstep " << kinEnergyHalfStep << G4endl;
88#endif
89
90 //--- rescale to energy lost at 1/2 step
91 kinEnergyEnd = theELossForExtrapolator->EnergyBeforeStep( kinEnergyHalfStep,
92 step_length,
93 aMaterial,
94 aParticleDef );
95 kinEnergyEnd = kinEnergyStart - (kinEnergyHalfStep - kinEnergyEnd );
96
97 // forward - energy decreased
98 } else {
99 kinEnergyEnd = theELossForExtrapolator->EnergyAfterStep( kinEnergyStart,
100 step_length,
101 aMaterial,
102 aParticleDef );
103 G4double kinEnergyHalfStep = (kinEnergyStart + kinEnergyEnd) * 0.5;
104#ifdef G4VERBOSE
106 G4cout << " G4ErrorEnergyLoss BCKD end " << kinEnergyEnd
107 << " halfstep " << kinEnergyHalfStep << G4endl;
108#endif
109
110 //--- rescale to energy lost at 1/2 step
111 kinEnergyEnd = theELossForExtrapolator->EnergyAfterStep( kinEnergyHalfStep,
112 step_length,
113 aMaterial,
114 aParticleDef );
115 kinEnergyEnd = kinEnergyStart - (kinEnergyHalfStep - kinEnergyEnd );
116 }
117
118 G4double edepo = kinEnergyEnd - kinEnergyStart;
119
120#ifdef G4VERBOSE
122 G4cout << "AlongStepDoIt Estart= " << kinEnergyStart << " Eend " << kinEnergyEnd
123 << " Ediff " << kinEnergyStart-kinEnergyEnd << " step= " << step_length
124 << " mate= " << aMaterial->GetName()
125 << " particle= " << aParticleDef->GetParticleName() << G4endl;
126#endif
127
131
132 aParticleChange.ProposeEnergy( kinEnergyEnd );
133
134 return &aParticleChange;
135}
@ G4ErrorMode_PropBackwards
double G4double
Definition: G4Types.hh:83
G4ParticleDefinition * GetDefinition() const
G4double EnergyBeforeStep(G4double kinEnergy, G4double step, const G4Material *, const G4ParticleDefinition *)
G4double EnergyAfterStep(G4double kinEnergy, G4double step, const G4Material *, const G4ParticleDefinition *)
static G4ErrorPropagatorData * GetErrorPropagatorData()
G4ErrorMode GetMode() const
const G4String & GetName() const
Definition: G4Material.hh:172
void Initialize(const G4Track &) override
void ProposeEnergy(G4double finalEnergy)
const G4String & GetParticleName() const
G4double GetStepLength() const
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
G4double GetKineticEnergy() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetNumberOfSecondaries(G4int totSecondaries)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:331

◆ GetContinuousStepLimit()

G4double G4ErrorEnergyLoss::GetContinuousStepLimit ( const G4Track aTrack,
G4double  ,
G4double  currentMinimumStep,
G4double  
)
overridevirtual

Implements G4VContinuousProcess.

Definition at line 139 of file G4ErrorEnergyLoss.cc.

141{
142 G4double ekin = aTrack.GetKineticEnergy();
143 const G4Material* mat = aTrack.GetMaterial();
144 const G4ParticleDefinition* part =
146 G4double range = theELossForExtrapolator->ComputeRange(ekin, part, mat);
147 G4double delta = std::max(range*theFractionLimit, theStepLimit);
148#ifdef G4VERBOSE
150 G4cout << " G4ErrorEnergyLoss: limiting Step " << delta
151 << " energy(GeV) " << ekin / CLHEP::GeV
152 << " for " << part->GetParticleName() << G4endl;
153 }
154#endif
155 return delta;
156}
G4double ComputeRange(G4double kinEnergy, const G4ParticleDefinition *, const G4Material *)

◆ GetStepLimit()

G4double G4ErrorEnergyLoss::GetStepLimit ( ) const
inline

Definition at line 69 of file G4ErrorEnergyLoss.hh.

69{ return theStepLimit; }

◆ IsApplicable()

G4bool G4ErrorEnergyLoss::IsApplicable ( const G4ParticleDefinition aParticleType)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 56 of file G4ErrorEnergyLoss.cc.

57{
58 return (aParticleType.GetPDGCharge() != 0);
59}
G4double GetPDGCharge() const

◆ operator=()

G4ErrorEnergyLoss & G4ErrorEnergyLoss::operator= ( const G4ErrorEnergyLoss right)
delete

◆ SetStepLimit()

void G4ErrorEnergyLoss::SetStepLimit ( G4double  val)
inline

Definition at line 70 of file G4ErrorEnergyLoss.hh.

70{ theStepLimit = val; }

Referenced by G4ErrorMessenger::SetNewValue().


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