Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VBiasingOperation.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//
27//
28// --------------------------------------------------------------------
29// Geant4 class header file
30//
31// Class Description:
32//
33// An abstract class to model the behavior of any type of biasing :
34// physics-based biasing (change of physics process behavior) or non-
35// physics-based one, like splitting, killing.
36//
37// o The change of behavior of a physics process can be:
38// - a change of the PostStep interaction probabilty, so-called
39// occurrence biasing
40// - a change in final state production
41// - both, provided above two are uncorrelated.
42// o The change of occurrence is driven by providing a biasing interaction
43// law (G4VBiasingInteractionLaw) that is used in place of the analog
44// exponential law.
45// This change of occurrence is controlled through many handles.
46// o The change in final state production is made through one single
47// method the user is fully responsible of.
48//
49// o Non-physics-based biasing is controlled by two methods : one to
50// specify where this biasing should happen, and one for generating
51// the related final-state.
52//
53// ----------------G4VBiasingOperation ----------------
54//
55// Author: M.Verderi (LLR), November 2013
56// - 05/11/13 : first implementation
57// - 07/11/14 : suppress DenyProcessPostStepDoIt(...) as redondant
58// and special case of ApplyFinalStateBiasing(...)
59// ---------------------------------------------------------------------
60
61
62
63#ifndef G4VBiasingOperation_hh
64#define G4VBiasingOperation_hh 1
65
66#include "globals.hh"
68class G4Track;
69class G4Step;
71class G4VProcess;
73#include "G4ForceCondition.hh"
74#include "G4GPILSelection.hh"
75
77public:
78 // ---------------
79 // -- Constructor:
80 // ---------------
81 //
82 // -- Constructor for biasing operations:
83 // -----------------------------------
84 // -- Operation is given a name.
86
87 // -- destructor:
88 virtual ~G4VBiasingOperation();
89
90public:
91 // -----------------------------
92 // -- Interface to sub-classes :
93 // -----------------------------
94 // --
95 // *************************************
96 // ** Methods for physics-based biasing:
97 // *************************************
98 // --
99 // ---- I. Biasing of the process occurrence:
100 // -----------------------------------------
101 // ---- The biasing of the process occurrence regards the occurrence of the PostStepDoIt
102 // ---- behavior. But the weight is manipulated by both AlongStep methods (weight for
103 // ---- non-interaction) and PostStep methods (weight for interaction). For this
104 // ---- reason, occurrence biasing is handled by both AlongStep and PostStep methods.
105 // ----
106 // ---- If the operation is returned to the G4BiasingProcessInterface process by the
107 // ---- ProposeOccurenceBiasingOperation(...)/GetProposedOccurenceBiasingOperation(...) method
108 // ---- of the biasing operator, all methods below will be called for this operation.
109 // ----
110 // ---- I.1) Methods called in at the PostStepGetPhysicalInteractionLength(...) level :
111 // ----
112 // ------ o Main and mandatory method for biasing of the PostStep process biasing occurrence :
113 // ------ - propose an interaction law to be substituted to the process that is biased
114 // ------ - the operation is told which is the G4BiasingProcessInterface calling it with
115 // ------ callingProcess argument.
116 // ------ - the returned law will have to have been sampled prior to be returned as it will be
117 // ------ asked for its GetSampledInteractionLength() by the callingProcess.
118 // ------ - the operation can propose a force condition in the PostStepGPIL (the passed value
119 // ------ to the operation is the one of the wrapped process, if proposeForceCondition is
120 // ------ unchanged, this same value will be used as the biasing foroce condition)
122 G4ForceCondition& /* proposeForceCondition */ ) = 0;
123 // ----
124 // ---- I.2) Methods called in at the AlongStepGetPhysicalInteractionLength(...) level :
125 // ----
126 // ------ o Operation can optionnally limit GPIL Along Step:
127 virtual G4double ProposeAlongStepLimit( const G4BiasingProcessInterface* /* callingProcess */ ) { return DBL_MAX; }
128 // ------ o Operation can propose a GPILSelection in the AlongStepGPIL
129 // ------ this selection superseeded the wrapped process selection
130 // ------ if the wrapped process exists, and if has along methods:
131 virtual G4GPILSelection ProposeGPILSelection( const G4GPILSelection wrappedProcessSelection )
132 {return wrappedProcessSelection;}
133
134 // ----
135 // ---- I.3) Methods called in at the AlongStepDoIt(...) level :
136 // ----
137 // ------ o Helper method to inform the operation of the move made in the along, and related non-interaction weight
138 // ------ applied to the primary track for this move:
139 virtual void AlongMoveBy( const G4BiasingProcessInterface* /* callingProcess */,
140 const G4Step* /* step */,
141 G4double /* weightForNonInteraction */ ) {}
142
143
144 // ---- II. Biasing of the process post step final state:
145 // ------------------------------------------------------
146 // ------ Mandatory method for biasing of the PostStepDoIt of the wrapped process
147 // ------ holds by the G4BiasingProcessInterface callingProcess.
148 // ------ User has full freedom for the particle change returned, and is reponsible for
149 // ------ the correctness of weights set to tracks.
150 // ------ The forcedBiasedFinalState should be left as is (ie false) in general. In this
151 // ------ way, if an occurrence biasing is also applied in the step, the weight correction
152 // ------ for it will be applied. If returned forceBiasedFinalState is returned true, then
153 // ------ the returned particle change will be returned as is to the stepping. Full
154 // ------ responsibility of the weight correctness is taken by the biasing operation.
155 // ------ The wrappedProcess can be accessed through the G4BiasingProcessInterface if needed.
156 // ------ This can be used in conjunction with an occurrence biasing, provided this final
157 // ------ state biasing is uncorrelated with the occurrence biasing (as single multiplication
158 // ------ of weights occur between these two biasings).
160 const G4Track* /* track */,
161 const G4Step* /* step */,
162 G4bool& /* forceBiasedFinalState */) = 0;
163
164
165 // ---- III. Biasing of the process along step final state:
166 // --------------------------------------------------------
167 // ---- Unprovided for now : requires significant developments.
168
169
170
171 // ***************************************************
172 // -- Methods for non-physics-based biasing operation:
173 // ***************************************************
174 // ----
175 // ---- If the operation is returned to the G4BiasingProcessInterface process by the
176 // ---- ProposeNonPhysicsBiasingOperation(...)/GetProposedNonPhysicsBiasingOperation(...) method
177 // ---- of the biasing operator, all methods below will be called for this operation.
178 // -----
179 // ---- 1) Method called in at the PostStepGetPhysicalInteractionLength(...) level :
180 // ----
181 // ---- o Return to the distance at which the operation should be applied, or may
182 // ---- play with the force condition flags.
183 virtual G4double DistanceToApplyOperation( const G4Track* /* track */,
184 G4double /* previousStepSize */,
185 G4ForceCondition* /* condition */) = 0;
186 // ----
187 // ---- 2) Method called in at the PostStepDoIt(...) level :
188 // ----
189 // ---- o Generate the final state for biasing (eg: splitting, killing, etc.)
191 const G4Step* /* step */) = 0;
192
193
194 // ----------------------------------------
195 // -- public interface and utility methods:
196 // ----------------------------------------
197public:
198 const G4String& GetName() const {return fName;}
199 std::size_t GetUniqueID() const {return fUniqueID;}
200
201
202private:
203 const G4String fName;
204 // -- better would be to have fUniqueID const, but pb on windows with constructor.
205 std::size_t fUniqueID;
206};
207
208#endif
G4ForceCondition
G4GPILSelection
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
Definition: G4Step.hh:62
virtual G4GPILSelection ProposeGPILSelection(const G4GPILSelection wrappedProcessSelection)
virtual void AlongMoveBy(const G4BiasingProcessInterface *, const G4Step *, G4double)
const G4String & GetName() const
virtual G4double DistanceToApplyOperation(const G4Track *, G4double, G4ForceCondition *)=0
std::size_t GetUniqueID() const
virtual G4VParticleChange * ApplyFinalStateBiasing(const G4BiasingProcessInterface *, const G4Track *, const G4Step *, G4bool &)=0
virtual G4VParticleChange * GenerateBiasingFinalState(const G4Track *, const G4Step *)=0
virtual const G4VBiasingInteractionLaw * ProvideOccurenceBiasingInteractionLaw(const G4BiasingProcessInterface *, G4ForceCondition &)=0
virtual G4double ProposeAlongStepLimit(const G4BiasingProcessInterface *)
#define DBL_MAX
Definition: templates.hh:62