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

#include <G4ContinuousGainOfEnergy.hh>

Inheritance diagram for G4ContinuousGainOfEnergy:

Public Member Functions

 G4ContinuousGainOfEnergy (const G4String &name="EnergyGain", G4ProcessType type=fElectromagnetic)
 ~G4ContinuousGainOfEnergy () override
G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &) override
void SetLossFluctuations (G4bool val)
void SetDirectEnergyLossProcess (G4VEnergyLossProcess *aProcess)
void SetDirectParticle (G4ParticleDefinition *p)
void ProcessDescription (std::ostream &) const override
void DumpInfo () const override
 G4ContinuousGainOfEnergy (G4ContinuousGainOfEnergy &)=delete
G4ContinuousGainOfEnergyoperator= (const G4ContinuousGainOfEnergy &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 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
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
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 &)

Protected Member Functions

G4double GetContinuousStepLimit (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety) override
Protected Member Functions inherited from G4VContinuousProcess
void SetGPILSelection (G4GPILSelection selection)
G4GPILSelection GetGPILSelection () const
Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
void ClearNumberOfInteractionLengthLeft ()

Additional Inherited Members

Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
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 51 of file G4ContinuousGainOfEnergy.hh.

Constructor & Destructor Documentation

◆ G4ContinuousGainOfEnergy() [1/2]

G4ContinuousGainOfEnergy::G4ContinuousGainOfEnergy ( const G4String & name = "EnergyGain",
G4ProcessType type = fElectromagnetic )
explicit

Definition at line 45 of file G4ContinuousGainOfEnergy.cc.

47 : G4VContinuousProcess(name, type)
48{}
G4VContinuousProcess(const G4String &aName, G4ProcessType aType=fNotDefined)

Referenced by G4ContinuousGainOfEnergy(), and operator=().

◆ ~G4ContinuousGainOfEnergy()

G4ContinuousGainOfEnergy::~G4ContinuousGainOfEnergy ( )
override

Definition at line 51 of file G4ContinuousGainOfEnergy.cc.

51{}

◆ G4ContinuousGainOfEnergy() [2/2]

G4ContinuousGainOfEnergy::G4ContinuousGainOfEnergy ( G4ContinuousGainOfEnergy & )
delete

Member Function Documentation

◆ AlongStepDoIt()

G4VParticleChange * G4ContinuousGainOfEnergy::AlongStepDoIt ( const G4Track & track,
const G4Step & step )
overridevirtual

Reimplemented from G4VContinuousProcess.

Definition at line 73 of file G4ContinuousGainOfEnergy.cc.

