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

#include <G4VParticleChange.hh>

+ Inheritance diagram for G4VParticleChange:

Public Member Functions

 G4VParticleChange ()
 
virtual ~G4VParticleChange ()=default
 
 G4VParticleChange (const G4VParticleChange &right)=delete
 
G4VParticleChangeoperator= (const G4VParticleChange &right)=delete
 
virtual G4StepUpdateStepForAtRest (G4Step *Step)
 
virtual G4StepUpdateStepForAlongStep (G4Step *Step)
 
virtual G4StepUpdateStepForPostStep (G4Step *Step)
 
virtual void Initialize (const G4Track &)
 
G4double GetTrueStepLength () const
 
void ProposeTrueStepLength (G4double truePathLength)
 
G4double GetLocalEnergyDeposit () const
 
void ProposeLocalEnergyDeposit (G4double anEnergyPart)
 
G4double GetNonIonizingEnergyDeposit () const
 
void ProposeNonIonizingEnergyDeposit (G4double anEnergyPart)
 
G4TrackStatus GetTrackStatus () const
 
void ProposeTrackStatus (G4TrackStatus status)
 
const G4TrackGetCurrentTrack () const
 
G4SteppingControl GetSteppingControl () const
 
void ProposeSteppingControl (G4SteppingControl StepControlFlag)
 
G4bool GetFirstStepInVolume () const
 
G4bool GetLastStepInVolume () const
 
void ProposeFirstStepInVolume (G4bool flag)
 
void ProposeLastStepInVolume (G4bool flag)
 
void Clear ()
 
void SetNumberOfSecondaries (G4int totSecondaries)
 
G4int GetNumberOfSecondaries () const
 
G4TrackGetSecondary (G4int anIndex) const
 
void AddSecondary (G4Track *aSecondary)
 
G4double GetWeight () const
 
G4double GetParentWeight () const
 
void ProposeWeight (G4double finalWeight)
 
void ProposeParentWeight (G4double finalWeight)
 
void SetSecondaryWeightByProcess (G4bool)
 
G4bool IsSecondaryWeightSetByProcess () const
 
void SetParentWeightByProcess (G4bool)
 
G4bool IsParentWeightSetByProcess () const
 
virtual void DumpInfo () const
 
void SetVerboseLevel (G4int vLevel)
 
G4int GetVerboseLevel () const
 
virtual G4bool CheckIt (const G4Track &)
 
void ClearDebugFlag ()
 
void SetDebugFlag ()
 
G4bool GetDebugFlag () const
 

Protected Member Functions

G4StepUpdateStepInfo (G4Step *Step)
 
void InitializeLocalEnergyDeposit ()
 
void InitializeSteppingControl ()
 
void InitializeParentWeight (const G4Track &)
 
void InitializeStatusChange (const G4Track &)
 
void InitializeSecondaries ()
 
void InitializeFromStep (const G4Step *)
 
G4double ComputeBeta (G4double kinEnergy)
 
G4bool CheckSecondary (G4Track &)
 
G4double GetAccuracyForWarning () const
 
G4double GetAccuracyForException () const
 

Protected Attributes

const G4TracktheCurrentTrack = nullptr
 
std::vector< G4Track * > theListOfSecondaries
 
G4TrackStatus theStatusChange = fAlive
 
G4SteppingControl theSteppingControlFlag = NormalCondition
 
G4double theLocalEnergyDeposit = 0.0
 
G4double theNonIonizingEnergyDeposit = 0.0
 
G4double theTrueStepLength = 0.0
 
G4double theParentWeight = 1.0
 
G4double theParentGlobalTime = 0.0
 
G4int theNumberOfSecondaries = 0
 
G4int theSizeOftheListOfSecondaries = 0
 
G4int verboseLevel = 1
 
G4int nError = 0
 
G4bool theFirstStepInVolume = false
 
G4bool theLastStepInVolume = false
 
G4bool isParentWeightProposed = false
 
G4bool fSetSecondaryWeightByProcess = false
 
G4bool debugFlag = false
 

Static Protected Attributes

static const G4double accuracyForWarning = 1.0e-9
 
static const G4double accuracyForException = 0.001
 
static const G4int maxError = 10
 

Detailed Description

Definition at line 68 of file G4VParticleChange.hh.

Constructor & Destructor Documentation

◆ G4VParticleChange() [1/2]

G4VParticleChange::G4VParticleChange ( )

Definition at line 40 of file G4VParticleChange.cc.

41{
42#ifdef G4VERBOSE
43 // activate CheckIt if in VERBOSE mode
44 debugFlag = true;
45#endif
46}

◆ ~G4VParticleChange()

virtual G4VParticleChange::~G4VParticleChange ( )
virtualdefault

