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

#include <G4BOptnForceCommonTruncatedExp.hh>

+ Inheritance diagram for G4BOptnForceCommonTruncatedExp:

Public Member Functions

 G4BOptnForceCommonTruncatedExp (const G4String &name)
 
virtual ~G4BOptnForceCommonTruncatedExp ()
 
virtual const G4VBiasingInteractionLawProvideOccurenceBiasingInteractionLaw (const G4BiasingProcessInterface *, G4ForceCondition &)
 
virtual G4double ProposeAlongStepLimit (const G4BiasingProcessInterface *)
 
virtual G4GPILSelection ProposeGPILSelection (const G4GPILSelection processSelection)
 
virtual G4VParticleChangeApplyFinalStateBiasing (const G4BiasingProcessInterface *, const G4Track *, const G4Step *, G4bool &)
 
virtual G4double DistanceToApplyOperation (const G4Track *, G4double, G4ForceCondition *)
 
virtual G4VParticleChangeGenerateBiasingFinalState (const G4Track *, const G4Step *)
 
G4ILawCommonTruncatedExpGetCommonTruncatedExpLaw ()
 
G4ILawForceFreeFlightGetForceFreeFlightLaw ()
 
void Initialize (const G4Track *)
 
void UpdateForStep (const G4Step *)
 
void Sample ()
 
const G4ThreeVectorGetInitialMomentum () const
 
G4double GetMaximumDistance () const
 
void ChooseProcessToApply ()
 
const G4VProcessGetProcessToApply () const
 
void AddCrossSection (const G4VProcess *, G4double)
 
std::size_t GetNumberOfSharing () const
 
void SetInteractionOccured (G4bool b)
 
G4bool GetInteractionOccured () const
 
- Public Member Functions inherited from G4VBiasingOperation
 G4VBiasingOperation (const G4String &name)
 
virtual ~G4VBiasingOperation ()=default
 
virtual void AlongMoveBy (const G4BiasingProcessInterface *, const G4Step *, G4double)
 
const G4StringGetName () const
 
std::size_t GetUniqueID () const
 

Detailed Description

Definition at line 56 of file G4BOptnForceCommonTruncatedExp.hh.

Constructor & Destructor Documentation

◆ G4BOptnForceCommonTruncatedExp()

G4BOptnForceCommonTruncatedExp::G4BOptnForceCommonTruncatedExp ( const G4String & name)

Definition at line 36 of file G4BOptnForceCommonTruncatedExp.cc.

39{
40 fCommonTruncatedExpLaw = new G4ILawCommonTruncatedExp("ExpLawForOperation"+name);
41 fForceFreeFlightLaw = new G4ILawForceFreeFlight ("FFFLawForOperation"+name);
42}
G4VBiasingOperation(const G4String &name)

◆ ~G4BOptnForceCommonTruncatedExp()

G4BOptnForceCommonTruncatedExp::~G4BOptnForceCommonTruncatedExp ( )
virtual

Definition at line 44 of file G4BOptnForceCommonTruncatedExp.cc.

45{
46 delete fCommonTruncatedExpLaw;
47 delete fForceFreeFlightLaw;
48}

Member Function Documentation

◆ AddCrossSection()

void G4BOptnForceCommonTruncatedExp::AddCrossSection ( const G4VProcess * process,
G4double crossSection )

Definition at line 114 of file G4BOptnForceCommonTruncatedExp.cc.

116{
117 fTotalCrossSection += crossSection;
118 fCrossSections[process] = crossSection;
119 fNumberOfSharing = fCrossSections.size();
120}

◆ ApplyFinalStateBiasing()

G4VParticleChange * G4BOptnForceCommonTruncatedExp::ApplyFinalStateBiasing ( const G4BiasingProcessInterface * callingProcess,
const G4Track * track,
const G4Step * step,
G4bool & forceFinalState )
virtual

Implements G4VBiasingOperation.

Definition at line 72 of file G4BOptnForceCommonTruncatedExp.cc.