75{
76 // Caution in this method the step length should be the true step length
77 // A problem is that this is computed by the multiple scattering that does
78 // not know the energy at the end of the adjoint step. This energy is used
79 // during the forward sim. Nothing we can really do against that at this
80 // time. This is inherent to the MS method
81
82 aParticleChange.Initialize(track);
83
84 // Get the actual (true) Step length
85 G4double length = step.GetStepLength();
86 G4double degain = 0.0;
87
88 // Compute this for weight change after continuous energy loss
89 G4double DEDX_before =
90 fDirectEnergyLossProcess->GetDEDX(fPreStepKinEnergy, fCurrentCouple);
91
92 // For the fluctuation we generate a new dynamic particle with energy
93 // = preEnergy+egain and then compute the fluctuation given in the direct
94 // case.
95 G4DynamicParticle* dynParticle = new G4DynamicParticle();
96 *dynParticle = *(track.GetDynamicParticle());
97 dynParticle->SetDefinition(fDirectPartDef);
98 G4double Tkin = dynParticle->GetKineticEnergy();
99
100 G4double dlength = length;
101 if(Tkin != fPreStepKinEnergy && fIsIon)
102 {
103 G4double chargeSqRatio = fCurrentModel->GetChargeSquareRatio(
104 fDirectPartDef, fCurrentMaterial, Tkin);
105 fDirectEnergyLossProcess->SetDynamicMassCharge(fMassRatio, chargeSqRatio);
106 }
107
108 G4double r = fDirectEnergyLossProcess->GetRange(Tkin, fCurrentCouple);
109 if(dlength <= fLinLossLimit * r)
110 {
111 degain = DEDX_before * dlength;
112 }
113 else
114 {
115 G4double x = r + dlength;
116 G4double E = fDirectEnergyLossProcess->GetKineticEnergy(x, fCurrentCouple);
117 if(fIsIon)
118 {
119 G4double chargeSqRatio = fCurrentModel->GetChargeSquareRatio(
120 fDirectPartDef, fCurrentMaterial, E);
121 fDirectEnergyLossProcess->SetDynamicMassCharge(fMassRatio, chargeSqRatio);
122 G4double x1 = fDirectEnergyLossProcess->GetRange(E, fCurrentCouple);
123
124 G4int ii = 0;
125 constexpr G4int iimax = 100;
126 while(std::abs(x - x1) > 0.01 * x)
127 {
128 E = fDirectEnergyLossProcess->GetKineticEnergy(x, fCurrentCouple);
129 chargeSqRatio = fCurrentModel->GetChargeSquareRatio(
130 fDirectPartDef, fCurrentMaterial, E);
131 fDirectEnergyLossProcess->SetDynamicMassCharge(fMassRatio,
132 chargeSqRatio);
133 x1 = fDirectEnergyLossProcess->GetRange(E, fCurrentCouple);
134 ++ii;
135 if(ii >= iimax)
136 {
137 break;
138 }
139 }
140 }
141
142 degain = E - Tkin;
143 }
144 G4double tmax = fCurrentModel->MaxSecondaryKinEnergy(dynParticle);
145 fCurrentTcut = std::min(fCurrentTcut, tmax);
146
147 dynParticle->SetKineticEnergy(Tkin + degain);
148
149 // Corrections, which cannot be tabulated for ions
150 fCurrentModel->CorrectionsAlongStep(fCurrentMaterial, fDirectPartDef, Tkin,
151 fCurrentTcut, dlength, degain);
152
153 // Sample fluctuations
154 G4double deltaE = 0.;
155 if(fLossFluctuationFlag)
156 {
157 deltaE = fCurrentModel->GetModelOfFluctuations()->SampleFluctuations(
158 fCurrentCouple, dynParticle, fCurrentTcut, tmax, dlength, degain)
159 - degain;
160 }
161
162 G4double egain = degain + deltaE;
163 if(egain <= 0.)
164 egain = degain;
165 Tkin += egain;
166 dynParticle->SetKineticEnergy(Tkin);
167
168 delete dynParticle;
169
170 if(fIsIon)
171 {
172 G4double chargeSqRatio = fCurrentModel->GetChargeSquareRatio(
173 fDirectPartDef, fCurrentMaterial, Tkin);
174 fDirectEnergyLossProcess->SetDynamicMassCharge(fMassRatio, chargeSqRatio);
175 }
176
177 G4double DEDX_after = fDirectEnergyLossProcess->GetDEDX(Tkin, fCurrentCouple);
178 G4double weight_correction = DEDX_after / DEDX_before;
179
180 aParticleChange.ProposeEnergy(Tkin);
181
182 // Caution!!! It is important to select the weight of the post_step_point
183 // as the current weight and not the weight of the track, as the weight of
184 // the track is changed after having applied all the along_step_do_it.
185
186 G4double new_weight =
187 weight_correction * step.GetPostStepPoint()->GetWeight();
188 aParticleChange.SetParentWeightByProcess(false);
189 aParticleChange.ProposeParentWeight(new_weight);
190
191 return &aParticleChange;
192}
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
G4double GetKineticEnergy() const
void SetKineticEnergy(G4double aEnergy)
G4double GetWeight() const
G4double GetStepLength() const
G4StepPoint * GetPostStepPoint() const
const G4DynamicParticle * GetDynamicParticle() const
G4ParticleChange aParticleChange

◆ DumpInfo()

void G4ContinuousGainOfEnergy::DumpInfo ( ) const
inlineoverridevirtual

Reimplemented from G4VProcess.

Definition at line 71 of file G4ContinuousGainOfEnergy.hh.

