Geant4-11
Public Member Functions | Protected Member Functions | Protected Attributes | Static Protected Attributes
G4VParticleChange Class Reference

#include <G4VParticleChange.hh>

Inheritance diagram for G4VParticleChange:
G4FastStep G4ParticleChange G4ParticleChangeForDecay G4ParticleChangeForGamma G4ParticleChangeForLoss G4ParticleChangeForMSC G4ParticleChangeForNothing G4ParticleChangeForOccurenceBiasing

Public Member Functions

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

Protected Member Functions

G4bool CheckSecondary (G4Track &)
 
 G4VParticleChange (const G4VParticleChange &right)
 
G4double GetAccuracyForException () const
 
G4double GetAccuracyForWarning () const
 
void InitializeLocalEnergyDeposit (const G4Track &)
 
void InitializeParentGlobalTime (const G4Track &)
 
void InitializeParentWeight (const G4Track &)
 
void InitializeSecondaries (const G4Track &)
 
void InitializeStatusChange (const G4Track &)
 
void InitializeStepInVolumeFlags (const G4Track &)
 
void InitializeSteppingControl (const G4Track &)
 
void InitializeTrueStepLength (const G4Track &)
 
G4VParticleChangeoperator= (const G4VParticleChange &right)
 
G4StepUpdateStepInfo (G4Step *Step)
 

Protected Attributes

G4bool debugFlag = false
 
G4bool fSetSecondaryWeightByProcess = false
 
G4bool isParentWeightProposed = false
 
G4bool theFirstStepInVolume = false
 
G4bool theLastStepInVolume = false
 
G4TrackFastVectortheListOfSecondaries = nullptr
 
G4double theLocalEnergyDeposit = 0.0
 
G4double theNonIonizingEnergyDeposit = 0.0
 
G4int theNumberOfSecondaries = 0
 
G4double theParentGlobalTime = 0.0
 
G4double theParentWeight = 1.0
 
G4int theSizeOftheListOfSecondaries = 0
 
G4TrackStatus theStatusChange = fAlive
 
G4SteppingControl theSteppingControlFlag = NormalCondition
 
G4double theTrueStepLength = 0.0
 
G4int verboseLevel = 1
 

Static Protected Attributes

static const G4double accuracyForException = 0.001
 
static const G4double accuracyForWarning = 1.0e-9
 

Detailed Description

Definition at line 70 of file G4VParticleChange.hh.

Constructor & Destructor Documentation

◆ G4VParticleChange() [1/2]

G4VParticleChange::G4VParticleChange ( )

Definition at line 42 of file G4VParticleChange.cc.

44{
45#ifdef G4VERBOSE
46 // activate CheckIt if in VERBOSE mode
47 debugFlag = true;
48#endif
50}
G4FastVector< G4Track, G4TrackFastVectorSize > G4TrackFastVector
const G4int G4TrackFastVectorSize
G4TrackFastVector * theListOfSecondaries

References debugFlag, and theListOfSecondaries.

◆ ~G4VParticleChange()

G4VParticleChange::~G4VParticleChange ( )
virtual

Definition at line 53 of file G4VParticleChange.cc.

54{
55 // check if tracks still exist in theListOfSecondaries
57 {
58 for(G4int index = 0; index < theNumberOfSecondaries; ++index)
59 {
60 delete(*theListOfSecondaries)[index];
61 }
62 }
64}
int G4int
Definition: G4Types.hh:85

References theListOfSecondaries, and theNumberOfSecondaries.

◆ G4VParticleChange() [2/2]

G4VParticleChange::G4VParticleChange ( const G4VParticleChange right)
protected

Definition at line 67 of file G4VParticleChange.cc.

79 , debugFlag(right.debugFlag)
80{
83 for(G4int index = 0; index < theNumberOfSecondaries; ++index)
84 {
85 G4Track* newTrack = new G4Track(*((*right.theListOfSecondaries)[index]));
86 theListOfSecondaries->SetElement(index, newTrack);
87 }
88}
void SetElement(G4int anIndex, Type *anElement)
Definition: G4FastVector.hh:72
G4TrackStatus theStatusChange
G4SteppingControl theSteppingControlFlag
G4double theNonIonizingEnergyDeposit