77{
78 if ( callingProcess->GetWrappedProcess() != fProcessToApply )
79 {
80 forceFinalState = true;
81 fDummyParticleChange.Initialize( *track );
82 return &fDummyParticleChange;
83 }
84 if ( fInteractionOccured )
85 {
86 forceFinalState = true;
87 fDummyParticleChange.Initialize( *track );
88 return &fDummyParticleChange;
89 }
90
91 // -- checks if process won the GPIL race:
92 G4double processGPIL = callingProcess->GetPostStepGPIL()
93 < callingProcess->GetAlongStepGPIL()
94 ? callingProcess->GetPostStepGPIL()
95 : callingProcess->GetAlongStepGPIL();
96 if ( processGPIL <= step->GetStepLength() )
97 {
98 // -- if process won, wrapped process produces the final state.
99 // -- In this case, the weight for occurrence biasing is applied
100 // -- by the callingProcess, at exit of present method. This is
101 // -- selected by "forceFinalState = false":
102 forceFinalState = false;
103 fInteractionOccured = true;
104 return callingProcess->GetWrappedProcess()->PostStepDoIt( *track, *step );
105 }
106 else
107 {
108 forceFinalState = true;
109 fDummyParticleChange.Initialize( *track );
110 return &fDummyParticleChange;
111 }
112}
double G4double
Definition G4Types.hh:83
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData)=0

◆ ChooseProcessToApply()

void G4BOptnForceCommonTruncatedExp::ChooseProcessToApply ( )

Definition at line 162 of file G4BOptnForceCommonTruncatedExp.cc.

163{
164 G4double sigmaRand = G4UniformRand() * fTotalCrossSection;
165 G4double sigmaSelect = 0.0;
166 for ( auto it = fCrossSections.cbegin(); it != fCrossSections.cend(); ++it)
167 {
168 sigmaSelect += (*it).second;
169 if ( sigmaRand <= sigmaSelect )
170 {
171 fProcessToApply = (*it).first;
172 break;
173 }
174 }
175}
#define G4UniformRand()
Definition Randomize.hh:52

Referenced by Sample().

◆ DistanceToApplyOperation()

virtual G4double G4BOptnForceCommonTruncatedExp::DistanceToApplyOperation ( const G4Track * ,
G4double ,
G4ForceCondition *  )
inlinevirtual

Implements G4VBiasingOperation.

Definition at line 78 of file G4BOptnForceCommonTruncatedExp.hh.

79 { return DBL_MAX; }
#define DBL_MAX
Definition templates.hh:62

◆ GenerateBiasingFinalState()

virtual G4VParticleChange * G4BOptnForceCommonTruncatedExp::GenerateBiasingFinalState ( const G4Track * ,
const G4Step *  )
inlinevirtual

Implements G4VBiasingOperation.

Definition at line 81 of file G4BOptnForceCommonTruncatedExp.hh.

81{ return nullptr; }

◆ GetCommonTruncatedExpLaw()

G4ILawCommonTruncatedExp * G4BOptnForceCommonTruncatedExp::GetCommonTruncatedExpLaw ( )
inline

Definition at line 86 of file G4BOptnForceCommonTruncatedExp.hh.

87 {
88 return fCommonTruncatedExpLaw;
89 }

◆ GetForceFreeFlightLaw()

G4ILawForceFreeFlight * G4BOptnForceCommonTruncatedExp::GetForceFreeFlightLaw ( )
inline

Definition at line 90 of file G4BOptnForceCommonTruncatedExp.hh.

91 {
92 return fForceFreeFlightLaw;
93 }

◆ GetInitialMomentum()

const G4ThreeVector & G4BOptnForceCommonTruncatedExp::GetInitialMomentum ( ) const
inline

Definition at line 98 of file G4BOptnForceCommonTruncatedExp.hh.

98{ return fInitialMomentum; }

◆ GetInteractionOccured()

G4bool G4BOptnForceCommonTruncatedExp::GetInteractionOccured ( ) const
inline

Definition at line 105 of file G4BOptnForceCommonTruncatedExp.hh.

105{ return fInteractionOccured; }

◆ GetMaximumDistance()

G4double G4BOptnForceCommonTruncatedExp::GetMaximumDistance ( ) const
inline

Definition at line 99 of file G4BOptnForceCommonTruncatedExp.hh.

99{ return fMaximumDistance; }

◆ GetNumberOfSharing()

std::size_t G4BOptnForceCommonTruncatedExp::GetNumberOfSharing ( ) const
inline

Definition at line 103 of file G4BOptnForceCommonTruncatedExp.hh.

103{ return fNumberOfSharing; }

◆ GetProcessToApply()

const G4VProcess * G4BOptnForceCommonTruncatedExp::GetProcessToApply ( ) const
inline

Definition at line 101 of file G4BOptnForceCommonTruncatedExp.hh.

101{ return fProcessToApply; }

