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

#include <G4BOptnForceFreeFlight.hh>

+ Inheritance diagram for G4BOptnForceFreeFlight:

Public Member Functions

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

Detailed Description

Definition at line 51 of file G4BOptnForceFreeFlight.hh.

Constructor & Destructor Documentation

◆ G4BOptnForceFreeFlight()

G4BOptnForceFreeFlight::G4BOptnForceFreeFlight ( const G4String & name)

Definition at line 33 of file G4BOptnForceFreeFlight.cc.

34 : G4VBiasingOperation( name )
35{
36 fForceFreeFlightInteractionLaw = new G4ILawForceFreeFlight("LawForOperation"+name);
37}
G4VBiasingOperation(const G4String &name)

◆ ~G4BOptnForceFreeFlight()

G4BOptnForceFreeFlight::~G4BOptnForceFreeFlight ( )
virtual

Definition at line 39 of file G4BOptnForceFreeFlight.cc.

40{
41 delete fForceFreeFlightInteractionLaw;
42}

Member Function Documentation

◆ AlongMoveBy()

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

Reimplemented from G4VBiasingOperation.

Definition at line 95 of file G4BOptnForceFreeFlight.cc.

98{
99 fCumulatedWeightChange *= weightChange;
100}

◆ ApplyFinalStateBiasing()

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

Implements G4VBiasingOperation.

Definition at line 53 of file G4BOptnForceFreeFlight.cc.

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

◆ DistanceToApplyOperation()

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

Implements G4VBiasingOperation.

Definition at line 73 of file G4BOptnForceFreeFlight.hh.

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

◆ GenerateBiasingFinalState()

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

Implements G4VBiasingOperation.

Definition at line 75 of file G4BOptnForceFreeFlight.hh.

76 { return nullptr; }

◆ GetForceFreeFlightLaw()

G4ILawForceFreeFlight * G4BOptnForceFreeFlight::GetForceFreeFlightLaw ( )
inline

Definition at line 81 of file G4BOptnForceFreeFlight.hh.

82 {
83 return fForceFreeFlightInteractionLaw;
84 }

◆ OperationComplete()

G4bool G4BOptnForceFreeFlight::OperationComplete ( ) const
inline

Definition at line 91 of file G4BOptnForceFreeFlight.hh.

92 {
93 return fOperationComplete;
94 }

◆ ProvideOccurenceBiasingInteractionLaw()

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

Implements G4VBiasingOperation.

Definition at line 44 of file G4BOptnForceFreeFlight.cc.

47{
48 fOperationComplete = false;
49 proposeForceCondition = Forced;
50 return fForceFreeFlightInteractionLaw;
51}
@ Forced

◆ ResetInitialTrackWeight()

void G4BOptnForceFreeFlight::ResetInitialTrackWeight ( G4double w)
inline

Definition at line 86 of file G4BOptnForceFreeFlight.hh.

87 {
88 fInitialTrackWeight = w;
89 fCumulatedWeightChange = 1.0;
90 }

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