References G4FastVector< Type, N >::SetElement(), theListOfSecondaries, and theNumberOfSecondaries.

Member Function Documentation

◆ AddSecondary()

void G4VParticleChange::AddSecondary ( G4Track aSecondary)

Definition at line 146 of file G4VParticleChange.cc.

147{
148 if(debugFlag)
149 CheckSecondary(*aTrack);
150
151 // add a secondary after size check
153 {
154 // set weight of secondary tracks
156 aTrack->SetWeight(theParentWeight);
159 }
160 else
161 {
162 delete aTrack;
163
164 G4Exception("G4VParticleChange::AddSecondary()", "TRACK101", JustWarning,
165 "Secondary buffer is full. The track is deleted!");
166 }
167}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
G4bool CheckSecondary(G4Track &)

References CheckSecondary(), debugFlag, fSetSecondaryWeightByProcess, G4Exception(), JustWarning, G4FastVector< Type, N >::SetElement(), G4Track::SetWeight(), theListOfSecondaries, theNumberOfSecondaries, theParentWeight, and theSizeOftheListOfSecondaries.

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

◆ CheckIt()

G4bool G4VParticleChange::CheckIt ( const G4Track )
virtual

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

Definition at line 299 of file G4VParticleChange.cc.

304{
305 G4bool exitWithError = false;
306 G4double accuracy;
307 static G4ThreadLocal G4int nError = 0;
308#ifdef G4VERBOSE
309 const G4int maxError = 30;
310#endif
311
312 // Energy deposit should not be negative
313 G4bool itsOKforEnergy = true;
314 accuracy = -1.0 * theLocalEnergyDeposit / MeV;
315 if(accuracy > accuracyForWarning)
316 {
317 itsOKforEnergy = false;
318 nError += 1;
319 exitWithError = (accuracy > accuracyForException);
320#ifdef G4VERBOSE
321 if(nError < maxError)
322 {
323 G4cout << " G4VParticleChange::CheckIt : ";
324 G4cout << "the energy deposit is negative !!"
325 << " Difference: " << accuracy << "[MeV] " << G4endl;
326 G4cout << aTrack.GetDefinition()->GetParticleName()
327 << " E=" << aTrack.GetKineticEnergy() / MeV
328 << " pos=" << aTrack.GetPosition().x() / m << ", "
329 << aTrack.GetPosition().y() / m << ", "
330 << aTrack.GetPosition().z() / m << G4endl;
331 }
332#endif
333 }
334
335 // true path length should not be negative
336 G4bool itsOKforStepLength = true;
337 accuracy = -1.0 * theTrueStepLength / mm;
338 if(accuracy > accuracyForWarning)
339 {
340 itsOKforStepLength = false;
341 nError += 1;
342 exitWithError = (accuracy > accuracyForException);
343#ifdef G4VERBOSE
344 if(nError < maxError)
345 {
346 G4cout << " G4VParticleChange::CheckIt : ";
347 G4cout << "the true step length is negative !!"
348 << " Difference: " << accuracy << "[MeV] " << G4endl;
349 G4cout << aTrack.GetDefinition()->GetParticleName()
350 << " E=" << aTrack.GetKineticEnergy() / MeV
351 << " pos=" << aTrack.GetPosition().x() / m << ", "
352 << aTrack.GetPosition().y() / m << ", "
353 << aTrack.GetPosition().z() / m << G4endl;
354 }
355#endif
356 }
357#ifdef G4VERBOSE
358 if(!itsOKforStepLength || !itsOKforEnergy)
359 {
360 // dump out information of this particle change
361 DumpInfo();
362 }
363#endif
364
365 // Exit with error
366 if(exitWithError)
367 {
368 G4Exception("G4VParticleChange::CheckIt()", "TRACK001", EventMustBeAborted,
369 "Step length and/or energy deposit was illegal");
370 }
371
372 // correction
373 if(!itsOKforStepLength)
374 {
375 theTrueStepLength = (1.e-12) * mm;
376 }
377 if(!itsOKforEnergy)
378 {
380 }
381 return (itsOKforStepLength && itsOKforEnergy);
382}
@ EventMustBeAborted
static constexpr double m
Definition: G4SIunits.hh:109
static constexpr double mm
Definition: G4SIunits.hh:95
static constexpr double MeV
Definition: G4SIunits.hh:200
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static const G4double accuracyForException
static const G4double accuracyForWarning
virtual void DumpInfo() const
#define G4ThreadLocal
Definition: tls.hh:77

