Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4BOptnForceFreeFlight Class Reference

#include <G4BOptnForceFreeFlight.hh>

Inheritance diagram for G4BOptnForceFreeFlight:
Collaboration diagram for G4BOptnForceFreeFlight:

Public Member Functions

 G4BOptnForceFreeFlight (G4String name)
 
virtual ~G4BOptnForceFreeFlight ()
 
virtual const
G4VBiasingInteractionLaw
ProvideOccurenceBiasingInteractionLaw (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 (G4String name)
 
virtual ~G4VBiasingOperation ()
 
virtual G4double ProposeAlongStepLimit (const G4BiasingProcessInterface *)
 
virtual G4GPILSelection ProposeGPILSelection (const G4GPILSelection wrappedProcessSelection)
 
const G4StringGetName () const
 
std::size_t GetUniqueID () const
 

Detailed Description

Definition at line 56 of file G4BOptnForceFreeFlight.hh.

Constructor & Destructor Documentation

G4BOptnForceFreeFlight::G4BOptnForceFreeFlight ( G4String  name)

Definition at line 33 of file G4BOptnForceFreeFlight.cc.

34  : G4VBiasingOperation ( name ),
35  fCumulatedWeightChange ( -1.0 ),
36  fInitialTrackWeight ( -1.0 ),
37  fOperationComplete ( true )
38 {
39  fForceFreeFlightInteractionLaw = new G4ILawForceFreeFlight("LawForOperation"+name);
40 }
G4VBiasingOperation(G4String name)
G4BOptnForceFreeFlight::~G4BOptnForceFreeFlight ( )
virtual

Definition at line 42 of file G4BOptnForceFreeFlight.cc.

43 {
44  if ( fForceFreeFlightInteractionLaw ) delete fForceFreeFlightInteractionLaw;
45 }

Member Function Documentation

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 }
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 
66  fParticleChange.Initialize( *track );
67  forceFinalState = true;
68  if ( step->GetPostStepPoint()->GetStepStatus() == fGeomBoundary )
69  {
70  // -- Sanity checks:
71  if ( fInitialTrackWeight <= DBL_MIN )
72  {
74  ed << " Initial track weight is null ! " << G4endl;
75  G4Exception(" G4BOptnForceFreeFlight::ApplyFinalStateBiasing(...)",
76  "BIAS.GEN.05",
78  ed);
79  }
80  if ( fCumulatedWeightChange <= DBL_MIN )
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 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4StepStatus GetStepStatus() const
void ProposeWeight(G4double finalWeight)
G4bool GetIsFirstPostStepDoItInterface(G4bool physOnly=true) const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
virtual void Initialize(const G4Track &)
G4StepPoint * GetPostStepPoint() const
#define DBL_MIN
Definition: templates.hh:75
G4double GetWeight() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

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

Implements G4VBiasingOperation.

Definition at line 72 of file G4BOptnForceFreeFlight.hh.

74  {return DBL_MAX;}
#define DBL_MAX
Definition: templates.hh:83
virtual G4VParticleChange* G4BOptnForceFreeFlight::GenerateBiasingFinalState ( const G4Track ,
const G4Step  
)
inlinevirtual

Implements G4VBiasingOperation.

Definition at line 75 of file G4BOptnForceFreeFlight.hh.

76  {return 0;}
G4ILawForceFreeFlight* G4BOptnForceFreeFlight::GetForceFreeFlightLaw ( )
inline

Definition at line 83 of file G4BOptnForceFreeFlight.hh.

83  {
84  return fForceFreeFlightInteractionLaw;
85  }
G4bool G4BOptnForceFreeFlight::OperationComplete ( ) const
inline

Definition at line 88 of file G4BOptnForceFreeFlight.hh.

88 { return fOperationComplete; }
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;
51  return fForceFreeFlightInteractionLaw;
52 }
void G4BOptnForceFreeFlight::ResetInitialTrackWeight ( G4double  w)
inline

Definition at line 87 of file G4BOptnForceFreeFlight.hh.

87 {fInitialTrackWeight = w; fCumulatedWeightChange = 1.0;}

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