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

#include <G4DynamicParticleIonisation.hh>

+ Inheritance diagram for G4DynamicParticleIonisation:

Public Member Functions

 G4DynamicParticleIonisation ()
 
 ~G4DynamicParticleIonisation () override
 
void BuildPhysicsTable (const G4ParticleDefinition &) override
 
G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection) override
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &) override
 
G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &) override
 
G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
G4double GetContinuousStepLimit (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety) override
 
void ProcessDescription (std::ostream &) const override
 
G4DynamicParticleIonisationoperator= (const G4DynamicParticleIonisation &right)=delete
 
 G4DynamicParticleIonisation (const G4DynamicParticleIonisation &)=delete
 
- Public Member Functions inherited from G4VContinuousDiscreteProcess
 G4VContinuousDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VContinuousDiscreteProcess (G4VContinuousDiscreteProcess &)
 
virtual ~G4VContinuousDiscreteProcess ()
 
G4VContinuousDiscreteProcessoperator= (const G4VContinuousDiscreteProcess &)=delete
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (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 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
 
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 G4VContinuousDiscreteProcess
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 65 of file G4DynamicParticleIonisation.hh.

Constructor & Destructor Documentation

◆ G4DynamicParticleIonisation() [1/2]

G4DynamicParticleIonisation::G4DynamicParticleIonisation ( )

Definition at line 67 of file G4DynamicParticleIonisation.cc.

68 : G4VContinuousDiscreteProcess("dynPartIoni")
69{
72 theElectron = G4Electron::Electron();
73
74 lManager = G4LossTableManager::Instance();
75 lManager->Register(this);
76
77 fUrban = new G4DynamicParticleFluctuation();
78
79 // define these flags only once
80 auto param = G4EmParameters::Instance();
81 fFluct = param->LossFluctuation();
82 fLinLimit = 5*param->LinearLossLimit();
83}
@ fDynamicIonisation
static G4Electron * Electron()
Definition G4Electron.cc:91
static G4EmParameters * Instance()
static G4LossTableManager * Instance()
G4VContinuousDiscreteProcess(const G4String &, G4ProcessType aType=fNotDefined)
void SetVerboseLevel(G4int value)
void SetProcessSubType(G4int)

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

◆ ~G4DynamicParticleIonisation()

G4DynamicParticleIonisation::~G4DynamicParticleIonisation ( )
override

Definition at line 87 of file G4DynamicParticleIonisation.cc.

88{
89 lManager->DeRegister(this);
90}

◆ G4DynamicParticleIonisation() [2/2]

G4DynamicParticleIonisation::G4DynamicParticleIonisation ( const G4DynamicParticleIonisation & )
delete

Member Function Documentation

◆ AlongStepDoIt()

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

Reimplemented from G4VContinuousDiscreteProcess.

Definition at line 190 of file G4DynamicParticleIonisation.cc.

192{
193 fParticleChange.InitializeForAlongStep(track);
194
195 // no energy loss
196 if (fCharge == 0.0) { return &fParticleChange; }
197
198 // stop low-energy object
199 if (fEkinPreStep <= fLowestEkin) {
200 fParticleChange.SetProposedKineticEnergy(0.0);
201 fParticleChange.ProposeLocalEnergyDeposit(fEkinPreStep);
202 return &fParticleChange;
203 }
204
205 G4double length = step.GetStepLength();
206 G4double dedxPre = ComputeDEDX(fEkinPreStep);
207 G4double eloss = dedxPre*length;
208 G4double ekinPostStep = fEkinPreStep - eloss;
209
210 // correction for large step if it is not the last step
211 if (fEkinPreStep*fLinLimit < eloss && ekinPostStep > fLowestEkin) {
212 G4double dedxPost = ComputeDEDX(ekinPostStep);
213 eloss = (eloss + dedxPost*length)*0.5;
214 }
215
216 // do not sample fluctuations at the last step
217 if (fFluct && fEkinPreStep > eloss) {
218 eloss = fUrban->SampleFluctuations(fCouple, track.GetDynamicParticle(),
219 fCut, fTmax, length, eloss);
220 }
221
222 ekinPostStep = fEkinPreStep - eloss;
223
224 // stop low-energy object
225 if (ekinPostStep <= fLowestEkin) {
226 fParticleChange.SetProposedKineticEnergy(0.0);
227 fParticleChange.ProposeLocalEnergyDeposit(fEkinPreStep);
228 } else {
229 fParticleChange.SetProposedKineticEnergy(ekinPostStep);
230 fParticleChange.ProposeLocalEnergyDeposit(eloss);
231 }
232 return &fParticleChange;
233}
double G4double
Definition G4Types.hh:83
G4double GetStepLength() const
const G4DynamicParticle * GetDynamicParticle() const

◆ AlongStepGetPhysicalInteractionLength()

G4double G4DynamicParticleIonisation::AlongStepGetPhysicalInteractionLength ( const G4Track & ,
G4double previousStepSize,
G4double currentMinimumStep,
G4double & currentSafety,
G4GPILSelection * selection )
overridevirtual

Reimplemented from G4VContinuousDiscreteProcess.

Definition at line 124 of file G4DynamicParticleIonisation.cc.

127{
128 *selection = CandidateForSelection;
129
130 // no step limit for the time being
131 return DBL_MAX;
132}
@ CandidateForSelection
#define DBL_MAX
Definition templates.hh:62

◆ BuildPhysicsTable()

void G4DynamicParticleIonisation::BuildPhysicsTable ( const G4ParticleDefinition & )
overridevirtual

Reimplemented from G4VProcess.

Definition at line 95 of file G4DynamicParticleIonisation.cc.

96{
97 auto theCoupleTable = G4ProductionCutsTable::GetProductionCutsTable();
98 fCuts = theCoupleTable->GetEnergyCutsVector(1);
99}
static G4ProductionCutsTable * GetProductionCutsTable()

◆ GetContinuousStepLimit()

G4double G4DynamicParticleIonisation::GetContinuousStepLimit ( const G4Track & track,
G4double previousStepSize,
G4double currentMinimumStep,
G4double & currentSafety )
overridevirtual

Implements G4VContinuousDiscreteProcess.

Definition at line 345 of file G4DynamicParticleIonisation.cc.

347{
348 return DBL_MAX;
349}

◆ GetMeanFreePath()

G4double G4DynamicParticleIonisation::GetMeanFreePath ( const G4Track & track,
G4double previousStepSize,
G4ForceCondition * condition )
overridevirtual

Implements G4VContinuousDiscreteProcess.

Definition at line 334 of file G4DynamicParticleIonisation.cc.

336{
337 // Note: this method is not used at run-time, so its implementation is simplified.
338 // It might be eventually refined later.
340 return DBL_MAX;
341}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced

◆ operator=()

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

◆ PostStepDoIt()

G4VParticleChange * G4DynamicParticleIonisation::PostStepDoIt ( const G4Track & track,
const G4Step &  )
overridevirtual

Reimplemented from G4VContinuousDiscreteProcess.

Definition at line 238 of file G4DynamicParticleIonisation.cc.

239{
241 fParticleChange.InitializeForPostStep(track);
242
243 auto dp = track.GetDynamicParticle();
244 G4double kinEnergy = dp->GetKineticEnergy();
245 const G4double totEnergy = kinEnergy + fMass;
246 const G4double beta2 = kinEnergy*(kinEnergy + 2.0*fMass)/(totEnergy*totEnergy);
247
248 G4double deltaKinEnergy, f;
249
250 CLHEP::HepRandomEngine* rndmEngineMod = G4Random::getTheEngine();
251 G4double rndm[2];
252
253 // sampling without nuclear size effect
254 do {
255 rndmEngineMod->flatArray(2, rndm);
256 deltaKinEnergy = fCut*fTmax/(fCut*(1.0 - rndm[0]) + fTmax*rndm[0]);
257 f = 1.0 - beta2*deltaKinEnergy/fTmax;
258 // Loop checking, 14-Aug-2024, Vladimir Ivanchenko
259 } while( rndm[1] > f);
260
261 G4double deltaMomentum =
262 std::sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*CLHEP::electron_mass_c2));
263 G4double cost = deltaKinEnergy * (totEnergy + CLHEP::electron_mass_c2) /
264 (deltaMomentum * dp->GetTotalMomentum());
265 cost = std::min(cost, 1.0);
266 const G4double sint = std::sqrt((1.0 - cost)*(1.0 + cost));
267 const G4double phi = CLHEP::twopi*rndmEngineMod->flat();
268
269 G4ThreeVector deltaDirection(sint*std::cos(phi), sint*std::sin(phi), cost);
270 deltaDirection.rotateUz(dp->GetMomentumDirection());
271
272 // create G4DynamicParticle object for delta ray
273 auto delta = new G4DynamicParticle(theElectron, deltaDirection, deltaKinEnergy);
274 auto t = new G4Track(delta, track.GetGlobalTime(), track.GetPosition());
275 t->SetTouchableHandle(track.GetTouchableHandle());
276 t->SetCreatorModelID(fSecID);
277 fParticleChange.AddSecondary(t);
278
279 // Change kinematics of primary particle
280 kinEnergy -= deltaKinEnergy;
281 G4ThreeVector finalP = dp->GetMomentum() - delta->GetMomentum();
282 finalP = finalP.unit();
283
284 fParticleChange.SetProposedKineticEnergy(kinEnergy);
285 fParticleChange.SetProposedMomentumDirection(finalP);
286 return &fParticleChange;
287}
CLHEP::Hep3Vector G4ThreeVector
Hep3Vector unit() const
virtual double flat()=0
virtual void flatArray(const int size, double *vect)=0
G4double GetKineticEnergy() const
const G4ThreeVector & GetPosition() const
G4double GetGlobalTime() const
const G4TouchableHandle & GetTouchableHandle() const
G4double theNumberOfInteractionLengthLeft