References accuracyForException, accuracyForWarning, DumpInfo(), EventMustBeAborted, G4cout, G4endl, G4Exception(), G4ThreadLocal, m, MeV, mm, theLocalEnergyDeposit, and theTrueStepLength.

Referenced by G4FastStep::CheckIt(), G4ParticleChange::CheckIt(), G4ParticleChangeForDecay::CheckIt(), G4ParticleChangeForGamma::CheckIt(), G4ParticleChangeForLoss::CheckIt(), and G4ParticleChangeForMSC::CheckIt().

◆ CheckSecondary()

G4bool G4VParticleChange::CheckSecondary ( G4Track aTrack)
protected

Definition at line 385 of file G4VParticleChange.cc.

386{
387 G4bool exitWithError = false;
388 G4double accuracy;
389 static G4ThreadLocal G4int nError = 0;
390#ifdef G4VERBOSE
391 const G4int maxError = 30;
392#endif
393
394 // MomentumDirection should be unit vector
395 G4bool itsOKforMomentum = true;
396 if(aTrack.GetKineticEnergy() > 0.)
397 {
398 accuracy = std::fabs((aTrack.GetMomentumDirection()).mag2() - 1.0);
399 if(accuracy > accuracyForWarning)
400 {
401 itsOKforMomentum = false;
402 nError += 1;
403 exitWithError = exitWithError || (accuracy > accuracyForException);
404#ifdef G4VERBOSE
405 if(nError < maxError)
406 {
407 G4cout << " G4VParticleChange::CheckSecondary : ";
408 G4cout << "the Momentum direction is not unit vector !! "
409 << " Difference: " << accuracy << G4endl;
411 << " E=" << aTrack.GetKineticEnergy() / MeV
412 << " pos=" << aTrack.GetPosition().x() / m << ", "
413 << aTrack.GetPosition().y() / m << ", "
414 << aTrack.GetPosition().z() / m << G4endl;
415 }
416#endif
417 }
418 }
419
420 // Kinetic Energy should not be negative
421 G4bool itsOKforEnergy = true;
422 accuracy = -1.0 * (aTrack.GetKineticEnergy()) / MeV;
423 if(accuracy > accuracyForWarning)
424 {
425 itsOKforEnergy = false;
426 nError += 1;
427 exitWithError = exitWithError || (accuracy > accuracyForException);
428#ifdef G4VERBOSE
429 if(nError < maxError)
430 {
431 G4cout << " G4VParticleChange::CheckSecondary : ";
432 G4cout << "the kinetic energy is negative !!"
433 << " Difference: " << accuracy << "[MeV] " << G4endl;
434 G4cout << " G4VParticleChange::CheckSecondary : ";
435 G4cout << "the global time of secondary is earlier than the parent !!"
436 << " Difference: " << accuracy << "[ns] " << G4endl;
438 << " E=" << aTrack.GetKineticEnergy() / MeV
439 << " pos=" << aTrack.GetPosition().x() / m << ", "
440 << aTrack.GetPosition().y() / m << ", "
441 << aTrack.GetPosition().z() / m << G4endl;
442 }
443#endif
444 }
445 // Check timing of secondaries
446 G4bool itsOKforTiming = true;
447
448 accuracy = (theParentGlobalTime - aTrack.GetGlobalTime()) / ns;
449 if(accuracy > accuracyForWarning)
450 {
451 itsOKforTiming = false;
452 nError += 1;
453 exitWithError = (accuracy > accuracyForException);
454#ifdef G4VERBOSE
455 if(nError < maxError)
456 {
457 G4cout << " G4VParticleChange::CheckSecondary : ";
458 G4cout
459 << "the global time of secondary goes back comapared to the parent !!"
460 << " Difference: " << accuracy << "[ns] " << G4endl;
462 << " E=" << aTrack.GetKineticEnergy() / MeV
463 << " pos=" << aTrack.GetPosition().x() / m << ", "
464 << aTrack.GetPosition().y() / m << ", "
465 << aTrack.GetPosition().z() / m
466 << " time=" << aTrack.GetGlobalTime() / ns
467 << " parent time=" << theParentGlobalTime / ns << G4endl;
468 }
469#endif
470 }
471
472 // Exit with error
473 if(exitWithError)
474 {
475 G4Exception("G4VParticleChange::CheckSecondary()", "TRACK001",
476 EventMustBeAborted, "Secondary with illegal energy/momentum ");
477 }
478
479 G4bool itsOK = itsOKforMomentum && itsOKforEnergy && itsOKforTiming;
480
481 // correction
482 if(!itsOKforMomentum)
483 {
484 G4double vmag = (aTrack.GetMomentumDirection()).mag();
485 aTrack.SetMomentumDirection((1. / vmag) * aTrack.GetMomentumDirection());
486 }
487 if(!itsOKforEnergy)
488 {
489 aTrack.SetKineticEnergy(0.0);
490 }
491
492 if(!itsOK)
493 {
494 this->DumpInfo();
495 }
496
497 return itsOK;
498}
double z() const
double x() const
double y() const
const G4String & GetParticleName() const
const G4ThreeVector & GetPosition() const
G4double GetGlobalTime() const
G4ParticleDefinition * GetDefinition() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
void SetKineticEnergy(const G4double aValue)
void SetMomentumDirection(const G4ThreeVector &aValue)
#define ns
Definition: xmlparse.cc:614

