Geant4-11
Public Member Functions | Private Attributes
G4BOptnForceFreeFlight Class Reference

#include <G4BOptnForceFreeFlight.hh>

Inheritance diagram for G4BOptnForceFreeFlight:
G4VBiasingOperation

Public Member Functions

virtual void AlongMoveBy (const G4BiasingProcessInterface *, const G4Step *, G4double)
 
virtual G4VParticleChangeApplyFinalStateBiasing (const G4BiasingProcessInterface *, const G4Track *, const G4Step *, G4bool &)
 
virtual G4double DistanceToApplyOperation (const G4Track *, G4double, G4ForceCondition *)
 
 G4BOptnForceFreeFlight (G4String name)
 
virtual G4VParticleChangeGenerateBiasingFinalState (const G4Track *, const G4Step *)
 
G4ILawForceFreeFlightGetForceFreeFlightLaw ()
 
const G4StringGetName () const
 
std::size_t GetUniqueID () const
 
G4bool OperationComplete () const
 
virtual G4double ProposeAlongStepLimit (const G4BiasingProcessInterface *)
 
virtual G4GPILSelection ProposeGPILSelection (const G4GPILSelection wrappedProcessSelection)
 
virtual const G4VBiasingInteractionLawProvideOccurenceBiasingInteractionLaw (const G4BiasingProcessInterface *, G4ForceCondition &)
 
void ResetInitialTrackWeight (G4double w)
 
virtual ~G4BOptnForceFreeFlight ()
 

Private Attributes

G4double fCumulatedWeightChange
 
G4ILawForceFreeFlightfForceFreeFlightInteractionLaw
 
G4double fInitialTrackWeight
 
const G4String fName
 
G4bool fOperationComplete
 
G4ParticleChange fParticleChange
 
std::size_t fUniqueID
 

Detailed Description

Definition at line 55 of file G4BOptnForceFreeFlight.hh.

Constructor & Destructor Documentation

◆ G4BOptnForceFreeFlight()

G4BOptnForceFreeFlight::G4BOptnForceFreeFlight ( G4String  name)

◆ ~G4BOptnForceFreeFlight()

G4BOptnForceFreeFlight::~G4BOptnForceFreeFlight ( )
virtual

Member Function Documentation

◆ AlongMoveBy()

void G4BOptnForceFreeFlight::AlongMoveBy ( const G4BiasingProcessInterface ,
const G4Step ,
G4double  weightChange 
)
virtual

Reimplemented from G4VBiasingOperation.

Definition at line 101 of file G4BOptnForceFreeFlight.cc.

102{
103 fCumulatedWeightChange *= weightChange;
104}

References fCumulatedWeightChange.

◆ ApplyFinalStateBiasing()

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

Implements G4VBiasingOperation.

Definition at line 55 of file G4BOptnForceFreeFlight.cc.

59{
60 // -- If the track is reaching the volume boundary, its free flight ends. In this case, its zero
61 // -- weight is brought back to non-zero value: its initial weight is restored by the first
62 // -- ApplyFinalStateBiasing operation called, and the weight for force free flight is applied
63 // -- is applied by each operation.
64 // -- If the track is not reaching the volume boundary, it zero weight flight continues.
65
67 forceFinalState = true;
69 {
70 // -- Sanity checks:
72 {
74 ed << " Initial track weight is null ! " << G4endl;
75 G4Exception(" G4BOptnForceFreeFlight::ApplyFinalStateBiasing(...)",
76 "BIAS.GEN.05",
78 ed);
79 }
81 {
83 ed << " Cumulated weight is null ! " << G4endl;
84 G4Exception(" G4BOptnForceFreeFlight::ApplyFinalStateBiasing(...)",
85 "BIAS.GEN.06",
87 ed);
88 }
89
90 G4double proposedWeight = track->GetWeight();
91 if ( callingProcess->GetIsFirstPostStepDoItInterface() ) proposedWeight = fCumulatedWeightChange * fInitialTrackWeight;
92 else proposedWeight *= fCumulatedWeightChange;
93 fParticleChange.ProposeWeight(proposedWeight);
94 fOperationComplete = true;
95 }
96
97 return &fParticleChange;
98}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
@ fGeomBoundary
Definition: G4StepStatus.hh:43
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
G4bool GetIsFirstPostStepDoItInterface(G4bool physOnly=true) const
virtual void Initialize(const G4Track &)
G4StepStatus GetStepStatus() const
G4StepPoint * GetPostStepPoint() const
G4double GetWeight() const
void ProposeWeight(G4double finalWeight)
#define DBL_MIN
Definition: templates.hh:54

