Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Member Functions
G4ContinuousGainOfEnergy Class Reference

#include <G4ContinuousGainOfEnergy.hh>

Inheritance diagram for G4ContinuousGainOfEnergy:
G4VContinuousProcess G4VProcess

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 &, 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 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 ()
 
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 ()
 

Additional Inherited Members

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

Constructor & Destructor Documentation

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;
51  lossFluctuationArePossible =true;
52  lossFluctuationFlag=true;
53  is_integral = false;
54 
55  //Will be properly set in SetDirectParticle()
56  IsIon=false;
57  massRatio =1.;
58  chargeSqRatio=1.;
59  preStepChargeSqRatio=1.;
60 
61  //Some initialization
62  currentCoupleIndex=9999999;
63  currentCutInRange=0.;
64  currentMaterialIndex=9999999;
65  currentTcut=0.;
66  preStepKinEnergy=0.;
67  preStepRange=0.;
68  preStepScaledKinEnergy=0.;
69 
70  currentCouple=0;
71 }
G4ContinuousGainOfEnergy::~G4ContinuousGainOfEnergy ( )
virtual

Definition at line 75 of file G4ContinuousGainOfEnergy.cc.

76 {
77 
78 }

Member Function Documentation

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

Reimplemented from G4VContinuousProcess.

Definition at line 115 of file G4ContinuousGainOfEnergy.cc.

