Geant4 10.7.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)
 
virtual ~G4ContinuousGainOfEnergy ()
 
void PreparePhysicsTable (const G4ParticleDefinition &)
 
void BuildPhysicsTable (const G4ParticleDefinition &)
 
G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
void SetLossFluctuations (G4bool val)
 
void SetIsIntegral (G4bool val)
 
void SetDirectEnergyLossProcess (G4VEnergyLossProcess *aProcess)
 
void SetDirectParticle (G4ParticleDefinition *p)
 
- 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 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 &)
 

Protected Member Functions

virtual G4double GetContinuousStepLimit (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)
 
- 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 ()
 

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 68 of file G4ContinuousGainOfEnergy.hh.

Constructor & Destructor Documentation

◆ G4ContinuousGainOfEnergy()

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

Definition at line 44 of file G4ContinuousGainOfEnergy.cc.

45 : G4VContinuousProcess(name, type)
46{
47
48
49 linLossLimit=0.05;
50 lossFluctuationArePossible =true;
51 lossFluctuationFlag=true;
52 is_integral = false;
53
54 //Will be properly set in SetDirectParticle()
55 IsIon=false;
56 massRatio =1.;
57 chargeSqRatio=1.;
58 preStepChargeSqRatio=1.;
59
60 //Some initialization
61 currentCoupleIndex=9999999;
62 currentCutInRange=0.;
63 currentMaterialIndex=9999999;
64 currentTcut=0.;
65 preStepKinEnergy=0.;
66 preStepRange=0.;
67 preStepScaledKinEnergy=0.;
68
69 currentCouple=0;
70}

◆ ~G4ContinuousGainOfEnergy()

G4ContinuousGainOfEnergy::~G4ContinuousGainOfEnergy ( )
virtual

Definition at line 74 of file G4ContinuousGainOfEnergy.cc.

75{
76
77}

Member Function Documentation

◆ AlongStepDoIt()

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

Reimplemented from G4VContinuousProcess.

Definition at line 114 of file G4ContinuousGainOfEnergy.cc.