◆ PostStepGetPhysicalInteractionLength()

G4double G4DynamicParticleIonisation::PostStepGetPhysicalInteractionLength ( const G4Track & track,
G4double previousStepSize,
G4ForceCondition * condition )
overridevirtual

Reimplemented from G4VContinuousDiscreteProcess.

Definition at line 136 of file G4DynamicParticleIonisation.cc.

139{
141 G4double x = DBL_MAX;
142 G4double xsec = 0.0;
143
144 PreStepInitialisation(track);
145
146 if (fCharge != 0.0) {
147 xsec = ComputeCrossSection(fEkinPreStep);
148 }
149
150 if (xsec <= 0.0) {
153
154 } else {
156
159
160 } else if(currentInteractionLength < DBL_MAX) {
161
162 // subtract NumberOfInteractionLengthLeft using previous step
164 previousStepSize/currentInteractionLength;
165
168 }
169 currentInteractionLength = 1.0/xsec;
171 }
172#ifdef G4VERBOSE
173 if (verboseLevel>2) {
174 G4cout << "G4DynamicParticleIonisation::PostStepGetPhysicalInteractionLength ";
175 G4cout << " Process: " << GetProcessName()
176 << " for unknown particle Mass(GeV)=" << fMass/CLHEP::GeV
177 << " charge=" << fCharge
178 << " Material " << fMaterial->GetName()
179 << " Ekin(MeV)=" << fEkinPreStep/CLHEP::MeV
180 << " MFP(cm)=" << currentInteractionLength/CLHEP::cm
181 << " ProposedLength(cm)=" << x/CLHEP::cm <<G4endl;
182 }
183#endif
184 return x;
185}
G4double G4Log(G4double x)
Definition G4Log.hh:227
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52
G4double currentInteractionLength
G4double theInitialNumberOfInteractionLength
G4int verboseLevel
const G4String & GetProcessName() const

◆ ProcessDescription()

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

Reimplemented from G4VProcess.

Definition at line 353 of file G4DynamicParticleIonisation.cc.

354{
355 out << "G4DynamicParticleIonisation: dynamic ionisation" << G4endl;
356}

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