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