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

#include <G4ParticleChangeForOccurenceBiasing.hh>

+ Inheritance diagram for G4ParticleChangeForOccurenceBiasing:

Public Member Functions

 G4ParticleChangeForOccurenceBiasing (G4String name)
 
 ~G4ParticleChangeForOccurenceBiasing ()
 
void SetOccurenceWeightForNonInteraction (G4double w)
 
G4double GetOccurenceWeightForNonInteraction () const
 
void SetOccurenceWeightForInteraction (G4double w)
 
G4double GetOccurenceWeightForInteraction () const
 
void SetWrappedParticleChange (G4VParticleChange *wpc)
 
G4VParticleChangeGetWrappedParticleChange () const
 
void StealSecondaries ()
 
virtual G4StepUpdateStepForAtRest (G4Step *step)
 
virtual G4StepUpdateStepForAlongStep (G4Step *step)
 
virtual G4StepUpdateStepForPostStep (G4Step *step)
 
const G4StringGetName () const
 
- Public Member Functions inherited from G4VParticleChange
 G4VParticleChange ()
 
virtual ~G4VParticleChange ()=default
 
 G4VParticleChange (const G4VParticleChange &right)=delete
 
G4VParticleChangeoperator= (const G4VParticleChange &right)=delete
 
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
 

Additional Inherited Members

- Protected Member Functions inherited from G4VParticleChange
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 inherited from G4VParticleChange
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 inherited from G4VParticleChange
static const G4double accuracyForWarning = 1.0e-9
 
static const G4double accuracyForException = 0.001
 
static const G4int maxError = 10
 

Detailed Description

Definition at line 47 of file G4ParticleChangeForOccurenceBiasing.hh.

Constructor & Destructor Documentation

◆ G4ParticleChangeForOccurenceBiasing()

G4ParticleChangeForOccurenceBiasing::G4ParticleChangeForOccurenceBiasing ( G4String name)

Definition at line 28 of file G4ParticleChangeForOccurenceBiasing.cc.

30 fName(name),
31 fWrappedParticleChange(0),
32 fOccurenceWeightForNonInteraction(-1.0),
33 fOccurenceWeightForInteraction(-1.0)
34{}

◆ ~G4ParticleChangeForOccurenceBiasing()

G4ParticleChangeForOccurenceBiasing::~G4ParticleChangeForOccurenceBiasing ( )

Definition at line 36 of file G4ParticleChangeForOccurenceBiasing.cc.

37{}

Member Function Documentation

◆ GetName()

const G4String & G4ParticleChangeForOccurenceBiasing::GetName ( ) const
inline

Definition at line 73 of file G4ParticleChangeForOccurenceBiasing.hh.

73{return fName;}

◆ GetOccurenceWeightForInteraction()

G4double G4ParticleChangeForOccurenceBiasing::GetOccurenceWeightForInteraction ( ) const
inline

Definition at line 56 of file G4ParticleChangeForOccurenceBiasing.hh.

56{return fOccurenceWeightForInteraction;}

◆ GetOccurenceWeightForNonInteraction()

G4double G4ParticleChangeForOccurenceBiasing::GetOccurenceWeightForNonInteraction ( ) const
inline

Definition at line 54 of file G4ParticleChangeForOccurenceBiasing.hh.

54{return fOccurenceWeightForNonInteraction;}

◆ GetWrappedParticleChange()

G4VParticleChange * G4ParticleChangeForOccurenceBiasing::GetWrappedParticleChange ( ) const
inline

Definition at line 61 of file G4ParticleChangeForOccurenceBiasing.hh.

61{return fWrappedParticleChange;}

◆ SetOccurenceWeightForInteraction()

void G4ParticleChangeForOccurenceBiasing::SetOccurenceWeightForInteraction ( G4double w)
inline

Definition at line 55 of file G4ParticleChangeForOccurenceBiasing.hh.

55{fOccurenceWeightForInteraction = w;}

Referenced by G4BiasingProcessInterface::PostStepDoIt().

