Geant4  10.02.p03
G4ContinuousGainOfEnergy Class Reference

#include <G4ContinuousGainOfEnergy.hh>

Inheritance diagram for G4ContinuousGainOfEnergy:
Collaboration diagram for G4ContinuousGainOfEnergy:

Public Member Functions

 G4ContinuousGainOfEnergy (const G4String &name="EnergyGain", G4ProcessType type=fElectromagnetic)
 
virtual ~G4ContinuousGainOfEnergy ()
 
void PreparePhysicsTable (const G4ParticleDefinition &)
 
void BuildPhysicsTable (const G4ParticleDefinition &)
 
G4VParticleChange * AlongStepDoIt (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 &, G4ProcessType aType=fNotDefined)
 
 G4VContinuousProcess (G4VContinuousProcess &)
 
virtual ~G4VContinuousProcess ()
 
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 G4VParticleChange * AtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChange * PostStepDoIt (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
 
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 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
 
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
void SetGPILSelection (G4GPILSelection selection)
 
G4GPILSelection GetGPILSelection () const
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Private Member Functions

void DefineMaterial (const G4MaterialCutsCouple *couple)
 
void SetDynamicMassCharge (const G4Track &track, G4double energy)
 
 G4ContinuousGainOfEnergy (G4ContinuousGainOfEnergy &)
 
G4ContinuousGainOfEnergyoperator= (const G4ContinuousGainOfEnergy &right)
 

Private Attributes

const G4MaterialcurrentMaterial
 
const G4MaterialCutsCouplecurrentCouple
 
size_t currentMaterialIndex
 
size_t currentCoupleIndex
 
G4double currentTcut
 
G4double currentCutInRange
 
G4double preStepKinEnergy
 
G4double linLossLimit
 
G4bool lossFluctuationFlag
 
G4bool lossFluctuationArePossible
 
G4VEnergyLossProcesstheDirectEnergyLossProcess
 
G4ParticleDefinitiontheDirectPartDef
 
G4bool is_integral
 
G4bool IsIon
 
G4double massRatio
 
G4double chargeSqRatio
 
G4VEmModelcurrentModel
 
G4double preStepChargeSqRatio
 
G4double preStepScaledKinEnergy
 
G4double preStepRange
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChange * pParticleChange
 
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 69 of file G4ContinuousGainOfEnergy.hh.

Constructor & Destructor Documentation

◆ G4ContinuousGainOfEnergy() [1/2]

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

Definition at line 45 of file G4ContinuousGainOfEnergy.cc.

46  : G4VContinuousProcess(name, type)
47 {
48 
49 
50  linLossLimit=0.05;
53  is_integral = false;
54 
55  //Will be properly set in SetDirectParticle()
56  IsIon=false;
57  massRatio =1.;
58  chargeSqRatio=1.;
60 
61  //Some initialization
62  currentCoupleIndex=9999999;
64  currentMaterialIndex=9999999;
65  currentTcut=0.;
67  preStepRange=0.;
69 
70  currentCouple=0;
71 }
const G4MaterialCutsCouple * currentCouple
Here is the caller graph for this function:

◆ ~G4ContinuousGainOfEnergy()

G4ContinuousGainOfEnergy::~G4ContinuousGainOfEnergy ( )
virtual

Definition at line 75 of file G4ContinuousGainOfEnergy.cc.

76 {
77 
78 }

◆ G4ContinuousGainOfEnergy() [2/2]

G4ContinuousGainOfEnergy::G4ContinuousGainOfEnergy ( G4ContinuousGainOfEnergy )
private

Member Function Documentation

◆ AlongStepDoIt()

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

Reimplemented from G4VContinuousProcess.

Definition at line 115 of file G4ContinuousGainOfEnergy.cc.

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

◆ BuildPhysicsTable()

void G4ContinuousGainOfEnergy::BuildPhysicsTable ( const G4ParticleDefinition )
virtual

Reimplemented from G4VProcess.

Definition at line 92 of file G4ContinuousGainOfEnergy.cc.

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

◆ DefineMaterial()

void G4ContinuousGainOfEnergy::DefineMaterial ( const G4MaterialCutsCouple couple)
inlineprivate

Definition at line 173 of file G4ContinuousGainOfEnergy.hh.

175 {
176  if(couple != currentCouple) {
177  currentCouple = couple;
178  currentMaterial = couple->GetMaterial();
179  currentCoupleIndex = couple->GetIndex();
181 
182  size_t idx=1;
183  const std::vector<G4double>* aVec = G4ProductionCutsTable::GetProductionCutsTable()->GetEnergyCutsVector(idx);
186  //G4cout<<"Define Material"<<G4endl;
187  //if(!meanFreePath) ResetNumberOfInteractionLengthLeft();
188  }
189 }
const G4Material * GetMaterial() const
G4ProductionCuts * GetProductionCuts() const
size_t GetIndex() const
Definition: G4Material.hh:262
const std::vector< G4double > * GetEnergyCutsVector(size_t pcIdx) const
G4ParticleDefinition * theDirectPartDef
G4double GetProductionCut(G4int index) const
const G4String & GetParticleName() const
const G4MaterialCutsCouple * currentCouple
static G4ProductionCutsTable * GetProductionCutsTable()
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetContinuousStepLimit()

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

Implements G4VContinuousProcess.

Definition at line 268 of file G4ContinuousGainOfEnergy.cc.

270 {
271  G4double x = DBL_MAX;
272  x=.1*mm;
273 
274 
275  DefineMaterial(track.GetMaterialCutsCouple());
276 
277  preStepKinEnergy = track.GetKineticEnergy();
278  preStepScaledKinEnergy = track.GetKineticEnergy()*massRatio;
280  G4double emax_model=currentModel->HighEnergyLimit();
281  if (IsIon) {
285  }
286 
287 
289  /*if (preStepKinEnergy< 0.05*MeV) maxE =2.*preStepKinEnergy;
290  else if (preStepKinEnergy< 0.1*MeV) maxE =1.5*preStepKinEnergy;
291  else if (preStepKinEnergy< 0.5*MeV) maxE =1.25*preStepKinEnergy;*/
292 
293  if (preStepKinEnergy < currentTcut) maxE = std::min(currentTcut,maxE);
294 
295  maxE=std::min(emax_model*1.001,maxE);
296 
298 
299  if (IsIon) {
302  }
303 
305 
307 
308 
309 
310  x=r1-preStepRange;
311  x=std::max(r1-preStepRange,0.001*mm);
312 
313  return x;
314 
315 
316 }
void SetDynamicMassCharge(G4double massratio, G4double charge2ratio)
G4ParticleDefinition * theDirectPartDef
G4VEnergyLossProcess * theDirectEnergyLossProcess
const G4MaterialCutsCouple * currentCouple
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
double maxE
Definition: plot_hist.C:8
virtual G4double GetChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:345
G4VEmModel * SelectModelForMaterial(G4double kinEnergy, size_t &idx) const
void DefineMaterial(const G4MaterialCutsCouple *couple)
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
static const double mm
Definition: G4SIunits.hh:114
G4double GetRange(G4double &kineticEnergy, const G4MaterialCutsCouple *)
Here is the call graph for this function:

◆ operator=()

G4ContinuousGainOfEnergy& G4ContinuousGainOfEnergy::operator= ( const G4ContinuousGainOfEnergy right)
private
Here is the caller graph for this function:

◆ PreparePhysicsTable()

void G4ContinuousGainOfEnergy::PreparePhysicsTable ( const G4ParticleDefinition )
virtual

Reimplemented from G4VProcess.

Definition at line 82 of file G4ContinuousGainOfEnergy.cc.

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

◆ SetDirectEnergyLossProcess()

void G4ContinuousGainOfEnergy::SetDirectEnergyLossProcess ( G4VEnergyLossProcess aProcess)
inline

Definition at line 112 of file G4ContinuousGainOfEnergy.hh.

112 {theDirectEnergyLossProcess=aProcess;};
G4VEnergyLossProcess * theDirectEnergyLossProcess
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetDirectParticle()

void G4ContinuousGainOfEnergy::SetDirectParticle ( G4ParticleDefinition p)

Definition at line 99 of file G4ContinuousGainOfEnergy.cc.

101  if (theDirectPartDef->GetParticleType()== "nucleus") {
102  IsIon=true;
105  chargeSqRatio=q*q;
106 
107 
108  }
109 
110 }
const G4String & GetParticleType() const
G4ParticleDefinition * theDirectPartDef
float proton_mass_c2
Definition: hepunit.py:275
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetDynamicMassCharge()

