Geant4 11.3.0
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//
26// G4BOptnForceFreeFlight
27// --------------------------------------------------------------------
31#include "G4Step.hh"
32
34 : G4VBiasingOperation( name )
35{
36 fForceFreeFlightInteractionLaw = new G4ILawForceFreeFlight("LawForOperation"+name);
37}
38
40{
41 delete fForceFreeFlightInteractionLaw;
42}
43
46 G4ForceCondition& proposeForceCondition )
47{
48 fOperationComplete = false;
49 proposeForceCondition = Forced;
50 return fForceFreeFlightInteractionLaw;
51}
52
55 const G4Track* track, const G4Step* step,
56 G4bool& forceFinalState)
57{
58 // -- If the track is reaching the volume boundary, its free flight ends. In this case, its zero
59 // -- weight is brought back to non-zero value: its initial weight is restored by the first
60 // -- ApplyFinalStateBiasing operation called, and the weight for force free flight is applied
61 // -- is applied by each operation.
62 // -- If the track is not reaching the volume boundary, it zero weight flight continues.
63
64 fParticleChange.Initialize( *track );
65 forceFinalState = true;
67 {
68 // -- Sanity checks:
69 if ( fInitialTrackWeight <= DBL_MIN )
70 {
72 ed << " Initial track weight is null ! " << G4endl;
73 G4Exception(" G4BOptnForceFreeFlight::ApplyFinalStateBiasing(...)",
74 "BIAS.GEN.05", JustWarning, ed);
75 }
76 if ( fCumulatedWeightChange <= DBL_MIN )
77 {
79 ed << " Cumulated weight is null ! " << G4endl;
80 G4Exception(" G4BOptnForceFreeFlight::ApplyFinalStateBiasing(...)",
81 "BIAS.GEN.06", JustWarning, ed);
82 }
83 G4double proposedWeight = track->GetWeight();
84 if ( callingProcess->GetIsFirstPostStepDoItInterface() )
85 proposedWeight = fCumulatedWeightChange * fInitialTrackWeight;
86 else
87 proposedWeight *= fCumulatedWeightChange;
88 fParticleChange.ProposeWeight(proposedWeight);
89 fOperationComplete = true;
90 }
91
92 return &fParticleChange;
93}
94
97 const G4Step*, G4double weightChange )
98{
99 fCumulatedWeightChange *= weightChange;
100}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
G4ForceCondition
@ Forced
@ fGeomBoundary
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
#define G4endl
Definition G4ios.hh:67
G4BOptnForceFreeFlight(const G4String &name)
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
G4StepStatus GetStepStatus() const
G4StepPoint * GetPostStepPoint() const
G4double GetWeight() const
G4VBiasingOperation(const G4String &name)
virtual void Initialize(const G4Track &)
#define DBL_MIN
Definition templates.hh:54