References accuracyForException, accuracyForWarning, DumpInfo(), EventMustBeAborted, G4cout, G4endl, G4Exception(), G4ThreadLocal, G4Track::GetDefinition(), G4Track::GetGlobalTime(), G4Track::GetKineticEnergy(), G4Track::GetMomentumDirection(), G4ParticleDefinition::GetParticleName(), G4Track::GetPosition(), m, MeV, ns, G4Track::SetKineticEnergy(), G4Track::SetMomentumDirection(), theParentGlobalTime, CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

Referenced by AddSecondary().

◆ Clear()

void G4VParticleChange::Clear ( )

◆ ClearDebugFlag()

void G4VParticleChange::ClearDebugFlag ( )

◆ DumpInfo()

void G4VParticleChange::DumpInfo ( ) const
virtual

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

Definition at line 232 of file G4VParticleChange.cc.

233{
234 // Show header
235 G4int olprc = G4cout.precision(3);
236 G4cout << " -----------------------------------------------" << G4endl;
237 G4cout << " G4ParticleChange Information " << std::setw(20) << G4endl;
238 G4cout << " -----------------------------------------------" << G4endl;
239
240 G4cout << " # of secondaries : " << std::setw(20)
242
244 {
245 G4cout << " Pointer to 2ndaries : " << std::setw(20)
246 << GetSecondary(0) << G4endl;
247 G4cout << " (Showed only 1st one)" << G4endl;
248 }
249 G4cout << " -----------------------------------------------" << G4endl;
250
251 G4cout << " Energy Deposit (MeV): " << std::setw(20)
253
254 G4cout << " Non-ionizing Energy Deposit (MeV): " << std::setw(11)
256
257 G4cout << " Track Status : " << std::setw(20);
259 {
260 G4cout << " Alive";
261 }
263 {
264 G4cout << " StopButAlive";
265 }
266 else if(theStatusChange == fStopAndKill)
267 {
268 G4cout << " StopAndKill";
269 }
271 {
272 G4cout << " KillTrackAndSecondaries";
273 }
274 else if(theStatusChange == fSuspend)
275 {
276 G4cout << " Suspend";
277 }
279 {
280 G4cout << " PostponeToNextEvent";
281 }
282 G4cout << G4endl;
283 G4cout << " True Path Length (mm) : " << std::setw(18)
285 G4cout << " Stepping Control : " << std::setw(20)
288 {
289 G4cout << " First step in volume" << G4endl;
290 }
292 {
293 G4cout << " Last step in volume" << G4endl;
294 }
295 G4cout.precision(olprc);
296}
@ fKillTrackAndSecondaries
@ fSuspend
@ fAlive
@ fStopAndKill
@ fStopButAlive
@ fPostponeToNextEvent
G4Track * GetSecondary(G4int anIndex) const