◆ G4VParticleChange() [2/2]

G4VParticleChange::G4VParticleChange ( const G4VParticleChange right)
delete

Member Function Documentation

◆ AddSecondary()

void G4VParticleChange::AddSecondary ( G4Track aSecondary)

Definition at line 49 of file G4VParticleChange.cc.

50{
51 if(debugFlag)
52 CheckSecondary(*aTrack);
53
55 aTrack->SetWeight(theParentWeight);
56
57 // add a secondary after size check
59 {
61 }
62 else
63 {
64 theListOfSecondaries.push_back(aTrack);
66 }
68}
std::vector< G4Track * > theListOfSecondaries
G4bool CheckSecondary(G4Track &)

Referenced by G4ParticleChangeForGamma::AddSecondary(), G4ParticleChange::AddSecondary(), G4eplusAnnihilation::AtRestDoIt(), G4FastStep::CreateSecondaryTrack(), G4RadioactiveDecay::DecayAnalog(), G4Decay::DecayIt(), G4UnknownDecay::DecayIt(), G4Radioactivation::DecayIt(), G4VEmProcess::PostStepDoIt(), G4VEnergyLossProcess::PostStepDoIt(), and G4ParticleChangeForOccurenceBiasing::StealSecondaries().

◆ CheckIt()

G4bool G4VParticleChange::CheckIt ( const G4Track aTrack)
virtual

Reimplemented in G4FastStep, and G4ParticleChangeForDecay.

Definition at line 221 of file G4VParticleChange.cc.

222{
223 G4bool isOK = true;
224
225 // Energy deposit should not be negative
226 if(theLocalEnergyDeposit < 0.0)
227 {
228 isOK = false;
229 ++nError;
230#ifdef G4VERBOSE
231 if(nError < maxError)
232 {
233 G4cout << " G4VParticleChange::CheckIt : ";
234 G4cout << "the energy deposit " << theLocalEnergyDeposit/MeV
235 << " MeV is negative !!" << G4endl;
236 }
237#endif
239 }
240
241 // true path length should not be negative
242 if(theTrueStepLength < 0.0)
243 {
244 isOK = false;
245 ++nError;
246#ifdef G4VERBOSE
247 if(nError < maxError)
248 {
249 G4cout << " G4VParticleChange::CheckIt : ";
250 G4cout << "true path length " << theTrueStepLength/mm
251 << " mm is negative !!" << G4endl;
252 }
253#endif
254 theTrueStepLength = (1.e-12) * mm;
255 }
256
257 if(!isOK)
258 {
259 if(nError < maxError)
260 {
261#ifdef G4VERBOSE
262 // dump out information of this particle change
263 DumpInfo();
264#endif
265 G4Exception("G4VParticleChange::CheckIt()", "TRACK001", JustWarning,
266 "Step length and/or energy deposit are illegal");
267 }
268 }
269 return isOK;
270}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
bool G4bool
Definition: G4Types.hh:86
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static const G4int maxError
virtual void DumpInfo() const

Referenced by G4FastStep::CheckIt(), G4ParticleChangeForDecay::CheckIt(), G4ParticleChangeForTransport::UpdateStepForAlongStep(), G4ParticleChange::UpdateStepForAlongStep(), G4ParticleChange::UpdateStepForAtRest(), and G4ParticleChange::UpdateStepForPostStep().

◆ CheckSecondary()

G4bool G4VParticleChange::CheckSecondary ( G4Track aTrack)
protected

Definition at line 273 of file G4VParticleChange.cc.