◆ Initialize()

void G4BOptnForceCommonTruncatedExp::Initialize ( const G4Track * track)

Definition at line 122 of file G4BOptnForceCommonTruncatedExp.cc.

123{
124 fCrossSections.clear();
125 fTotalCrossSection = 0.0;
126 fNumberOfSharing = 0;
127 fProcessToApply = 0;
128 fInteractionOccured = false;
129 fInitialMomentum = track->GetMomentum();
130
131 G4VSolid* currentSolid = track->GetVolume()->GetLogicalVolume()->GetSolid();
133 GetNavigatorForTracking()->
134 GetGlobalToLocalTransform()).TransformPoint(track->GetPosition());
136 GetNavigatorForTracking()->
137 GetGlobalToLocalTransform()).TransformAxis(track->GetMomentumDirection());
138 fMaximumDistance = currentSolid->DistanceToOut(localPosition, localDirection);
139 if ( fMaximumDistance <= DBL_MIN ) { fMaximumDistance = 0.0; }
140 fCommonTruncatedExpLaw->SetMaximumDistance( fMaximumDistance );
141}
CLHEP::Hep3Vector G4ThreeVector
G4VSolid * GetSolid() const
G4VPhysicalVolume * GetVolume() const
const G4ThreeVector & GetPosition() const
G4ThreeVector GetMomentum() const
const G4ThreeVector & GetMomentumDirection() const
static G4TransportationManager * GetTransportationManager()
G4LogicalVolume * GetLogicalVolume() const
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const =0
#define DBL_MIN
Definition templates.hh:54

◆ ProposeAlongStepLimit()

virtual G4double G4BOptnForceCommonTruncatedExp::ProposeAlongStepLimit ( const G4BiasingProcessInterface * )
inlinevirtual

Reimplemented from G4VBiasingOperation.

Definition at line 70 of file G4BOptnForceCommonTruncatedExp.hh.

70{ return DBL_MAX; }

◆ ProposeGPILSelection()

G4GPILSelection G4BOptnForceCommonTruncatedExp::ProposeGPILSelection ( const G4GPILSelection processSelection)
virtual

Reimplemented from G4VBiasingOperation.

Definition at line 66 of file G4BOptnForceCommonTruncatedExp.cc.

68{
70}
@ NotCandidateForSelection

◆ ProvideOccurenceBiasingInteractionLaw()

const G4VBiasingInteractionLaw * G4BOptnForceCommonTruncatedExp::ProvideOccurenceBiasingInteractionLaw ( const G4BiasingProcessInterface * callingProcess,
G4ForceCondition & proposeForceCondition )
virtual

Implements G4VBiasingOperation.

Definition at line 50 of file G4BOptnForceCommonTruncatedExp.cc.

53{
54 if ( callingProcess->GetWrappedProcess() == fProcessToApply )
55 {
56 proposeForceCondition = Forced;
57 return fCommonTruncatedExpLaw;
58 }
59 else
60 {
61 proposeForceCondition = Forced;
62 return fForceFreeFlightLaw;
63 }
64}
@ Forced

◆ Sample()

void G4BOptnForceCommonTruncatedExp::Sample ( )

Definition at line 154 of file G4BOptnForceCommonTruncatedExp.cc.

155{
156 fCommonTruncatedExpLaw->SetForceCrossSection( fTotalCrossSection );
157 fCommonTruncatedExpLaw->Sample();
159 fCommonTruncatedExpLaw->SetSelectedProcessXSfraction(fCrossSections[fProcessToApply] / fTotalCrossSection);
160}

◆ SetInteractionOccured()

void G4BOptnForceCommonTruncatedExp::SetInteractionOccured ( G4bool b)
inline

Definition at line 104 of file G4BOptnForceCommonTruncatedExp.hh.

104{ fInteractionOccured = b; }

◆ UpdateForStep()

void G4BOptnForceCommonTruncatedExp::UpdateForStep ( const G4Step * step)

Definition at line 143 of file G4BOptnForceCommonTruncatedExp.cc.

144{
145 fCrossSections.clear();
146 fTotalCrossSection = 0.0;
147 fNumberOfSharing = 0;
148 fProcessToApply = 0;
149
150 fCommonTruncatedExpLaw->UpdateForStep( step->GetStepLength() );
151 fMaximumDistance = fCommonTruncatedExpLaw->GetMaximumDistance();
152}
G4double GetStepLength() const

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