References fAlive, fKillTrackAndSecondaries, fPostponeToNextEvent, fStopAndKill, fStopButAlive, fSuspend, G4cout, G4endl, GetSecondary(), MeV, mm, theFirstStepInVolume, theLastStepInVolume, theLocalEnergyDeposit, theNonIonizingEnergyDeposit, theNumberOfSecondaries, theStatusChange, theSteppingControlFlag, and theTrueStepLength.

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

◆ GetAccuracyForException()

G4double G4VParticleChange::GetAccuracyForException ( ) const
protected

Definition at line 507 of file G4VParticleChange.cc.

508{
510}

References accuracyForException.

Referenced by G4FastStep::CheckIt().

◆ GetAccuracyForWarning()

G4double G4VParticleChange::GetAccuracyForWarning ( ) const
protected

Definition at line 501 of file G4VParticleChange.cc.

502{
503 return accuracyForWarning;
504}

References accuracyForWarning.

Referenced by G4FastStep::CheckIt().

◆ GetDebugFlag()

G4bool G4VParticleChange::GetDebugFlag ( ) const

◆ GetFirstStepInVolume()

G4bool G4VParticleChange::GetFirstStepInVolume ( ) const

◆ GetLastStepInVolume()

G4bool G4VParticleChange::GetLastStepInVolume ( ) const

◆ GetLocalEnergyDeposit()

G4double G4VParticleChange::GetLocalEnergyDeposit ( ) const

◆ GetNonIonizingEnergyDeposit()

G4double G4VParticleChange::GetNonIonizingEnergyDeposit ( ) const

◆ GetNumberOfSecondaries()

G4int G4VParticleChange::GetNumberOfSecondaries ( ) const

◆ GetParentWeight()

G4double G4VParticleChange::GetParentWeight ( ) const

◆ GetSecondary()

G4Track * G4VParticleChange::GetSecondary ( G4int  anIndex) const

◆ GetSteppingControl()

G4SteppingControl G4VParticleChange::GetSteppingControl ( ) const

◆ GetTrackStatus()

G4TrackStatus G4VParticleChange::GetTrackStatus ( ) const

◆ GetTrueStepLength()

G4double G4VParticleChange::GetTrueStepLength ( ) const

◆ GetVerboseLevel()

G4int G4VParticleChange::GetVerboseLevel ( ) const

◆ GetWeight()

G4double G4VParticleChange::GetWeight ( ) const

◆ Initialize()

virtual void G4VParticleChange::Initialize ( const G4Track )
virtual

◆ InitializeLocalEnergyDeposit()

void G4VParticleChange::InitializeLocalEnergyDeposit ( const G4Track )
protected

◆ InitializeParentGlobalTime()

void G4VParticleChange::InitializeParentGlobalTime ( const G4Track )
protected

◆ InitializeParentWeight()

void G4VParticleChange::InitializeParentWeight ( const G4Track )
protected

◆ InitializeSecondaries()

void G4VParticleChange::InitializeSecondaries ( const G4Track )
protected

◆ InitializeStatusChange()

void G4VParticleChange::InitializeStatusChange ( const G4Track )
protected

◆ InitializeStepInVolumeFlags()

void G4VParticleChange::InitializeStepInVolumeFlags ( const G4Track )
protected