274{
275 G4bool isOK = true;
276
277 // MomentumDirection should be unit vector
278 G4double ekin = aTrack.GetKineticEnergy();
279 auto dir = aTrack.GetMomentumDirection();
280 G4double accuracy = std::abs(dir.mag2() - 1.0);
281 if(accuracy > accuracyForWarning)
282 {
283 isOK = false;
284 ++nError;
285#ifdef G4VERBOSE
286 if(nError < maxError)
287 {
288 G4String mname = aTrack.GetCreatorModelName();
289 G4cout << " G4VParticleChange::CheckSecondary : " << G4endl;
290 G4cout << " the momentum direction " << dir
291 << " is not unit vector !!" << G4endl;
292 G4cout << " Difference=" << accuracy
293 << " Ekin(MeV)=" << ekin/MeV
294 << " " << aTrack.GetParticleDefinition()->GetParticleName()
295 << " created by " << mname
296 << G4endl;
297 }
298#endif
299 aTrack.SetMomentumDirection(dir.unit());
300 }
301
302 // Kinetic Energy should not be negative
303 if(ekin < 0.0)
304 {
305 isOK = false;
306 ++nError;
307#ifdef G4VERBOSE
308 if(nError < maxError)
309 {
310 G4String mname = aTrack.GetCreatorModelName();
311 G4cout << " G4VParticleChange::CheckSecondary : " << G4endl;
312 G4cout << " Ekin(MeV)=" << ekin << " is negative !! "
314 << " created by " << mname
315 << G4endl;
316 }
317#endif
318 aTrack.SetKineticEnergy(0.0);
319 }
320
321 // Check timing of secondaries
322 G4double time = aTrack.GetGlobalTime();
323 if(time < theParentGlobalTime)
324 {
325 isOK = false;
326 ++nError;
327#ifdef G4VERBOSE
328 if(nError < maxError)
329 {
330 G4String mname = aTrack.GetCreatorModelName();
331 G4cout << " G4VParticleChange::CheckSecondary : " << G4endl;
332 G4cout << " The global time of secondary goes back compared to the parent !!" << G4endl;
333 G4cout << " ParentTime(ns)=" << theParentGlobalTime/ns
334 << " SecondaryTime(ns)= " << time/ns
335 << " Difference(ns)=" << (theParentGlobalTime - time)/ns
336 << G4endl;
337 G4cout << " Ekin(MeV)=" << ekin
339 << " created by " << mname << G4endl;
340 }
341#endif
343 }
344
345 // Exit with error
346 if(!isOK)
347 {
348 if(nError < maxError)
349 {
350#ifdef G4VERBOSE
351 DumpInfo();
352#endif
353 G4Exception("G4VParticleChange::CheckSecondary()", "TRACK001",
354 JustWarning, "Secondary with illegal time and/or energy and/or momentum");
355 }
356 }
357 return isOK;
358}
double G4double
Definition: G4Types.hh:83
const G4String & GetParticleName() const
const G4ParticleDefinition * GetParticleDefinition() const
const G4String GetCreatorModelName() const
G4double GetGlobalTime() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
void SetKineticEnergy(const G4double aValue)
void SetMomentumDirection(const G4ThreeVector &aValue)
void SetGlobalTime(const G4double aValue)
static const G4double accuracyForWarning
#define ns(x)
Definition: xmltok.c:1649

Referenced by AddSecondary().

◆ Clear()

◆ ClearDebugFlag()

void G4VParticleChange::ClearDebugFlag ( )
inline

◆ ComputeBeta()

G4double G4VParticleChange::ComputeBeta ( G4double  kinEnergy)
inlineprotected

◆ DumpInfo()

void G4VParticleChange::DumpInfo ( ) const
virtual

Reimplemented in G4FastStep, G4ParticleChangeForDecay, G4ParticleChangeForLoss, G4ParticleChangeForTransport, G4ParticleChange, and G4ParticleChangeForGamma.

Definition at line 133 of file G4VParticleChange.cc.

134{
135 auto vol = theCurrentTrack->GetVolume();
136 G4String vname = (nullptr == vol) ? "" : vol->GetName();
137 G4long olprc = G4cout.precision(8);
138 G4cout << " -----------------------------------------------" << G4endl;
139 G4cout << " G4VParticleChange Information " << G4endl;
140 G4cout << " TrackID : " << theCurrentTrack->GetTrackID()
141 << G4endl;
142 G4cout << " ParentID : " << theCurrentTrack->GetParentID()
143 << G4endl;
144 G4cout << " Particle : "
146 << G4endl;
147 G4cout << " Kinetic energy (MeV): "
149 G4cout << " Position (mm) : "
151 G4cout << " Direction : "
153 G4cout << " PhysicsVolume : " << vname << G4endl;
154 G4cout << " Material : "
156 G4cout << " -----------------------------------------------" << G4endl;
157
158 G4cout << " # of secondaries : " << std::setw(20)
160
161 G4cout << " -----------------------------------------------" << G4endl;
162
163 G4cout << " Energy Deposit (MeV): " << std::setw(20)
164 << theLocalEnergyDeposit / MeV << G4endl;
165
166 G4cout << " NIEL Energy Deposit (MeV): " << std::setw(20)
168
169 G4cout << " Track Status : " << std::setw(20);
171 {
172 G4cout << " Alive";
173 }
175 {
176 G4cout << " StopButAlive";
177 }
178 else if(theStatusChange == fStopAndKill)
179 {
180 G4cout << " StopAndKill";
181 }
183 {
184 G4cout << " KillTrackAndSecondaries";
185 }
186 else if(theStatusChange == fSuspend)
187 {
188 G4cout << " Suspend";
189 }
191 {
192 G4cout << " PostponeToNextEvent";
193 }
194 G4cout << G4endl;
195 G4cout << " TruePathLength (mm) : " << std::setw(20)
196 << theTrueStepLength / mm << G4endl;
197 G4cout << " Stepping Control : " << std::setw(20)
200 {
201 G4cout << " First step in volume" << G4endl;
202 }
204 {
205 G4cout << " Last step in volume" << G4endl;
206 }
207
208#ifdef G4VERBOSE
209 if(nError == maxError)
210 {
211 G4cout << " -----------------------------------------------" << G4endl;
212 G4cout << " G4VParticleChange warnings closed " << G4endl;
213 G4cout << " -----------------------------------------------" << G4endl;
214 }
215#endif
216
217 G4cout.precision(olprc);
218}
@ fKillTrackAndSecondaries
@ fSuspend
@ fAlive
@ fStopAndKill
@ fStopButAlive
@ fPostponeToNextEvent
long G4long
Definition: G4Types.hh:87
const G4String & GetName() const
Definition: G4Material.hh:172
G4int GetTrackID() const
G4VPhysicalVolume * GetVolume() const
const G4ThreeVector & GetPosition() const
G4Material * GetMaterial() const
G4int GetParentID() const
G4TrackStatus theStatusChange
G4SteppingControl theSteppingControlFlag
G4double theNonIonizingEnergyDeposit
const G4Track * theCurrentTrack

