Geant4  10.02.p02
G4BOptnForceFreeFlight.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
27 #include "G4ILawForceFreeFlight.hh"
29 #include "G4Step.hh"
30 
31 
32 
34  : G4VBiasingOperation ( name ),
35  fCumulatedWeightChange ( -1.0 ),
36  fInitialTrackWeight ( -1.0 ),
37  fOperationComplete ( true )
38 {
39  fForceFreeFlightInteractionLaw = new G4ILawForceFreeFlight("LawForOperation"+name);
40 }
41 
43 {
45 }
46 
48 {
49  fOperationComplete = false;
50  proposeForceCondition = Forced;
52 }
53 
55 {
56  // -- force free flight always deny process to apply its doit.
57  // -- if reaching boundary, track is restored with non-zero weight
59  {
61  ed << " Initial track weight is null ! " << G4endl;
62  G4Exception(" G4BOptnForceFreeFlight::DenyProcessPostStepDoIt(...)",
63  "BIAS.GEN.05",
65  ed);
66  }
68  {
70  ed << " Cumulated weight is null ! " << G4endl;
71  G4Exception(" G4BOptnForceFreeFlight::DenyProcessPostStepDoIt(...)",
72  "BIAS.GEN.06",
74  ed);
75  }
76  if ( step->GetPostStepPoint()->GetStepStatus() == fGeomBoundary )
77  {
78  if ( proposedWeight <= DBL_MIN ) proposedWeight = fCumulatedWeightChange * fInitialTrackWeight;
79  else proposedWeight *= fCumulatedWeightChange;
80  }
81 
82  return true;
83 }
84 
86  const G4Track* track,
87  const G4Step* step,
88  G4bool& forceFinalState)
89 {
90  // -- If the track is reaching the volume boundary, its free flight ends. In this case, its zero
91  // -- weight is brought back to non-zero value: its initial weight is restored by the first
92  // -- ApplyFinalStateBiasing operation called, and the weight for force free flight is applied
93  // -- is applied by each operation.
94  // -- If the track is not reaching the volume boundary, it zero weight flight continues.
95 
96  fParticleChange.Initialize( *track );
97  forceFinalState = true;
98  if ( step->GetPostStepPoint()->GetStepStatus() == fGeomBoundary )
99  {
100  // -- Sanity checks:
101  if ( fInitialTrackWeight <= DBL_MIN )
102  {
104  ed << " Initial track weight is null ! " << G4endl;
105  G4Exception(" G4BOptnForceFreeFlight::ApplyFinalStateBiasing(...)",
106  "BIAS.GEN.05",
107  JustWarning,
108  ed);
109  }
111  {
113  ed << " Cumulated weight is null ! " << G4endl;
114  G4Exception(" G4BOptnForceFreeFlight::ApplyFinalStateBiasing(...)",
115  "BIAS.GEN.06",
116  JustWarning,
117  ed);
118  }
119 
120  G4double proposedWeight = track->GetWeight();
121  if ( callingProcess->GetIsFirstPostStepDoItInterface() ) proposedWeight = fCumulatedWeightChange * fInitialTrackWeight;
122  else proposedWeight *= fCumulatedWeightChange;
123  fParticleChange.ProposeWeight(proposedWeight);
124  fOperationComplete = true;
125  }
126 
127  return &fParticleChange;
128 }
129 
130 
132 {
133  fCumulatedWeightChange *= weightChange;
134 }
virtual const G4VBiasingInteractionLaw * ProvideOccurenceBiasingInteractionLaw(const G4BiasingProcessInterface *, G4ForceCondition &)
virtual G4VParticleChange * ApplyFinalStateBiasing(const G4BiasingProcessInterface *, const G4Track *, const G4Step *, G4bool &)
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4ILawForceFreeFlight * fForceFreeFlightInteractionLaw
G4String name
Definition: TRTMaterials.hh:40
G4StepStatus GetStepStatus() const
void ProposeWeight(G4double finalWeight)
bool G4bool
Definition: G4Types.hh:79
G4bool GetIsFirstPostStepDoItInterface(G4bool physOnly=true) const
Definition: G4Step.hh:76
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
virtual void Initialize(const G4Track &)
virtual G4bool DenyProcessPostStepDoIt(const G4BiasingProcessInterface *, const G4Track *, const G4Step *, G4double &)
G4StepPoint * GetPostStepPoint() const
#define DBL_MIN
Definition: templates.hh:75
G4double GetWeight() const
#define G4endl
Definition: G4ios.hh:61
virtual void AlongMoveBy(const G4BiasingProcessInterface *, const G4Step *, G4double)
double G4double
Definition: G4Types.hh:76
G4ForceCondition