References DBL_MIN, fCumulatedWeightChange, fGeomBoundary, fInitialTrackWeight, fOperationComplete, fParticleChange, G4endl, G4Exception(), G4BiasingProcessInterface::GetIsFirstPostStepDoItInterface(), G4Step::GetPostStepPoint(), G4StepPoint::GetStepStatus(), G4Track::GetWeight(), G4ParticleChange::Initialize(), JustWarning, and G4VParticleChange::ProposeWeight().

◆ DistanceToApplyOperation()

virtual G4double G4BOptnForceFreeFlight::DistanceToApplyOperation ( const G4Track ,
G4double  ,
G4ForceCondition  
)
inlinevirtual

Implements G4VBiasingOperation.

Definition at line 71 of file G4BOptnForceFreeFlight.hh.

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

References DBL_MAX.

◆ GenerateBiasingFinalState()

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

Implements G4VBiasingOperation.

Definition at line 74 of file G4BOptnForceFreeFlight.hh.

75 {return 0;}

◆ GetForceFreeFlightLaw()

G4ILawForceFreeFlight * G4BOptnForceFreeFlight::GetForceFreeFlightLaw ( )
inline

Definition at line 82 of file G4BOptnForceFreeFlight.hh.

82 {
84 }

References fForceFreeFlightInteractionLaw.

◆ GetName()

const G4String & G4VBiasingOperation::GetName ( ) const
inlineinherited

Definition at line 198 of file G4VBiasingOperation.hh.

198{return fName;}

References G4VBiasingOperation::fName.

◆ GetUniqueID()

std::size_t G4VBiasingOperation::GetUniqueID ( ) const
inlineinherited

Definition at line 199 of file G4VBiasingOperation.hh.

199{return fUniqueID;}

References G4VBiasingOperation::fUniqueID.

◆ OperationComplete()

G4bool G4BOptnForceFreeFlight::OperationComplete ( ) const
inline

Definition at line 87 of file G4BOptnForceFreeFlight.hh.

87{ return fOperationComplete; }

References fOperationComplete.

◆ ProposeAlongStepLimit()

virtual G4double G4VBiasingOperation::ProposeAlongStepLimit ( const G4BiasingProcessInterface )
inlinevirtualinherited

Reimplemented in G4BOptnForceCommonTruncatedExp.

Definition at line 127 of file G4VBiasingOperation.hh.

127{ return DBL_MAX; }

References DBL_MAX.

Referenced by G4BiasingProcessInterface::AlongStepGetPhysicalInteractionLength().

◆ ProposeGPILSelection()

virtual G4GPILSelection G4VBiasingOperation::ProposeGPILSelection ( const G4GPILSelection  wrappedProcessSelection)
inlinevirtualinherited

Reimplemented in G4BOptnForceCommonTruncatedExp.

Definition at line 131 of file G4VBiasingOperation.hh.

132 {return wrappedProcessSelection;}

Referenced by G4BiasingProcessInterface::AlongStepGetPhysicalInteractionLength().

◆ ProvideOccurenceBiasingInteractionLaw()

const G4VBiasingInteractionLaw * G4BOptnForceFreeFlight::ProvideOccurenceBiasingInteractionLaw ( const G4BiasingProcessInterface ,
G4ForceCondition proposeForceCondition 
)
virtual

Implements G4VBiasingOperation.

Definition at line 47 of file G4BOptnForceFreeFlight.cc.

48{
49 fOperationComplete = false;
50 proposeForceCondition = Forced;
52}
@ Forced

References fForceFreeFlightInteractionLaw, fOperationComplete, and Forced.

◆ ResetInitialTrackWeight()

void G4BOptnForceFreeFlight::ResetInitialTrackWeight ( G4double  w)
inline

Field Documentation

◆ fCumulatedWeightChange

G4double G4BOptnForceFreeFlight::fCumulatedWeightChange
private

◆ fForceFreeFlightInteractionLaw

G4ILawForceFreeFlight* G4BOptnForceFreeFlight::fForceFreeFlightInteractionLaw
private

◆ fInitialTrackWeight

G4double G4BOptnForceFreeFlight::fInitialTrackWeight
private

Definition at line 92 of file G4BOptnForceFreeFlight.hh.

Referenced by ApplyFinalStateBiasing(), and ResetInitialTrackWeight().

◆ fName

const G4String G4VBiasingOperation::fName
privateinherited

Definition at line 203 of file G4VBiasingOperation.hh.

Referenced by G4VBiasingOperation::GetName().

◆ fOperationComplete

G4bool G4BOptnForceFreeFlight::fOperationComplete
private

◆ fParticleChange

G4ParticleChange G4BOptnForceFreeFlight::fParticleChange
private

Definition at line 93 of file G4BOptnForceFreeFlight.hh.

Referenced by ApplyFinalStateBiasing().

◆ fUniqueID

std::size_t G4VBiasingOperation::fUniqueID
privateinherited

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