Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4BOptnForceFreeFlight Class Reference

#include <G4BOptnForceFreeFlight.hh>

Inheritance diagram for G4BOptnForceFreeFlight:
G4VBiasingOperation

Public Member Functions

 G4BOptnForceFreeFlight (G4String name)
 
virtual ~G4BOptnForceFreeFlight ()
 
virtual const
G4VBiasingInteractionLaw
ProvideOccurenceBiasingInteractionLaw (const G4BiasingProcessInterface *)
 
virtual G4ForceCondition ProposeForceCondition (const G4ForceCondition)
 
virtual G4bool DenyProcessPostStepDoIt (const G4BiasingProcessInterface *, const G4Track *, const G4Step *, G4double &)
 
virtual void AlongMoveBy (const G4BiasingProcessInterface *, const G4Step *, G4double)
 
virtual G4VParticleChangeApplyFinalStateBiasing (const G4BiasingProcessInterface *, const G4Track *, const G4Step *)
 
virtual G4double DistanceToApplyOperation (const G4Track *, G4double, G4ForceCondition *)
 
virtual G4VParticleChangeGenerateBiasingFinalState (const G4Track *, const G4Step *)
 
G4ILawForceFreeFlightGetForceFreeFlightLaw ()
 
void ResetInitialTrackWeight (G4double w)
 
- 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 54 of file G4BOptnForceFreeFlight.hh.

Constructor & Destructor Documentation

G4BOptnForceFreeFlight::G4BOptnForceFreeFlight ( G4String  name)

Definition at line 32 of file G4BOptnForceFreeFlight.cc.

33  : G4VBiasingOperation(name)
34 {
35  fForceFreeFlightInteractionLaw = new G4ILawForceFreeFlight("LawForOperation"+name);
36 }
G4VBiasingOperation(G4String name)
G4BOptnForceFreeFlight::~G4BOptnForceFreeFlight ( )
virtual

Definition at line 38 of file G4BOptnForceFreeFlight.cc.

39 {}

Member Function Documentation

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

Reimplemented from G4VBiasingOperation.

Definition at line 77 of file G4BOptnForceFreeFlight.cc.

78 {
79  fCumulatedWeightChange *= weightChange;
80 }
virtual G4VParticleChange* G4BOptnForceFreeFlight::ApplyFinalStateBiasing ( const G4BiasingProcessInterface ,
const G4Track ,
const G4Step  
)
inlinevirtual

Implements G4VBiasingOperation.

Definition at line 71 of file G4BOptnForceFreeFlight.hh.

73  {return 0;}
G4bool G4BOptnForceFreeFlight::DenyProcessPostStepDoIt ( const G4BiasingProcessInterface ,
const G4Track ,
const G4Step step,
G4double proposedWeight 
)
virtual

Reimplemented from G4VBiasingOperation.

Definition at line 46 of file G4BOptnForceFreeFlight.cc.

References DBL_MIN, fGeomBoundary, G4endl, G4Exception(), G4Step::GetPostStepPoint(), G4StepPoint::GetStepStatus(), and JustWarning.

47 {
48  // -- force free flight always deny process to apply its doit.
49  // -- if reaching boundary, track is restored with non-zero weight
50  if ( fInitialTrackWeight <= DBL_MIN )
51  {
53  ed << " Initial track weight is null ! " << G4endl;
54  G4Exception(" G4BOptnForceFreeFlight::DenyProcessPostStepDoIt(...)",
55  "BIAS.GEN.05",
57  ed);
58  }
59  if ( fCumulatedWeightChange <= DBL_MIN )
60  {
62  ed << " Cumulated weight is null ! " << G4endl;
63  G4Exception(" G4BOptnForceFreeFlight::DenyProcessPostStepDoIt(...)",
64  "BIAS.GEN.06",
66  ed);
67  }
68  if ( step->GetPostStepPoint()->GetStepStatus() == fGeomBoundary )
69  {
70  if ( proposedWeight <= DBL_MIN ) proposedWeight = fCumulatedWeightChange * fInitialTrackWeight;
71  else proposedWeight *= fCumulatedWeightChange;
72  }
73 
74  return true;
75 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4StepStatus GetStepStatus() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4StepPoint * GetPostStepPoint() const
#define DBL_MIN
Definition: templates.hh:75
#define G4endl
Definition: G4ios.hh:61
virtual G4double G4BOptnForceFreeFlight::DistanceToApplyOperation ( const G4Track ,
G4double  ,
G4ForceCondition  
)
inlinevirtual

Implements G4VBiasingOperation.

Definition at line 74 of file G4BOptnForceFreeFlight.hh.

References DBL_MAX.

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

Implements G4VBiasingOperation.

Definition at line 77 of file G4BOptnForceFreeFlight.hh.

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

Definition at line 85 of file G4BOptnForceFreeFlight.hh.

85  {
86  return fForceFreeFlightInteractionLaw;
87  }
virtual G4ForceCondition G4BOptnForceFreeFlight::ProposeForceCondition ( const G4ForceCondition  )
inlinevirtual

Reimplemented from G4VBiasingOperation.

Definition at line 66 of file G4BOptnForceFreeFlight.hh.

References Forced.

66 {return Forced;}
const G4VBiasingInteractionLaw * G4BOptnForceFreeFlight::ProvideOccurenceBiasingInteractionLaw ( const G4BiasingProcessInterface )
virtual

Implements G4VBiasingOperation.

Definition at line 41 of file G4BOptnForceFreeFlight.cc.

42 {
43  return fForceFreeFlightInteractionLaw;
44 }
void G4BOptnForceFreeFlight::ResetInitialTrackWeight ( G4double  w)
inline

Definition at line 89 of file G4BOptnForceFreeFlight.hh.

89 {fInitialTrackWeight = w; fCumulatedWeightChange = 1.0;}

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