Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4BOptnForceCommonTruncatedExp Class Reference

#include <G4BOptnForceCommonTruncatedExp.hh>

Inheritance diagram for G4BOptnForceCommonTruncatedExp:
Collaboration diagram for G4BOptnForceCommonTruncatedExp:

Public Member Functions

 G4BOptnForceCommonTruncatedExp (G4String name)
 
virtual ~G4BOptnForceCommonTruncatedExp ()
 
virtual const
G4VBiasingInteractionLaw
ProvideOccurenceBiasingInteractionLaw (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)
 
size_t GetNumberOfSharing () const
 
void SetInteractionOccured (G4bool b)
 
G4bool GetInteractionOccured () const
 
- Public Member Functions inherited from G4VBiasingOperation
 G4VBiasingOperation (G4String name)
 
virtual ~G4VBiasingOperation ()
 
virtual void AlongMoveBy (const G4BiasingProcessInterface *, const G4Step *, G4double)
 
const G4StringGetName () const
 
std::size_t GetUniqueID () const
 

Detailed Description

Definition at line 60 of file G4BOptnForceCommonTruncatedExp.hh.

Constructor & Destructor Documentation

G4BOptnForceCommonTruncatedExp::G4BOptnForceCommonTruncatedExp ( G4String  name)

Definition at line 34 of file G4BOptnForceCommonTruncatedExp.cc.

35  : G4VBiasingOperation(name),
36  fNumberOfSharing(0),
37  fProcessToApply(nullptr),
38  fInteractionOccured(false),
39  fMaximumDistance(-1.0)
40 {
41  fCommonTruncatedExpLaw = new G4ILawCommonTruncatedExp("ExpLawForOperation"+name);
42  fForceFreeFlightLaw = new G4ILawForceFreeFlight ("FFFLawForOperation"+name);
43 
44  fTotalCrossSection = 0.0;
45 }
G4VBiasingOperation(G4String name)
G4BOptnForceCommonTruncatedExp::~G4BOptnForceCommonTruncatedExp ( )
virtual

Definition at line 47 of file G4BOptnForceCommonTruncatedExp.cc.

48 {
49  if ( fCommonTruncatedExpLaw ) delete fCommonTruncatedExpLaw;
50  if ( fForceFreeFlightLaw ) delete fForceFreeFlightLaw;
51 }

Member Function Documentation

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

Definition at line 116 of file G4BOptnForceCommonTruncatedExp.cc.

117 {
118  fTotalCrossSection += crossSection;
119  fCrossSections[process] = crossSection;
120  fNumberOfSharing = fCrossSections.size();
121 }
G4VParticleChange * G4BOptnForceCommonTruncatedExp::ApplyFinalStateBiasing ( const G4BiasingProcessInterface callingProcess,
const G4Track track,
const G4Step step,
G4bool forceFinalState 
)
virtual

Implements G4VBiasingOperation.

Definition at line 76 of file G4BOptnForceCommonTruncatedExp.cc.

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

Here is the call graph for this function:

void G4BOptnForceCommonTruncatedExp::ChooseProcessToApply ( )

Definition at line 167 of file G4BOptnForceCommonTruncatedExp.cc.

168 {
169  G4double sigmaRand = G4UniformRand() * fTotalCrossSection;
170  G4double sigmaSelect = 0.0;
171  for ( std::map< const G4VProcess*, G4double>::const_iterator it = fCrossSections.begin();
172  it != fCrossSections.end();
173  it++)
174  {
175  sigmaSelect += (*it).second;
176  if ( sigmaRand <= sigmaSelect )
177  {
178  fProcessToApply = (*it).first;
179  break;
180  }
181  }
182 }
#define G4UniformRand()
Definition: Randomize.hh:97
double G4double
Definition: G4Types.hh:76

Here is the caller graph for this function:

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

Implements G4VBiasingOperation.

Definition at line 79 of file G4BOptnForceCommonTruncatedExp.hh.

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

Implements G4VBiasingOperation.

Definition at line 82 of file G4BOptnForceCommonTruncatedExp.hh.