Referenced by CheckIt(), CheckSecondary(), G4FastStep::DumpInfo(), G4ParticleChangeForDecay::DumpInfo(), G4ParticleChangeForLoss::DumpInfo(), G4ParticleChange::DumpInfo(), G4ParticleChangeForGamma::DumpInfo(), G4ITSteppingVerbose::VerboseParticleChange(), G4SteppingVerbose::VerboseParticleChange(), and G4SteppingVerboseWithUnits::VerboseParticleChange().

◆ GetAccuracyForException()

G4double G4VParticleChange::GetAccuracyForException ( ) const
protected

Definition at line 367 of file G4VParticleChange.cc.

368{
370}
static const G4double accuracyForException

Referenced by G4FastStep::CheckIt().

◆ GetAccuracyForWarning()

G4double G4VParticleChange::GetAccuracyForWarning ( ) const
protected

Definition at line 361 of file G4VParticleChange.cc.

362{
363 return accuracyForWarning;
364}

Referenced by G4FastStep::CheckIt().

◆ GetCurrentTrack()

◆ GetDebugFlag()

G4bool G4VParticleChange::GetDebugFlag ( ) const
inline

◆ GetFirstStepInVolume()

G4bool G4VParticleChange::GetFirstStepInVolume ( ) const
inline

◆ GetLastStepInVolume()

G4bool G4VParticleChange::GetLastStepInVolume ( ) const
inline

◆ GetLocalEnergyDeposit()

G4double G4VParticleChange::GetLocalEnergyDeposit ( ) const
inline

◆ GetNonIonizingEnergyDeposit()

G4double G4VParticleChange::GetNonIonizingEnergyDeposit ( ) const
inline

◆ GetNumberOfSecondaries()

◆ GetParentWeight()

◆ GetSecondary()

◆ GetSteppingControl()

G4SteppingControl G4VParticleChange::GetSteppingControl ( ) const
inline

◆ GetTrackStatus()

◆ GetTrueStepLength()

G4double G4VParticleChange::GetTrueStepLength ( ) const
inline

◆ GetVerboseLevel()

G4int G4VParticleChange::GetVerboseLevel ( ) const
inline

◆ GetWeight()

G4double G4VParticleChange::GetWeight ( ) const
inline

◆ Initialize()

◆ InitializeFromStep()

void G4VParticleChange::InitializeFromStep ( const G4Step )
inlineprotected

◆ InitializeLocalEnergyDeposit()

void G4VParticleChange::InitializeLocalEnergyDeposit ( )
inlineprotected

◆ InitializeParentWeight()

void G4VParticleChange::InitializeParentWeight ( const G4Track )
inlineprotected

◆ InitializeSecondaries()

void G4VParticleChange::InitializeSecondaries ( )
inlineprotected

◆ InitializeStatusChange()

void G4VParticleChange::InitializeStatusChange ( const G4Track )
inlineprotected

◆ InitializeSteppingControl()

void G4VParticleChange::InitializeSteppingControl ( )
inlineprotected

◆ IsParentWeightSetByProcess()

G4bool G4VParticleChange::IsParentWeightSetByProcess ( ) const

Definition at line 376 of file G4VParticleChange.cc.

376{ return true; }

◆ IsSecondaryWeightSetByProcess()

G4bool G4VParticleChange::IsSecondaryWeightSetByProcess ( ) const
inline

◆ operator=()

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

