Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4BOptnForceFreeFlight.hh
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//
28// Class Description:
29//
30// A G4VBiasingOperation physics-based biasing operation.
31// If forces the physics process to not act on the track.
32// In this implementation (meant for the ForceCollision
33// operator) the free flight is done under zero weight for
34// the track, and the action is meant to accumulate the weight
35// change for making this uninteracting flight,
36// cumulatedWeightChange.
37// When the track reaches the current volume boundary, its
38// weight is restored with value: initialWeight * cumulatedWeightChange
39//
40// Author: Marc Verderi, November 2013.
41// --------------------------------------------------------------------
42#ifndef G4BOptnForceFreeFlight_hh
43#define G4BOptnForceFreeFlight_hh 1
44
46#include "G4ForceCondition.hh"
47#include "G4ParticleChange.hh" // Should add a dedicated "weight change only"
48 // particle change
50
52{
53 public:
54
55 // -- Constructor :
57 // -- destructor:
59
60 // -- Methods from G4VBiasingOperation interface:
61 // -------------------------------------------
62 // -- Used:
63 virtual const G4VBiasingInteractionLaw*
66 virtual void AlongMoveBy( const G4BiasingProcessInterface*,
67 const G4Step*, G4double );
68 virtual G4VParticleChange*
70 const G4Track*, const G4Step*, G4bool& );
71
72 // -- Unused:
76 const G4Step* ) { return nullptr; }
77
78 // -- Additional methods, specific to this class:
79 // ----------------------------------------------
80 // -- return concrete type of interaction law:
82 {
83 return fForceFreeFlightInteractionLaw;
84 }
85 // -- initialization for weight:
87 {
88 fInitialTrackWeight = w;
89 fCumulatedWeightChange = 1.0;
90 }
92 {
93 return fOperationComplete;
94 }
95
96 private:
97
98 G4ILawForceFreeFlight* fForceFreeFlightInteractionLaw = nullptr;
99 G4double fCumulatedWeightChange = -1.0,
100 fInitialTrackWeight = -1.0;
101 G4ParticleChange fParticleChange;
102 G4bool fOperationComplete = true;
103};
104
105#endif
G4ForceCondition
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
G4BOptnForceFreeFlight(const G4String &name)
G4ILawForceFreeFlight * GetForceFreeFlightLaw()
virtual const G4VBiasingInteractionLaw * ProvideOccurenceBiasingInteractionLaw(const G4BiasingProcessInterface *, G4ForceCondition &)
virtual G4VParticleChange * GenerateBiasingFinalState(const G4Track *, const G4Step *)
virtual G4double DistanceToApplyOperation(const G4Track *, G4double, G4ForceCondition *)
virtual void AlongMoveBy(const G4BiasingProcessInterface *, const G4Step *, G4double)
void ResetInitialTrackWeight(G4double w)
virtual G4VParticleChange * ApplyFinalStateBiasing(const G4BiasingProcessInterface *, const G4Track *, const G4Step *, G4bool &)
G4VBiasingOperation(const G4String &name)
#define DBL_MAX
Definition templates.hh:62