References G4VProcess::aParticleChange, G4VEmModel::CorrectionsAlongStep(), G4VEmModel::GetChargeSquareRatio(), G4VEnergyLossProcess::GetDEDX(), G4Track::GetDynamicParticle(), G4DynamicParticle::GetKineticEnergy(), G4VEnergyLossProcess::GetKineticEnergy(), G4VEmModel::GetModelOfFluctuations(), G4Step::GetPostStepPoint(), G4VEnergyLossProcess::GetRange(), G4Step::GetStepLength(), G4StepPoint::GetWeight(), G4ParticleChange::Initialize(), G4VEmModel::MaxSecondaryKinEnergy(), G4INCL::Math::min(), n, G4ParticleChange::ProposeEnergy(), G4VParticleChange::ProposeParentWeight(), G4DynamicParticle::SetDefinition(), G4VEnergyLossProcess::SetDynamicMassCharge(), G4DynamicParticle::SetKineticEnergy(), G4VParticleChange::SetParentWeightByProcess(), and test::x.

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 
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  //-------------------------------------------------------------
137  G4double DEDX_before = theDirectEnergyLossProcess->GetDEDX(preStepKinEnergy, currentCouple);
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) {
156  chargeSqRatio = currentModel->GetChargeSquareRatio(theDirectPartDef,currentMaterial,Tkin);
157  theDirectEnergyLossProcess->SetDynamicMassCharge(massRatio,chargeSqRatio);
158 
159  }
160 
161  G4double r = theDirectEnergyLossProcess->GetRange(Tkin, currentCouple);
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);
168  G4double E = theDirectEnergyLossProcess->GetKineticEnergy(x,currentCouple);
169  if (IsIon){
170  chargeSqRatio = currentModel->GetChargeSquareRatio(theDirectPartDef,currentMaterial,E);
171  theDirectEnergyLossProcess->SetDynamicMassCharge(massRatio,chargeSqRatio);
172  G4double x1= theDirectEnergyLossProcess->GetRange(E, currentCouple);
173  while (std::abs(x-x1)>0.01*x) {
174  E = theDirectEnergyLossProcess->GetKineticEnergy(x,currentCouple);
175  chargeSqRatio = currentModel->GetChargeSquareRatio(theDirectPartDef,currentMaterial,E);
176  theDirectEnergyLossProcess->SetDynamicMassCharge(massRatio,chargeSqRatio);
177  x1= theDirectEnergyLossProcess->GetRange(E, currentCouple);
178 
179  }
180  }
181 
182  degain=E-Tkin;
183 
184 
185 
186  }
187  //G4cout<<degain<<G4endl;
188  G4double tmax = currentModel->MaxSecondaryKinEnergy(dynParticle);
189  tmax = std::min(tmax,currentTcut);
190 
191 
192  dynParticle->SetKineticEnergy(Tkin+degain);
193 
194  // Corrections, which cannot be tabulated for ions
195  //----------------------------------------
196  G4double esecdep=0;//not used in most models
197  currentModel->CorrectionsAlongStep(currentCouple, dynParticle, degain,esecdep, dlength);
198 
199  // Sample fluctuations
200  //-------------------
201 
202 
203  G4double deltaE =0.;
204  if (lossFluctuationFlag ) {
205  deltaE = currentModel->GetModelOfFluctuations()->
206  SampleFluctuations(currentCouple,dynParticle,tmax,dlength,degain)-degain;
207  }
208 
209  G4double egain=degain+deltaE;
210  if (egain <=0) egain=degain;
211  Tkin+=egain;
212  dynParticle->SetKineticEnergy(Tkin);
213  }
214 
215 
216 
217 
218 
219  delete dynParticle;
220 
221  if (IsIon){
222  chargeSqRatio = currentModel->GetChargeSquareRatio(theDirectPartDef,currentMaterial,Tkin);
223  theDirectEnergyLossProcess->SetDynamicMassCharge(massRatio,chargeSqRatio);
224 
225  }
226 
227  G4double DEDX_after = theDirectEnergyLossProcess->GetDEDX(Tkin, currentCouple);
228 
229 
230  G4double weight_correction=DEDX_after/DEDX_before;
231 
232 
234 
235 
236  //Caution!!!
237  // It is important to select the weight of the post_step_point
238  // as the current weight and not the weight of the track, as t
239  // the weight of the track is changed after having applied all
240  // the along_step_do_it.
241 
242  // G4double new_weight=weight_correction*track.GetWeight(); //old
243  G4double new_weight=weight_correction*step.GetPostStepPoint()->GetWeight();
246 
247 
248  return &aParticleChange;
249 
250 }
G4double GetKineticEnergy(G4double &range, const G4MaterialCutsCouple *)
G4double MaxSecondaryKinEnergy(const G4DynamicParticle *dynParticle)
Definition: G4VEmModel.hh:448
G4double GetWeight() const
virtual void CorrectionsAlongStep(const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
Definition: G4VEmModel.cc:339
void SetDynamicMassCharge(G4double massratio, G4double charge2ratio)
G4double GetKineticEnergy() const
G4double GetStepLength() const
const G4DynamicParticle * GetDynamicParticle() const
G4double GetDEDX(G4double &kineticEnergy, const G4MaterialCutsCouple *)
void ProposeParentWeight(G4double finalWeight)
G4VEmFluctuationModel * GetModelOfFluctuations()
Definition: G4VEmModel.hh:571
void SetParentWeightByProcess(G4bool)
const G4int n
void SetKineticEnergy(G4double aEnergy)
virtual G4double GetChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:322
virtual void Initialize(const G4Track &)
G4StepPoint * GetPostStepPoint() const
void ProposeEnergy(G4double finalEnergy)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
double G4double
Definition: G4Types.hh:76
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
G4double GetRange(G4double &kineticEnergy, const G4MaterialCutsCouple *)
void G4ContinuousGainOfEnergy::BuildPhysicsTable ( const G4ParticleDefinition )
virtual

Reimplemented from G4VProcess.

Definition at line 92 of file G4ContinuousGainOfEnergy.cc.

93 {//theDirectEnergyLossProcess->BuildPhysicsTable(part);
94 ;
95 }
G4double G4ContinuousGainOfEnergy::GetContinuousStepLimit ( const G4Track track,
G4double  previousStepSize,
G4double  currentMinimumStep,
G4double currentSafety 
)
protectedvirtual

Implements G4VContinuousProcess.

Definition at line 263 of file G4ContinuousGainOfEnergy.cc.

References DBL_MAX, G4VEmModel::GetChargeSquareRatio(), G4Track::GetKineticEnergy(), G4Track::GetMaterialCutsCouple(), G4VEnergyLossProcess::GetRange(), G4VEmModel::HighEnergyLimit(), G4INCL::Math::max(), G4INCL::Math::min(), python.hepunit::mm, G4VEnergyLossProcess::SelectModelForMaterial(), G4VEnergyLossProcess::SetDynamicMassCharge(), and test::x.

265 {
266  G4double x = DBL_MAX;
267  x=.1*mm;
268 
269 
270  DefineMaterial(track.GetMaterialCutsCouple());
271 
272  preStepKinEnergy = track.GetKineticEnergy();
273  preStepScaledKinEnergy = track.GetKineticEnergy()*massRatio;
274  currentModel = theDirectEnergyLossProcess->SelectModelForMaterial(preStepScaledKinEnergy,currentCoupleIndex);
275  G4double emax_model=currentModel->HighEnergyLimit();
276  if (IsIon) {
277  chargeSqRatio = currentModel->GetChargeSquareRatio(theDirectPartDef,currentMaterial,preStepKinEnergy);
278  preStepChargeSqRatio = chargeSqRatio;
279  theDirectEnergyLossProcess->SetDynamicMassCharge(massRatio,preStepChargeSqRatio);
280  }
281 
282 
283  G4double maxE =1.1*preStepKinEnergy;
284  /*if (preStepKinEnergy< 0.05*MeV) maxE =2.*preStepKinEnergy;
285  else if (preStepKinEnergy< 0.1*MeV) maxE =1.5*preStepKinEnergy;
286  else if (preStepKinEnergy< 0.5*MeV) maxE =1.25*preStepKinEnergy;*/
287 
288  if (preStepKinEnergy < currentTcut) maxE = std::min(currentTcut,maxE);
289 
290  maxE=std::min(emax_model*1.001,maxE);
291 
292  preStepRange = theDirectEnergyLossProcess->GetRange(preStepKinEnergy, currentCouple);
293 
294  if (IsIon) {
295  G4double chargeSqRatioAtEmax = currentModel->GetChargeSquareRatio(theDirectPartDef,currentMaterial,maxE);
296  theDirectEnergyLossProcess->SetDynamicMassCharge(massRatio,chargeSqRatioAtEmax);
297  }
298 
299  G4double r1 = theDirectEnergyLossProcess->GetRange(maxE, currentCouple);
300 
301  if (IsIon) theDirectEnergyLossProcess->SetDynamicMassCharge(massRatio,preStepChargeSqRatio);
302 
303 
304 
305  x=r1-preStepRange;
306  x=std::max(r1-preStepRange,0.001*mm);
307 
308  return x;
309 
310 
311 }
void SetDynamicMassCharge(G4double massratio, G4double charge2ratio)
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:592
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
G4double GetKineticEnergy() const
virtual G4double GetChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:322
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
G4VEmModel * SelectModelForMaterial(G4double kinEnergy, size_t &idx) const
G4double GetRange(G4double &kineticEnergy, const G4MaterialCutsCouple *)
void G4ContinuousGainOfEnergy::PreparePhysicsTable ( const G4ParticleDefinition )
virtual

Reimplemented from G4VProcess.

Definition at line 82 of file G4ContinuousGainOfEnergy.cc.

84 {//theDirectEnergyLossProcess->PreparePhysicsTable(part);
85 
86 ;
87 }
void G4ContinuousGainOfEnergy::SetDirectEnergyLossProcess ( G4VEnergyLossProcess aProcess)
inline

Definition at line 112 of file G4ContinuousGainOfEnergy.hh.

Referenced by G4AdjointPhysicsList::ConstructEM().

112 {theDirectEnergyLossProcess=aProcess;};
void G4ContinuousGainOfEnergy::SetDirectParticle ( G4ParticleDefinition p)

Definition at line 99 of file G4ContinuousGainOfEnergy.cc.

References G4ParticleDefinition::GetParticleType(), G4ParticleDefinition::GetPDGCharge(), G4ParticleDefinition::GetPDGMass(), and python.hepunit::proton_mass_c2.

Referenced by G4AdjointPhysicsList::ConstructEM().

100 {theDirectPartDef=p;
101  if (theDirectPartDef->GetParticleType()== "nucleus") {
102  IsIon=true;
103  massRatio = proton_mass_c2/theDirectPartDef->GetPDGMass();
104  G4double q=theDirectPartDef->GetPDGCharge();
105  chargeSqRatio=q*q;
106 
107 
108  }
109 
110 }
const char * p
Definition: xmltok.h:285
float proton_mass_c2
Definition: hepunit.py:275
const G4String & GetParticleType() const
G4double GetPDGMass() const
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
void G4ContinuousGainOfEnergy::SetIsIntegral ( G4bool  val)
inline

Definition at line 110 of file G4ContinuousGainOfEnergy.hh.

110 {is_integral= val;}
void G4ContinuousGainOfEnergy::SetLossFluctuations ( G4bool  val)

Definition at line 253 of file G4ContinuousGainOfEnergy.cc.

Referenced by G4AdjointPhysicsList::ConstructEM().

254 {
255  if(val && !lossFluctuationArePossible) return;
256  lossFluctuationFlag = val;
257 }

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