◆ ProposeFirstStepInVolume()

void G4VParticleChange::ProposeFirstStepInVolume ( G4bool  flag)
inline

◆ ProposeLastStepInVolume()

void G4VParticleChange::ProposeLastStepInVolume ( G4bool  flag)
inline

◆ ProposeLocalEnergyDeposit()

void G4VParticleChange::ProposeLocalEnergyDeposit ( G4double  anEnergyPart)
inline

Referenced by G4VEnergyLossProcess::AlongStepDoIt(), G4ErrorEnergyLoss::AlongStepDoIt(), G4NuclearStopping::AlongStepDoIt(), G4hImpactIonisation::AlongStepDoIt(), G4HadronStoppingProcess::AtRestDoIt(), G4MuonMinusAtomicCapture::AtRestDoIt(), G4eplusAnnihilation::AtRestDoIt(), G4RadioactiveDecay::DecayAnalog(), G4DNAMolecularDissociation::DecayIt(), G4Decay::DecayIt(), G4UnknownDecay::DecayIt(), G4Radioactivation::DecayIt(), G4RadioactiveDecay::DecayIt(), G4HadronicProcess::FillResult(), G4SpecialCuts::PostStepDoIt(), G4UserSpecialCuts::PostStepDoIt(), G4LowECapture::PostStepDoIt(), G4VEmProcess::PostStepDoIt(), G4VEnergyLossProcess::PostStepDoIt(), G4NeutrinoElectronProcess::PostStepDoIt(), G4UCNBoundaryProcess::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4TauNeutrinoNucleusProcess::PostStepDoIt(), G4OpAbsorption::PostStepDoIt(), G4OpBoundaryProcess::PostStepDoIt(), G4hImpactIonisation::PostStepDoIt(), G4SynchrotronRadiationInMat::PostStepDoIt(), G4DNABornExcitationModel1::SampleSecondaries(), G4DNABornExcitationModel2::SampleSecondaries(), G4DNABornIonisationModel1::SampleSecondaries(), G4DNABornIonisationModel2::SampleSecondaries(), G4DNACPA100ElasticModel::SampleSecondaries(), G4DNACPA100ExcitationModel::SampleSecondaries(), G4DNACPA100IonisationModel::SampleSecondaries(), G4DNADingfelderChargeIncreaseModel::SampleSecondaries(), G4DNADiracRMatrixExcitationModel::SampleSecondaries(), G4DNAELSEPAElasticModel::SampleSecondaries(), G4DNAEmfietzoglouExcitationModel::SampleSecondaries(), G4DNAEmfietzoglouIonisationModel::SampleSecondaries(), G4DNAIonElasticModel::SampleSecondaries(), G4DNAMeltonAttachmentModel::SampleSecondaries(), G4DNAMillerGreenExcitationModel::SampleSecondaries(), G4DNAQuinnPlasmonExcitationModel::SampleSecondaries(), G4DNARelativisticIonisationModel::SampleSecondaries(), G4DNARuddIonisationExtendedModel::SampleSecondaries(), G4DNARuddIonisationModel::SampleSecondaries(), G4DNASancheExcitationModel::SampleSecondaries(), G4DNATransformElectronModel::SampleSecondaries(), G4BoldyshevTripletModel::SampleSecondaries(), G4eSingleCoulombScatteringModel::SampleSecondaries(), G4IonCoulombScatteringModel::SampleSecondaries(), G4PAIPhotModel::SampleSecondaries(), G4DNADingfelderChargeDecreaseModel::SampleSecondaries(), G4DNARPWBAExcitationModel::SampleSecondaries(), G4DNARPWBAIonisationModel::SampleSecondaries(), G4JAEAElasticScatteringModel::SampleSecondaries(), G4JAEAPolarizedElasticScatteringModel::SampleSecondaries(), G4LivermoreComptonModel::SampleSecondaries(), G4LivermoreIonisationModel::SampleSecondaries(), G4LivermorePhotoElectricModel::SampleSecondaries(), G4LivermorePolarizedComptonModel::SampleSecondaries(), G4LivermorePolarizedRayleighModel::SampleSecondaries(), G4LowEPComptonModel::SampleSecondaries(), G4LowEPPolarizedComptonModel::SampleSecondaries(), G4MicroElecElasticModel::SampleSecondaries(), G4MicroElecElasticModel_new::SampleSecondaries(), G4MicroElecInelasticModel::SampleSecondaries(), G4MicroElecInelasticModel_new::SampleSecondaries(), G4PenelopeBremsstrahlungModel::SampleSecondaries(), G4PenelopeComptonModel::SampleSecondaries(), G4PenelopeGammaConversionModel::SampleSecondaries(), G4PenelopeIonisationModel::SampleSecondaries(), G4PenelopePhotoElectricModel::SampleSecondaries(), G4PenelopeRayleighModel::SampleSecondaries(), G4PenelopeRayleighModelMI::SampleSecondaries(), G4PolarizedComptonModel::SampleSecondaries(), G4eCoulombScatteringModel::SampleSecondaries(), G4hCoulombScatteringModel::SampleSecondaries(), G4KleinNishinaCompton::SampleSecondaries(), G4KleinNishinaModel::SampleSecondaries(), G4PEEffectFluoModel::SampleSecondaries(), G4LEPTSAttachmentModel::SampleSecondaries(), G4LEPTSDissociationModel::SampleSecondaries(), G4LEPTSElasticModel::SampleSecondaries(), G4LEPTSExcitationModel::SampleSecondaries(), G4LEPTSIonisationModel::SampleSecondaries(), G4LEPTSPositroniumModel::SampleSecondaries(), G4LEPTSRotExcitationModel::SampleSecondaries(), G4LEPTSVibExcitationModel::SampleSecondaries(), G4DNAPTBElasticModel::SampleSecondaries(), G4DNAPTBExcitationModel::SampleSecondaries(), and G4DNAPTBIonisationModel::SampleSecondaries().