◆ InitializeSteppingControl()

void G4VParticleChange::InitializeSteppingControl ( const G4Track )
protected

◆ InitializeTrueStepLength()

void G4VParticleChange::InitializeTrueStepLength ( const G4Track )
protected

◆ IsParentWeightSetByProcess()

G4bool G4VParticleChange::IsParentWeightSetByProcess ( ) const

Definition at line 516 of file G4VParticleChange.cc.

516{ return true; }

◆ IsSecondaryWeightSetByProcess()

G4bool G4VParticleChange::IsSecondaryWeightSetByProcess ( ) const

◆ operator!=()

G4bool G4VParticleChange::operator!= ( const G4VParticleChange right) const

Definition at line 140 of file G4VParticleChange.cc.

141{
142 return (this != (G4VParticleChange*) &right);
143}

◆ operator=()

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

Definition at line 91 of file G4VParticleChange.cc.

92{
93 if(this != &right)
94 {
96 {
97 for(G4int index = 0; index < theNumberOfSecondaries; ++index)
98 {
99 if((*theListOfSecondaries)[index])
100 delete((*theListOfSecondaries)[index]);
101 }
102 }
104
107 for(G4int index = 0; index < theNumberOfSecondaries; ++index)
108 {
109 G4Track* newTrack = new G4Track(*((*right.theListOfSecondaries)[index]));
110 theListOfSecondaries->SetElement(index, newTrack);
111 }
117
120
124
126
128 debugFlag = right.debugFlag;
129 }
130 return *this;
131}

References debugFlag, fSetSecondaryWeightByProcess, isParentWeightProposed, G4FastVector< Type, N >::SetElement(), theFirstStepInVolume, theLastStepInVolume, theListOfSecondaries, theLocalEnergyDeposit, theNonIonizingEnergyDeposit, theNumberOfSecondaries, theParentGlobalTime, theParentWeight, theStatusChange, theSteppingControlFlag, theTrueStepLength, and verboseLevel.

Referenced by G4FastStep::operator=().

◆ operator==()

G4bool G4VParticleChange::operator== ( const G4VParticleChange right) const

Definition at line 134 of file G4VParticleChange.cc.

135{
136 return (this == (G4VParticleChange*) &right);
137}

◆ ProposeFirstStepInVolume()

void G4VParticleChange::ProposeFirstStepInVolume ( G4bool  flag)

◆ ProposeLastStepInVolume()

void G4VParticleChange::ProposeLastStepInVolume ( G4bool  flag)

◆ ProposeLocalEnergyDeposit()

void G4VParticleChange::ProposeLocalEnergyDeposit ( G4double  anEnergyPart)

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(), G4MuonicAtomDecay::DecayIt(), G4OpBoundaryProcess::DoAbsorption(), G4HadronicProcess::FillResult(), G4MuonicAtomDecay::FillResult(), G4SpecialCuts::PostStepDoIt(), G4UserSpecialCuts::PostStepDoIt(), G4LowECapture::PostStepDoIt(), G4VEmProcess::PostStepDoIt(), G4VEnergyLossProcess::PostStepDoIt(), G4NeutrinoElectronProcess::PostStepDoIt(), G4UCNBoundaryProcess::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4OpAbsorption::PostStepDoIt(), G4OpBoundaryProcess::PostStepDoIt(), G4hImpactIonisation::PostStepDoIt(), G4SynchrotronRadiationInMat::PostStepDoIt(), G4DNABornExcitationModel1::SampleSecondaries(), G4DNABornExcitationModel2::SampleSecondaries(), G4DNABornIonisationModel1::SampleSecondaries(), G4DNABornIonisationModel2::SampleSecondaries(), G4DNACPA100ElasticModel::SampleSecondaries(), G4DNACPA100ExcitationModel::SampleSecondaries(), G4DNACPA100IonisationModel::SampleSecondaries(), G4DNADingfelderChargeDecreaseModel::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(), 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()

void G4VParticleChange::ProposeNonIonizingEnergyDeposit ( G4double  anEnergyPart)