void G4ContinuousGainOfEnergy::SetDynamicMassCharge ( const G4Track &  track,
G4double  energy 
)
private

Definition at line 321 of file G4ContinuousGainOfEnergy.cc.

322 {
323 
326 }
static G4LossTableManager * Instance()
void SetDynamicMassCharge(G4double massratio, G4double charge2ratio)
G4double EffectiveChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4ParticleDefinition * theDirectPartDef
G4VEnergyLossProcess * theDirectEnergyLossProcess
double energy
Definition: plottest35.C:25
G4EmCorrections * EmCorrections()
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetIsIntegral()

void G4ContinuousGainOfEnergy::SetIsIntegral ( G4bool  val)
inline

Definition at line 110 of file G4ContinuousGainOfEnergy.hh.

◆ SetLossFluctuations()

void G4ContinuousGainOfEnergy::SetLossFluctuations ( G4bool  val)

Definition at line 258 of file G4ContinuousGainOfEnergy.cc.

259 {
260  if(val && !lossFluctuationArePossible) return;
261  lossFluctuationFlag = val;
262 }
Here is the caller graph for this function:

Member Data Documentation

◆ chargeSqRatio

G4double G4ContinuousGainOfEnergy::chargeSqRatio
private

Definition at line 159 of file G4ContinuousGainOfEnergy.hh.