◆ ProposeNonIonizingEnergyDeposit()

◆ ProposeParentWeight()

◆ ProposeSteppingControl()

void G4VParticleChange::ProposeSteppingControl ( G4SteppingControl  StepControlFlag)
inline

◆ ProposeTrackStatus()

void G4VParticleChange::ProposeTrackStatus ( G4TrackStatus  status)
inline

Referenced by G4BiasingProcessInterface::AlongStepDoIt(), G4ITTransportation::AlongStepDoIt(), G4Transportation::AlongStepDoIt(), G4hImpactIonisation::AlongStepDoIt(), G4BOptnLeadingParticle::ApplyFinalStateBiasing(), G4HadronStoppingProcess::AtRestDoIt(), G4MuonMinusAtomicCapture::AtRestDoIt(), G4RadioactiveDecay::DecayAnalog(), G4DNAMolecularDissociation::DecayIt(), G4Decay::DecayIt(), G4UnknownDecay::DecayIt(), G4Radioactivation::DecayIt(), G4RadioactiveDecay::DecayIt(), G4DNABrownianTransportation::Diffusion(), G4HadronicProcess::FillResult(), G4FastStep::KillPrimaryTrack(), G4ImportanceProcess::KillTrack(), G4WeightWindowProcess::KillTrack(), G4SpecialCuts::PostStepDoIt(), G4UserSpecialCuts::PostStepDoIt(), G4WeightCutOffProcess::PostStepDoIt(), G4DNASecondOrderReaction::PostStepDoIt(), G4FastSimulationManagerProcess::PostStepDoIt(), G4PhononDownconversion::PostStepDoIt(), G4PhononReflection::PostStepDoIt(), G4PhononScattering::PostStepDoIt(), G4NeutronKiller::PostStepDoIt(), G4DNAScavengerProcess::PostStepDoIt(), G4LowECapture::PostStepDoIt(), G4VEmProcess::PostStepDoIt(), G4VEnergyLossProcess::PostStepDoIt(), G4NeutronGeneralProcess::PostStepDoIt(), G4NeutrinoElectronProcess::PostStepDoIt(), G4UCNAbsorption::PostStepDoIt(), G4UCNBoundaryProcess::PostStepDoIt(), G4UCNLoss::PostStepDoIt(), G4AnnihiToMuPair::PostStepDoIt(), G4GammaConversionToMuons::PostStepDoIt(), G4Cerenkov::PostStepDoIt(), G4Scintillation::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4TauNeutrinoNucleusProcess::PostStepDoIt(), G4OpAbsorption::PostStepDoIt(), G4OpBoundaryProcess::PostStepDoIt(), G4OpWLS::PostStepDoIt(), G4OpWLS2::PostStepDoIt(), G4ITTransportation::PostStepDoIt(), G4DNAPolyNucleotideReactionProcess::PostStepDoIt(), G4BiasingProcessInterface::PostStepDoIt(), G4hImpactIonisation::PostStepDoIt(), G4SynchrotronRadiationInMat::PostStepDoIt(), G4CoupledTransportation::PostStepDoIt(), G4Transportation::PostStepDoIt(), G4AdjointBremsstrahlungModel::RapidSampleSecondaries(), G4AdjointComptonModel::RapidSampleSecondaries(), G4AdjointhIonisationModel::RapidSampleSecondaries(), G4AdjointBremsstrahlungModel::SampleSecondaries(), G4AdjointComptonModel::SampleSecondaries(), G4AdjointeIonisationModel::SampleSecondaries(), G4AdjointhIonisationModel::SampleSecondaries(), G4AdjointIonIonisationModel::SampleSecondaries(), G4AdjointPhotoElectricModel::SampleSecondaries(), G4LivermoreBremsstrahlungModel::SampleSecondaries(), G4eBremParametrizedModel::SampleSecondaries(), G4eBremsstrahlungRelModel::SampleSecondaries(), G4SeltzerBergerModel::SampleSecondaries(), G4DNADingfelderChargeIncreaseModel::SampleSecondaries(), G4DNAELSEPAElasticModel::SampleSecondaries(), G4DNAIonElasticModel::SampleSecondaries(), G4DNAMeltonAttachmentModel::SampleSecondaries(), G4DNARuddIonisationExtendedModel::SampleSecondaries(), G4DNARuddIonisationModel::SampleSecondaries(), G4DNATransformElectronModel::SampleSecondaries(), G4BoldyshevTripletModel::SampleSecondaries(), G4PolarizedAnnihilationModel::SampleSecondaries(), G4DNADingfelderChargeDecreaseModel::SampleSecondaries(), G4JAEAElasticScatteringModel::SampleSecondaries(), G4JAEAPolarizedElasticScatteringModel::SampleSecondaries(), G4LivermoreComptonModel::SampleSecondaries(), G4LivermoreNuclearGammaConversionModel::SampleSecondaries(), G4LivermorePhotoElectricModel::SampleSecondaries(), G4LivermorePolarizedComptonModel::SampleSecondaries(), G4LivermorePolarizedGammaConversionModel::SampleSecondaries(), G4LivermorePolarizedRayleighModel::SampleSecondaries(), G4LowEPComptonModel::SampleSecondaries(), G4LowEPPolarizedComptonModel::SampleSecondaries(), G4MicroElecElasticModel::SampleSecondaries(), G4MicroElecElasticModel_new::SampleSecondaries(), G4PenelopeAnnihilationModel::SampleSecondaries(), G4PenelopeComptonModel::SampleSecondaries(), G4PenelopeGammaConversionModel::SampleSecondaries(), G4PenelopePhotoElectricModel::SampleSecondaries(), G4PenelopeRayleighModel::SampleSecondaries(), G4PenelopeRayleighModelMI::SampleSecondaries(), G4MuBremsstrahlungModel::SampleSecondaries(), G4MuonToMuonPairProductionModel::SampleSecondaries(), G4MuPairProductionModel::SampleSecondaries(), G4PolarizedComptonModel::SampleSecondaries(), G4BetheHeitlerModel::SampleSecondaries(), G4eeToTwoGammaModel::SampleSecondaries(), G4KleinNishinaCompton::SampleSecondaries(), G4KleinNishinaModel::SampleSecondaries(), G4PairProductionRelModel::SampleSecondaries(), G4PEEffectFluoModel::SampleSecondaries(), G4eplusTo2GammaOKVIModel::SampleSecondaries(), G4eplusTo3GammaOKVIModel::SampleSecondaries(), G4eeToHadronsMultiModel::SampleSecondaries(), G4LEPTSAttachmentModel::SampleSecondaries(), G4LEPTSElasticModel::SampleSecondaries(), G4LEPTSPositroniumModel::SampleSecondaries(), G4DNAPTBElasticModel::SampleSecondaries(), and G4BetheHeitler5DModel::SampleSecondaries().