◆ ProposeParentWeight()

void G4VParticleChange::ProposeParentWeight ( G4double  finalWeight)

◆ ProposeSteppingControl()

void G4VParticleChange::ProposeSteppingControl ( G4SteppingControl  StepControlFlag)

◆ ProposeTrackStatus()

void G4VParticleChange::ProposeTrackStatus ( G4TrackStatus  status)

Referenced by G4BiasingProcessInterface::AlongStepDoIt(), G4ITTransportation::AlongStepDoIt(), G4CoupledTransportation::AlongStepDoIt(), G4Transportation::AlongStepDoIt(), G4hImpactIonisation::AlongStepDoIt(), G4BOptnLeadingParticle::ApplyFinalStateBiasing(), G4HadronStoppingProcess::AtRestDoIt(), G4MuonMinusAtomicCapture::AtRestDoIt(), G4RadioactiveDecay::DecayAnalog(), G4DNAMolecularDissociation::DecayIt(), G4Decay::DecayIt(), G4UnknownDecay::DecayIt(), G4Radioactivation::DecayIt(), G4RadioactiveDecay::DecayIt(), G4MuonicAtomDecay::DecayIt(), G4DNABrownianTransportation::Diffusion(), G4OpBoundaryProcess::DoAbsorption(), G4HadronicProcess::FillResult(), G4MuonicAtomDecay::FillResult(), G4FastStep::KillPrimaryTrack(), G4ImportanceProcess::KillTrack(), G4WeightWindowProcess::KillTrack(), G4DNAElectronHoleRecombination::MakeReaction(), G4SpecialCuts::PostStepDoIt(), G4WeightCutOffProcess::PostStepDoIt(), G4DNASecondOrderReaction::PostStepDoIt(), G4FastSimulationManagerProcess::PostStepDoIt(), G4PhononDownconversion::PostStepDoIt(), G4PhononReflection::PostStepDoIt(), G4PhononScattering::PostStepDoIt(), G4NeutronKiller::PostStepDoIt(), G4UserSpecialCuts::PostStepDoIt(), G4DNAScavengerProcess::PostStepDoIt(), G4LowECapture::PostStepDoIt(), G4VEmProcess::PostStepDoIt(), G4VEnergyLossProcess::PostStepDoIt(), G4NeutrinoElectronProcess::PostStepDoIt(), G4UCNAbsorption::PostStepDoIt(), G4UCNBoundaryProcess::PostStepDoIt(), G4UCNLoss::PostStepDoIt(), G4AnnihiToMuPair::PostStepDoIt(), G4GammaConversionToMuons::PostStepDoIt(), G4Cerenkov::PostStepDoIt(), G4Scintillation::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4OpAbsorption::PostStepDoIt(), G4OpBoundaryProcess::PostStepDoIt(), G4OpWLS::PostStepDoIt(), G4OpWLS2::PostStepDoIt(), G4ITTransportation::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(), G4DNADingfelderChargeDecreaseModel::SampleSecondaries(), G4DNADingfelderChargeIncreaseModel::SampleSecondaries(), G4DNAELSEPAElasticModel::SampleSecondaries(), G4DNAIonElasticModel::SampleSecondaries(), G4DNAMeltonAttachmentModel::SampleSecondaries(), G4DNARuddIonisationExtendedModel::SampleSecondaries(), G4DNARuddIonisationModel::SampleSecondaries(), G4DNATransformElectronModel::SampleSecondaries(), G4BoldyshevTripletModel::SampleSecondaries(), G4PolarizedAnnihilationModel::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(), 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()

void G4VParticleChange::ProposeTrueStepLength ( G4double  truePathLength)

◆ ProposeWeight()

void G4VParticleChange::ProposeWeight ( G4double  finalWeight)

◆ SetDebugFlag()

void G4VParticleChange::SetDebugFlag ( )

◆ SetNumberOfSecondaries()

void G4VParticleChange::SetNumberOfSecondaries ( G4int  totSecondaries)

◆ SetParentWeightByProcess()