83  {return 0;}
G4ILawCommonTruncatedExp* G4BOptnForceCommonTruncatedExp::GetCommonTruncatedExpLaw ( )
inline

Definition at line 89 of file G4BOptnForceCommonTruncatedExp.hh.

90  {
91  return fCommonTruncatedExpLaw;
92  }
G4ILawForceFreeFlight* G4BOptnForceCommonTruncatedExp::GetForceFreeFlightLaw ( )
inline

Definition at line 93 of file G4BOptnForceCommonTruncatedExp.hh.

94  {
95  return fForceFreeFlightLaw;
96  }
const G4ThreeVector& G4BOptnForceCommonTruncatedExp::GetInitialMomentum ( ) const
inline

Definition at line 101 of file G4BOptnForceCommonTruncatedExp.hh.

101 { return fInitialMomentum; }
G4bool G4BOptnForceCommonTruncatedExp::GetInteractionOccured ( ) const
inline

Definition at line 108 of file G4BOptnForceCommonTruncatedExp.hh.

108 { return fInteractionOccured; }

Here is the caller graph for this function:

G4double G4BOptnForceCommonTruncatedExp::GetMaximumDistance ( ) const
inline

Definition at line 102 of file G4BOptnForceCommonTruncatedExp.hh.

102 { return fMaximumDistance; }
size_t G4BOptnForceCommonTruncatedExp::GetNumberOfSharing ( ) const
inline

Definition at line 106 of file G4BOptnForceCommonTruncatedExp.hh.

106 { return fNumberOfSharing;}
const G4VProcess* G4BOptnForceCommonTruncatedExp::GetProcessToApply ( ) const
inline

Definition at line 104 of file G4BOptnForceCommonTruncatedExp.hh.

104 { return fProcessToApply; }
void G4BOptnForceCommonTruncatedExp::Initialize ( const G4Track track)

Definition at line 124 of file G4BOptnForceCommonTruncatedExp.cc.

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

Here is the call graph for this function:

virtual G4double G4BOptnForceCommonTruncatedExp::ProposeAlongStepLimit ( const G4BiasingProcessInterface )
inlinevirtual

Reimplemented from G4VBiasingOperation.

Definition at line 72 of file G4BOptnForceCommonTruncatedExp.hh.

72 { return DBL_MAX; }
#define DBL_MAX
Definition: templates.hh:83
G4GPILSelection G4BOptnForceCommonTruncatedExp::ProposeGPILSelection ( const G4GPILSelection  processSelection)
virtual

Reimplemented from G4VBiasingOperation.

Definition at line 70 of file G4BOptnForceCommonTruncatedExp.cc.

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

Implements G4VBiasingOperation.

Definition at line 54 of file G4BOptnForceCommonTruncatedExp.cc.

56 {
57  if ( callingProcess->GetWrappedProcess() == fProcessToApply )
58  {
59  proposeForceCondition = Forced;
60  return fCommonTruncatedExpLaw;
61  }
62  else
63  {
64  proposeForceCondition = Forced;
65  return fForceFreeFlightLaw;
66  }
67 }
G4VProcess * GetWrappedProcess() const

Here is the call graph for this function:

void G4BOptnForceCommonTruncatedExp::Sample ( )

Definition at line 158 of file G4BOptnForceCommonTruncatedExp.cc.

159 {
160  fCommonTruncatedExpLaw->SetForceCrossSection( fTotalCrossSection );
161  fCommonTruncatedExpLaw->Sample();
163  fCommonTruncatedExpLaw->SetSelectedProcessXSfraction(fCrossSections[fProcessToApply] / fTotalCrossSection);
164 }
void SetSelectedProcessXSfraction(G4double fXS)

Here is the call graph for this function:

void G4BOptnForceCommonTruncatedExp::SetInteractionOccured ( G4bool  b)
inline

Definition at line 107 of file G4BOptnForceCommonTruncatedExp.hh.

107 { fInteractionOccured = b; }
void G4BOptnForceCommonTruncatedExp::UpdateForStep ( const G4Step step)

Definition at line 146 of file G4BOptnForceCommonTruncatedExp.cc.

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

Here is the call graph for this function:


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