◆ ProposeTrueStepLength()

◆ ProposeWeight()

◆ SetDebugFlag()

void G4VParticleChange::SetDebugFlag ( )
inline

◆ SetNumberOfSecondaries()

◆ SetParentWeightByProcess()

◆ SetSecondaryWeightByProcess()

◆ SetVerboseLevel()

void G4VParticleChange::SetVerboseLevel ( G4int  vLevel)
inline

◆ UpdateStepForAlongStep()

G4Step * G4VParticleChange::UpdateStepForAlongStep ( G4Step Step)
virtual

Reimplemented in G4ParticleChangeForNothing, G4ParticleChangeForOccurenceBiasing, G4ParticleChangeForLoss, G4ParticleChangeForMSC, G4ParticleChangeForTransport, and G4ParticleChange.

Definition at line 110 of file G4VParticleChange.cc.

111{
113 {
114 G4double initialWeight = Step->GetPreStepPoint()->GetWeight();
115 G4double currentWeight = Step->GetPostStepPoint()->GetWeight();
116 G4double finalWeight = (theParentWeight / initialWeight) * currentWeight;
117 Step->GetPostStepPoint()->SetWeight(finalWeight);
118 }
119 return UpdateStepInfo(Step);
120}
void SetWeight(G4double aValue)
G4double GetWeight() const
G4StepPoint * GetPreStepPoint() const
G4StepPoint * GetPostStepPoint() const
G4Step * UpdateStepInfo(G4Step *Step)