void G4VParticleChange::SetParentWeightByProcess ( G4bool  )

◆ SetSecondaryWeightByProcess()

void G4VParticleChange::SetSecondaryWeightByProcess ( G4bool  )

◆ SetVerboseLevel()

void G4VParticleChange::SetVerboseLevel ( G4int  vLevel)

◆ UpdateStepForAlongStep()

G4Step * G4VParticleChange::UpdateStepForAlongStep ( G4Step Step)
virtual

◆ UpdateStepForAtRest()

G4Step * G4VParticleChange::UpdateStepForAtRest ( G4Step Step)
virtual

◆ UpdateStepForPostStep()

G4Step * G4VParticleChange::UpdateStepForPostStep ( G4Step Step)
virtual

◆ UpdateStepInfo()

G4Step * G4VParticleChange::UpdateStepInfo ( G4Step Step)
protected

Definition at line 170 of file G4VParticleChange.cc.

171{
172 // Update the G4Step specific attributes
173 pStep->SetStepLength(theTrueStepLength);
174 pStep->AddTotalEnergyDeposit(theLocalEnergyDeposit);
175 pStep->AddNonIonizingEnergyDeposit(theNonIonizingEnergyDeposit);
176 pStep->SetControlFlag(theSteppingControlFlag);
177
179 {
180 pStep->SetFirstStepFlag();
181 }
182 else
183 {
184 pStep->ClearFirstStepFlag();
185 }
187 {
188 pStep->SetLastStepFlag();
189 }
190 else
191 {
192 pStep->ClearLastStepFlag();
193 }
194
195 return pStep;
196}

References G4Step::AddNonIonizingEnergyDeposit(), G4Step::AddTotalEnergyDeposit(), G4Step::ClearFirstStepFlag(), G4Step::ClearLastStepFlag(), G4Step::SetControlFlag(), G4Step::SetFirstStepFlag(), G4Step::SetLastStepFlag(), G4Step::SetStepLength(), theFirstStepInVolume, theLastStepInVolume, theLocalEnergyDeposit, theNonIonizingEnergyDeposit, theSteppingControlFlag, and theTrueStepLength.

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

Field Documentation

◆ accuracyForException

const G4double G4VParticleChange::accuracyForException = 0.001
staticprotected

◆ accuracyForWarning

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

◆ debugFlag

G4bool G4VParticleChange::debugFlag = false
protected

◆ fSetSecondaryWeightByProcess

G4bool G4VParticleChange::fSetSecondaryWeightByProcess = false
protected

◆ isParentWeightProposed

G4bool G4VParticleChange::isParentWeightProposed = false
protected

◆ theFirstStepInVolume

G4bool G4VParticleChange::theFirstStepInVolume = false
protected

Definition at line 277 of file G4VParticleChange.hh.

Referenced by DumpInfo(), operator=(), and UpdateStepInfo().

◆ theLastStepInVolume

G4bool G4VParticleChange::theLastStepInVolume = false
protected

Definition at line 278 of file G4VParticleChange.hh.

Referenced by DumpInfo(), operator=(), and UpdateStepInfo().

◆ theListOfSecondaries

G4TrackFastVector* G4VParticleChange::theListOfSecondaries = nullptr
protected

◆ theLocalEnergyDeposit

G4double G4VParticleChange::theLocalEnergyDeposit = 0.0
protected

◆ theNonIonizingEnergyDeposit

G4double G4VParticleChange::theNonIonizingEnergyDeposit = 0.0
protected

◆ theNumberOfSecondaries

G4int G4VParticleChange::theNumberOfSecondaries = 0
protected

◆ theParentGlobalTime

G4double G4VParticleChange::theParentGlobalTime = 0.0
protected

Definition at line 264 of file G4VParticleChange.hh.

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

◆ theParentWeight

G4double G4VParticleChange::theParentWeight = 1.0
protected

◆ theSizeOftheListOfSecondaries

G4int G4VParticleChange::theSizeOftheListOfSecondaries = 0
protected

◆ 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

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