Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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//
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{
44 if ( fForceFreeFlightInteractionLaw ) delete fForceFreeFlightInteractionLaw;
45}
46
48{
49 fOperationComplete = false;
50 proposeForceCondition = Forced;
51 return fForceFreeFlightInteractionLaw;
52}
53
54
56 const G4Track* track,
57 const G4Step* step,
58 G4bool& forceFinalState)
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;
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}
99
100
102{
103 fCumulatedWeightChange *= weightChange;
104}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4ForceCondition
@ Forced
@ fGeomBoundary
Definition: G4StepStatus.hh:43
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
#define G4endl
Definition: G4ios.hh:57
virtual const G4VBiasingInteractionLaw * ProvideOccurenceBiasingInteractionLaw(const G4BiasingProcessInterface *, G4ForceCondition &)
virtual void AlongMoveBy(const G4BiasingProcessInterface *, const G4Step *, G4double)
virtual G4VParticleChange * ApplyFinalStateBiasing(const G4BiasingProcessInterface *, const G4Track *, const G4Step *, G4bool &)
G4bool GetIsFirstPostStepDoItInterface(G4bool physOnly=true) const
void Initialize(const G4Track &) override
G4StepStatus GetStepStatus() const
Definition: G4Step.hh:62
G4StepPoint * GetPostStepPoint() const
G4double GetWeight() const
void ProposeWeight(G4double finalWeight)
#define DBL_MIN
Definition: templates.hh:54