Referenced by G4ITStepProcessor::InvokeAlongStepDoItProcs(), and G4ParticleChangeForOccurenceBiasing::UpdateStepForAlongStep().

◆ UpdateStepForAtRest()

G4Step * G4VParticleChange::UpdateStepForAtRest ( G4Step Step)
virtual

◆ UpdateStepForPostStep()

◆ UpdateStepInfo()

G4Step * G4VParticleChange::UpdateStepInfo ( G4Step Step)
protected

Definition at line 71 of file G4VParticleChange.cc.

72{
73 // Update the G4Step specific attributes
74 pStep->SetStepLength(theTrueStepLength);
75 pStep->AddTotalEnergyDeposit(theLocalEnergyDeposit);
76 pStep->AddNonIonizingEnergyDeposit(theNonIonizingEnergyDeposit);
77 pStep->SetControlFlag(theSteppingControlFlag);
78
80 {
81 pStep->SetFirstStepFlag();
82 }
83 else
84 {
85 pStep->ClearFirstStepFlag();
86 }
88 {
89 pStep->SetLastStepFlag();
90 }
91 else
92 {
93 pStep->ClearLastStepFlag();
94 }
95
96 return pStep;
97}

Referenced by UpdateStepForAlongStep(), G4FastStep::UpdateStepForAtRest(), UpdateStepForAtRest(), G4ParticleChangeForDecay::UpdateStepForAtRest(), G4FastStep::UpdateStepForPostStep(), UpdateStepForPostStep(), and G4ParticleChangeForDecay::UpdateStepForPostStep().

Member Data Documentation

◆ accuracyForException

const G4double G4VParticleChange::accuracyForException = 0.001
staticprotected

Definition at line 228 of file G4VParticleChange.hh.

Referenced by GetAccuracyForException().

◆ accuracyForWarning

const G4double G4VParticleChange::accuracyForWarning = 1.0e-9
staticprotected

Definition at line 227 of file G4VParticleChange.hh.

Referenced by CheckSecondary(), and GetAccuracyForWarning().

◆ debugFlag

◆ fSetSecondaryWeightByProcess

G4bool G4VParticleChange::fSetSecondaryWeightByProcess = false
protected

Definition at line 283 of file G4VParticleChange.hh.

Referenced by AddSecondary().

◆ isParentWeightProposed

◆ maxError

const G4int G4VParticleChange::maxError = 10
staticprotected

◆ nError

G4int G4VParticleChange::nError = 0
protected

◆ theCurrentTrack

◆ theFirstStepInVolume

G4bool G4VParticleChange::theFirstStepInVolume = false
protected

Definition at line 277 of file G4VParticleChange.hh.

Referenced by DumpInfo(), and UpdateStepInfo().

◆ theLastStepInVolume

G4bool G4VParticleChange::theLastStepInVolume = false
protected

Definition at line 278 of file G4VParticleChange.hh.

Referenced by DumpInfo(), and UpdateStepInfo().

◆ theListOfSecondaries

std::vector<G4Track*> G4VParticleChange::theListOfSecondaries
protected

Definition at line 234 of file G4VParticleChange.hh.

Referenced by AddSecondary().

◆ theLocalEnergyDeposit

◆ theNonIonizingEnergyDeposit

◆ theNumberOfSecondaries

G4int G4VParticleChange::theNumberOfSecondaries = 0
protected

◆ theParentGlobalTime

G4double G4VParticleChange::theParentGlobalTime = 0.0
protected

Definition at line 262 of file G4VParticleChange.hh.

Referenced by CheckSecondary().

◆ theParentWeight

◆ theSizeOftheListOfSecondaries

G4int G4VParticleChange::theSizeOftheListOfSecondaries = 0
protected

Definition at line 269 of file G4VParticleChange.hh.

Referenced by AddSecondary().

◆ theStatusChange

G4TrackStatus G4VParticleChange::theStatusChange = fAlive
protected

◆ theSteppingControlFlag

G4SteppingControl G4VParticleChange::theSteppingControlFlag = NormalCondition
protected

◆ theTrueStepLength

G4double G4VParticleChange::theTrueStepLength = 0.0
protected

◆ verboseLevel

G4int G4VParticleChange::verboseLevel = 1
protected

Definition at line 272 of file G4VParticleChange.hh.

Referenced by G4FastStep::G4FastStep(), and G4FastStep::~G4FastStep().


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