G4GLOB_DLL std::ostream G4cout
void ProcessDescription(std::ostream &) const override

◆ GetContinuousStepLimit()

G4double G4ContinuousGainOfEnergy::GetContinuousStepLimit ( const G4Track & track,
G4double previousStepSize,
G4double currentMinimumStep,
G4double & currentSafety )
overrideprotectedvirtual

Implements G4VContinuousProcess.

Definition at line 203 of file G4ContinuousGainOfEnergy.cc.

206{
207 DefineMaterial(track.GetMaterialCutsCouple());
208
209 fPreStepKinEnergy = track.GetKineticEnergy();
210 fCurrentModel = fDirectEnergyLossProcess->SelectModelForMaterial(
211 track.GetKineticEnergy() * fMassRatio, fCurrentCoupleIndex);
212 G4double emax_model = fCurrentModel->HighEnergyLimit();
213 G4double preStepChargeSqRatio = 0.;
214 if(fIsIon)
215 {
216 G4double chargeSqRatio = fCurrentModel->GetChargeSquareRatio(
217 fDirectPartDef, fCurrentMaterial, fPreStepKinEnergy);
218 preStepChargeSqRatio = chargeSqRatio;
219 fDirectEnergyLossProcess->SetDynamicMassCharge(fMassRatio,
220 preStepChargeSqRatio);
221 }
222
223 G4double maxE = 1.1 * fPreStepKinEnergy;
224
225 if(fPreStepKinEnergy < fCurrentTcut)
226 maxE = std::min(fCurrentTcut, maxE);
227
228 maxE = std::min(emax_model * 1.001, maxE);
229
230 G4double preStepRange =
231 fDirectEnergyLossProcess->GetRange(fPreStepKinEnergy, fCurrentCouple);
232
233 if(fIsIon)
234 {
235 G4double chargeSqRatioAtEmax = fCurrentModel->GetChargeSquareRatio(
236 fDirectPartDef, fCurrentMaterial, maxE);
237 fDirectEnergyLossProcess->SetDynamicMassCharge(fMassRatio,
238 chargeSqRatioAtEmax);
239 }
240
241 G4double r1 = fDirectEnergyLossProcess->GetRange(maxE, fCurrentCouple);
242
243 if(fIsIon)
244 fDirectEnergyLossProcess->SetDynamicMassCharge(fMassRatio,
245 preStepChargeSqRatio);
246
247 return std::max(r1 - preStepRange, 0.001 * mm);
248}
G4double GetKineticEnergy() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const

◆ operator=()

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

◆ ProcessDescription()

void G4ContinuousGainOfEnergy::ProcessDescription ( std::ostream & out) const
overridevirtual

Reimplemented from G4VProcess.

Definition at line 54 of file G4ContinuousGainOfEnergy.cc.

55{
56 out << "Continuous process acting on adjoint particles to compute the "
57 "continuous gain of energy of charged particles when they are "
58 "tracked back.\n";
59}

Referenced by DumpInfo().

◆ SetDirectEnergyLossProcess()

void G4ContinuousGainOfEnergy::SetDirectEnergyLossProcess ( G4VEnergyLossProcess * aProcess)
inline

Definition at line 63 of file G4ContinuousGainOfEnergy.hh.

64 {
65 fDirectEnergyLossProcess = aProcess;
66 };

◆ SetDirectParticle()

void G4ContinuousGainOfEnergy::SetDirectParticle ( G4ParticleDefinition * p)

Definition at line 62 of file G4ContinuousGainOfEnergy.cc.

63{
64 fDirectPartDef = p;
65 if(fDirectPartDef->GetParticleType() == "nucleus")
66 {
67 fIsIon = true;
68 fMassRatio = proton_mass_c2 / fDirectPartDef->GetPDGMass();
69 }
70}

◆ SetLossFluctuations()

void G4ContinuousGainOfEnergy::SetLossFluctuations ( G4bool val)

Definition at line 195 of file G4ContinuousGainOfEnergy.cc.

196{
197 if(val && !fLossFluctuationArePossible)
198 return;
199 fLossFluctuationFlag = val;
200}

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