◆ SetOccurenceWeightForNonInteraction()

void G4ParticleChangeForOccurenceBiasing::SetOccurenceWeightForNonInteraction ( G4double w)
inline

Definition at line 53 of file G4ParticleChangeForOccurenceBiasing.hh.

53{fOccurenceWeightForNonInteraction = w;}

Referenced by G4BiasingProcessInterface::AlongStepDoIt().

◆ SetWrappedParticleChange()

void G4ParticleChangeForOccurenceBiasing::SetWrappedParticleChange ( G4VParticleChange * wpc)

Definition at line 39 of file G4ParticleChangeForOccurenceBiasing.cc.

40{
41 fWrappedParticleChange = wpc;
42}

Referenced by G4BiasingProcessInterface::AlongStepDoIt(), and G4BiasingProcessInterface::PostStepDoIt().

◆ StealSecondaries()

void G4ParticleChangeForOccurenceBiasing::StealSecondaries ( )

Definition at line 44 of file G4ParticleChangeForOccurenceBiasing.cc.

45{
46 SetNumberOfSecondaries( fWrappedParticleChange->GetNumberOfSecondaries() );
47 for (G4int isecond = 0; isecond < fWrappedParticleChange->GetNumberOfSecondaries(); isecond++)
48 {
49 G4Track* secondary = fWrappedParticleChange->GetSecondary(isecond);
50 secondary->SetWeight ( secondary->GetWeight() * fOccurenceWeightForInteraction );
51 AddSecondary( secondary );
52 }
53 fWrappedParticleChange->Clear();
54}
int G4int
Definition G4Types.hh:85
G4double GetWeight() const
void SetWeight(G4double aValue)
G4int GetNumberOfSecondaries() const
void AddSecondary(G4Track *aSecondary)
void SetNumberOfSecondaries(G4int totSecondaries)
G4Track * GetSecondary(G4int anIndex) const

Referenced by G4BiasingProcessInterface::PostStepDoIt().

◆ UpdateStepForAlongStep()

G4Step * G4ParticleChangeForOccurenceBiasing::UpdateStepForAlongStep ( G4Step * step)
virtual

Reimplemented from G4VParticleChange.

Definition at line 62 of file G4ParticleChangeForOccurenceBiasing.cc.

63{
64 // -- make particle change of wrapped process to apply its changes:
65 if ( fWrappedParticleChange ) fWrappedParticleChange->UpdateStepForAlongStep( step );
66
67 // -- multiply parent weight by weight due to occurrence biasing:
68 G4StepPoint* postStepPoint = step->GetPostStepPoint();
69 postStepPoint->SetWeight( postStepPoint->GetWeight() * fOccurenceWeightForNonInteraction );
70
71 return step;
72}
void SetWeight(G4double aValue)
G4double GetWeight() const
G4StepPoint * GetPostStepPoint() const
virtual G4Step * UpdateStepForAlongStep(G4Step *Step)

◆ UpdateStepForAtRest()

G4Step * G4ParticleChangeForOccurenceBiasing::UpdateStepForAtRest ( G4Step * step)
virtual

Reimplemented from G4VParticleChange.

Definition at line 57 of file G4ParticleChangeForOccurenceBiasing.cc.

58{
59 return step;
60}

◆ UpdateStepForPostStep()

G4Step * G4ParticleChangeForOccurenceBiasing::UpdateStepForPostStep ( G4Step * step)
virtual

Reimplemented from G4VParticleChange.

Definition at line 74 of file G4ParticleChangeForOccurenceBiasing.cc.

75{
76 // -- let make first wrapped process to apply its changes:
77 fWrappedParticleChange->UpdateStepForPostStep(step);
78 // -- then apply weight correction due to occurrence biasing:
79 G4StepPoint* postStepPoint = step->GetPostStepPoint();
80 postStepPoint->SetWeight( postStepPoint->GetWeight() * fOccurenceWeightForInteraction );
81
82 return step;
83}
virtual G4Step * UpdateStepForPostStep(G4Step *Step)

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