◆ currentCouple

const G4MaterialCutsCouple* G4ContinuousGainOfEnergy::currentCouple
private

Definition at line 136 of file G4ContinuousGainOfEnergy.hh.

◆ currentCoupleIndex

size_t G4ContinuousGainOfEnergy::currentCoupleIndex
private

Definition at line 138 of file G4ContinuousGainOfEnergy.hh.

◆ currentCutInRange

G4double G4ContinuousGainOfEnergy::currentCutInRange
private

Definition at line 140 of file G4ContinuousGainOfEnergy.hh.

◆ currentMaterial

const G4Material* G4ContinuousGainOfEnergy::currentMaterial
private

Definition at line 135 of file G4ContinuousGainOfEnergy.hh.

◆ currentMaterialIndex

size_t G4ContinuousGainOfEnergy::currentMaterialIndex
private

Definition at line 137 of file G4ContinuousGainOfEnergy.hh.

◆ currentModel

G4VEmModel* G4ContinuousGainOfEnergy::currentModel
private

Definition at line 160 of file G4ContinuousGainOfEnergy.hh.

◆ currentTcut

G4double G4ContinuousGainOfEnergy::currentTcut
private

Definition at line 139 of file G4ContinuousGainOfEnergy.hh.

◆ is_integral

G4bool G4ContinuousGainOfEnergy::is_integral
private

Definition at line 153 of file G4ContinuousGainOfEnergy.hh.

◆ IsIon

G4bool G4ContinuousGainOfEnergy::IsIon
private

Definition at line 157 of file G4ContinuousGainOfEnergy.hh.

◆ linLossLimit

G4double G4ContinuousGainOfEnergy::linLossLimit
private

Definition at line 145 of file G4ContinuousGainOfEnergy.hh.

◆ lossFluctuationArePossible

G4bool G4ContinuousGainOfEnergy::lossFluctuationArePossible
private

Definition at line 147 of file G4ContinuousGainOfEnergy.hh.

◆ lossFluctuationFlag

G4bool G4ContinuousGainOfEnergy::lossFluctuationFlag
private

Definition at line 146 of file G4ContinuousGainOfEnergy.hh.

◆ massRatio

G4double G4ContinuousGainOfEnergy::massRatio
private

Definition at line 158 of file G4ContinuousGainOfEnergy.hh.

◆ preStepChargeSqRatio

G4double G4ContinuousGainOfEnergy::preStepChargeSqRatio
private

Definition at line 161 of file G4ContinuousGainOfEnergy.hh.

◆ preStepKinEnergy

G4double G4ContinuousGainOfEnergy::preStepKinEnergy
private

Definition at line 141 of file G4ContinuousGainOfEnergy.hh.

◆ preStepRange

G4double G4ContinuousGainOfEnergy::preStepRange
private

Definition at line 163 of file G4ContinuousGainOfEnergy.hh.

◆ preStepScaledKinEnergy

G4double G4ContinuousGainOfEnergy::preStepScaledKinEnergy
private

Definition at line 162 of file G4ContinuousGainOfEnergy.hh.

◆ theDirectEnergyLossProcess

G4VEnergyLossProcess* G4ContinuousGainOfEnergy::theDirectEnergyLossProcess
private

Definition at line 149 of file G4ContinuousGainOfEnergy.hh.

◆ theDirectPartDef

G4ParticleDefinition* G4ContinuousGainOfEnergy::theDirectPartDef
private

Definition at line 150 of file G4ContinuousGainOfEnergy.hh.


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