116{
117
118 //Caution in this method the step length should be the true step length
119 // A problem is that this is compute by the multiple scattering that does not know the energy at the end of the adjoint step. This energy is used during the
120 //Forward sim. Nothing we can really do against that at this time. This is inherent to the MS method
121 //
122
123
124
126
127 // Get the actual (true) Step length
128 //----------------------------------
129 G4double length = step.GetStepLength();
130 G4double degain = 0.0;
131
132
133
134 // Compute this for weight change after continuous energy loss
135 //-------------------------------------------------------------
136 G4double DEDX_before = theDirectEnergyLossProcess->GetDEDX(preStepKinEnergy, currentCouple);
137
138
139
140 // For the fluctuation we generate a new dynamic particle with energy =preEnergy+egain
141 // and then compute the fluctuation given in the direct case.
142 //-----------------------------------------------------------------------
143 G4DynamicParticle* dynParticle = new G4DynamicParticle();
144 *dynParticle = *(track.GetDynamicParticle());
145 dynParticle->SetDefinition(theDirectPartDef);
146 G4double Tkin = dynParticle->GetKineticEnergy();
147
148
149 size_t n=1;
150 if (is_integral ) n=10;
151 n=1;
152 G4double dlength= length/n;
153 for (size_t i=0;i<n;i++) {
154 if (Tkin != preStepKinEnergy && IsIon) {
155 chargeSqRatio = currentModel->GetChargeSquareRatio(theDirectPartDef,currentMaterial,Tkin);
156 theDirectEnergyLossProcess->SetDynamicMassCharge(massRatio,chargeSqRatio);
157
158 }
159
160 G4double r = theDirectEnergyLossProcess->GetRange(Tkin, currentCouple);
161 if( dlength <= linLossLimit * r ) {
162 degain = DEDX_before*dlength;
163 }
164 else {
165 G4double x = r + dlength;
166 //degain = theDirectEnergyLossProcess->GetKineticEnergy(x,currentCouple) - theDirectEnergyLossProcess->GetKineticEnergy(r,currentCouple);
167 G4double E = theDirectEnergyLossProcess->GetKineticEnergy(x,currentCouple);
168 if (IsIon){
169 chargeSqRatio = currentModel->GetChargeSquareRatio(theDirectPartDef,currentMaterial,E);
170 theDirectEnergyLossProcess->SetDynamicMassCharge(massRatio,chargeSqRatio);
171 G4double x1= theDirectEnergyLossProcess->GetRange(E, currentCouple);
172
173 // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
174 G4int ii=0;
175 const G4int iimax = 100;
176 while (std::abs(x-x1)>0.01*x) {
177 E = theDirectEnergyLossProcess->GetKineticEnergy(x,currentCouple);
178 chargeSqRatio = currentModel->GetChargeSquareRatio(theDirectPartDef,currentMaterial,E);
179 theDirectEnergyLossProcess->SetDynamicMassCharge(massRatio,chargeSqRatio);
180 x1= theDirectEnergyLossProcess->GetRange(E, currentCouple);
181 ++ii;
182 if(ii >= iimax) { break; }
183 }
184 }
185
186 degain=E-Tkin;
187
188
189
190 }
191 //G4cout<<degain<<G4endl;
192 G4double tmax = currentModel->MaxSecondaryKinEnergy(dynParticle);
193 tmax = std::min(tmax,currentTcut);
194
195
196 dynParticle->SetKineticEnergy(Tkin+degain);
197
198 // Corrections, which cannot be tabulated for ions
199 //----------------------------------------
200 G4double esecdep=0;//not used in most models
201 currentModel->CorrectionsAlongStep(currentCouple, dynParticle, degain,esecdep, dlength);
202
203 // Sample fluctuations
204 //-------------------
205
206
207 G4double deltaE =0.;
208 if (lossFluctuationFlag ) {
209 deltaE = currentModel->GetModelOfFluctuations()->
210 SampleFluctuations(currentCouple,dynParticle,tmax,dlength,degain)-degain;
211 }
212
213 G4double egain=degain+deltaE;
214 if (egain <=0) egain=degain;
215 Tkin+=egain;
216 dynParticle->SetKineticEnergy(Tkin);
217 }
218
219
220
221
222
223 delete dynParticle;
224
225 if (IsIon){
226 chargeSqRatio = currentModel->GetChargeSquareRatio(theDirectPartDef,currentMaterial,Tkin);
227 theDirectEnergyLossProcess->SetDynamicMassCharge(massRatio,chargeSqRatio);
228
229 }
230
231 G4double DEDX_after = theDirectEnergyLossProcess->GetDEDX(Tkin, currentCouple);
232
233
234 G4double weight_correction=DEDX_after/DEDX_before;
235
236
238
239
240 //Caution!!!
241 // It is important to select the weight of the post_step_point
242 // as the current weight and not the weight of the track, as t
243 // the weight of the track is changed after having applied all
244 // the along_step_do_it.
245
246 // G4double new_weight=weight_correction*track.GetWeight(); //old
247 G4double new_weight=weight_correction*step.GetPostStepPoint()->GetWeight();
250
251
252 return &aParticleChange;
253
254}
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
G4double GetKineticEnergy() const
void SetKineticEnergy(G4double aEnergy)
void ProposeEnergy(G4double finalEnergy)
virtual void Initialize(const G4Track &)
G4double GetWeight() const
G4double GetStepLength() const
G4StepPoint * GetPostStepPoint() const
const G4DynamicParticle * GetDynamicParticle() const
virtual void CorrectionsAlongStep(const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
Definition: G4VEmModel.cc:408
G4VEmFluctuationModel * GetModelOfFluctuations()
Definition: G4VEmModel.hh:604
virtual G4double GetChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:391
G4double MaxSecondaryKinEnergy(const G4DynamicParticle *dynParticle)
Definition: G4VEmModel.hh:510
G4double GetKineticEnergy(G4double range, const G4MaterialCutsCouple *)
void SetDynamicMassCharge(G4double massratio, G4double charge2ratio)
G4double GetDEDX(G4double kineticEnergy, const G4MaterialCutsCouple *)
G4double GetRange(G4double kineticEnergy, const G4MaterialCutsCouple *)
void SetParentWeightByProcess(G4bool)
void ProposeParentWeight(G4double finalWeight)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:327

◆ BuildPhysicsTable()

void G4ContinuousGainOfEnergy::BuildPhysicsTable ( const G4ParticleDefinition )
virtual

Reimplemented from G4VProcess.

Definition at line 91 of file G4ContinuousGainOfEnergy.cc.

92{//theDirectEnergyLossProcess->BuildPhysicsTable(part);
93;
94}

◆ GetContinuousStepLimit()

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

Implements G4VContinuousProcess.

Definition at line 267 of file G4ContinuousGainOfEnergy.cc.

269{
270 G4double x = DBL_MAX;
271 x=.1*mm;
272
273
274 DefineMaterial(track.GetMaterialCutsCouple());
275
276 preStepKinEnergy = track.GetKineticEnergy();
277 preStepScaledKinEnergy = track.GetKineticEnergy()*massRatio;
278 currentModel = theDirectEnergyLossProcess->SelectModelForMaterial(preStepScaledKinEnergy,currentCoupleIndex);
279 G4double emax_model=currentModel->HighEnergyLimit();
280 if (IsIon) {
281 chargeSqRatio = currentModel->GetChargeSquareRatio(theDirectPartDef,currentMaterial,preStepKinEnergy);
282 preStepChargeSqRatio = chargeSqRatio;
283 theDirectEnergyLossProcess->SetDynamicMassCharge(massRatio,preStepChargeSqRatio);
284 }
285
286
287 G4double maxE =1.1*preStepKinEnergy;
288 /*if (preStepKinEnergy< 0.05*MeV) maxE =2.*preStepKinEnergy;
289 else if (preStepKinEnergy< 0.1*MeV) maxE =1.5*preStepKinEnergy;
290 else if (preStepKinEnergy< 0.5*MeV) maxE =1.25*preStepKinEnergy;*/
291
292 if (preStepKinEnergy < currentTcut) maxE = std::min(currentTcut,maxE);
293
294 maxE=std::min(emax_model*1.001,maxE);
295
296 preStepRange = theDirectEnergyLossProcess->GetRange(preStepKinEnergy, currentCouple);
297
298 if (IsIon) {
299 G4double chargeSqRatioAtEmax = currentModel->GetChargeSquareRatio(theDirectPartDef,currentMaterial,maxE);
300 theDirectEnergyLossProcess->SetDynamicMassCharge(massRatio,chargeSqRatioAtEmax);
301 }
302
303 G4double r1 = theDirectEnergyLossProcess->GetRange(maxE, currentCouple);
304
305 if (IsIon) theDirectEnergyLossProcess->SetDynamicMassCharge(massRatio,preStepChargeSqRatio);
306
307
308
309 x=r1-preStepRange;
310 x=std::max(r1-preStepRange,0.001*mm);
311
312 return x;
313
314
315}
G4double GetKineticEnergy() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:645
G4VEmModel * SelectModelForMaterial(G4double kinEnergy, size_t &idx) const
#define DBL_MAX
Definition: templates.hh:62

◆ PreparePhysicsTable()

void G4ContinuousGainOfEnergy::PreparePhysicsTable ( const G4ParticleDefinition )
virtual

Reimplemented from G4VProcess.

Definition at line 81 of file G4ContinuousGainOfEnergy.cc.

83{//theDirectEnergyLossProcess->PreparePhysicsTable(part);
84
85;
86}

◆ SetDirectEnergyLossProcess()

void G4ContinuousGainOfEnergy::SetDirectEnergyLossProcess ( G4VEnergyLossProcess aProcess)
inline

Definition at line 111 of file G4ContinuousGainOfEnergy.hh.

111{theDirectEnergyLossProcess=aProcess;};

◆ SetDirectParticle()

void G4ContinuousGainOfEnergy::SetDirectParticle ( G4ParticleDefinition p)

Definition at line 98 of file G4ContinuousGainOfEnergy.cc.

99{theDirectPartDef=p;
100 if (theDirectPartDef->GetParticleType()== "nucleus") {
101 IsIon=true;
102 massRatio = proton_mass_c2/theDirectPartDef->GetPDGMass();
103 G4double q=theDirectPartDef->GetPDGCharge();
104 chargeSqRatio=q*q;
105
106
107 }
108
109}
const G4String & GetParticleType() const
G4double GetPDGCharge() const

◆ SetIsIntegral()

void G4ContinuousGainOfEnergy::SetIsIntegral ( G4bool  val)
inline

Definition at line 109 of file G4ContinuousGainOfEnergy.hh.

109{is_integral= val;}

◆ SetLossFluctuations()

void G4ContinuousGainOfEnergy::SetLossFluctuations ( G4bool  val)

Definition at line 257 of file G4ContinuousGainOfEnergy.cc.

258{
259 if(val && !lossFluctuationArePossible) return;
260 lossFluctuationFlag = val;
261}

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