Geant4_10
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 // ---------------------------------------------------------------------
59 
60 
61 
62 #ifndef G4VBiasingOperation_hh
63 #define G4VBiasingOperation_hh 1
64 
65 #include "globals.hh"
66 class G4VParticleChange;
67 class G4Track;
68 class G4Step;
70 class G4VProcess;
72 #include "G4ForceCondition.hh"
73 #include "G4GPILSelection.hh"
74 
76 public:
77  // ---------------
78  // -- Constructor:
79  // ---------------
80  //
81  // -- Constructor for biasing operations:
82  // -----------------------------------
83  // -- Operation is given a name.
85 
86  // -- destructor:
87  virtual ~G4VBiasingOperation();
88 
89 public:
90  // -----------------------------
91  // -- Interface to sub-classes :
92  // -----------------------------
93  // --
94  // *************************************
95  // ** Methods for physics-based biasing:
96  // *************************************
97  // --
98  // ---- I. Biasing of the process occurence:
99  // -----------------------------------------
100  // ---- The biasing of the process occurence regards the occurence of the PostStepDoIt
101  // ---- behavior. But the weight is manipulated by both AlongStep methods (weight for
102  // ---- non-interaction) and PostStep methods (weight for interaction). For this
103  // ---- reason, occurence biasing is handled by both AlongStep and PostStep methods.
104  // ----
105  // ---- If the operation is returned to the G4BiasingProcessInterface process by the
106  // ---- ProposeOccurenceBiasingOperation(...)/GetProposedOccurenceBiasingOperation(...) method
107  // ---- of the biasing operator, all methods below will be called for this operation.
108  // ----
109  // ---- I.1) Methods called in at the PostStepGetPhysicalInteractionLength(...) level :
110  // ----
111  // ------ o Main and mandatory method for biasing of the PostStep process biasing occurence :
112  // ------ - propose an interaction law to be substituted to the process that is biased
113  // ------ - the operation is told which is the G4BiasingProcessInterface calling it with
114  // ------ callingProcess argument.
115  // ------ - the returned law will have to have been sampled prior to be returned as it will be
116  // ------ asked for its GetSampledInteractionLength() by the callingProcess.
117  virtual const G4VBiasingInteractionLaw* ProvideOccurenceBiasingInteractionLaw( const G4BiasingProcessInterface* /* callingProcess */ ) = 0;
118  // ------ o Additionnal method for PostStep occurence:
119  // ------ - operation can propose a force condition in the PostStepGPIL
120  // ------ - this condition superseedes the wrapped process condition (taken as default condition)
121  // ------ - it is called after above ProvideOccurenceBiasingInteractionLaw(...) by the same calling process.
122  virtual G4ForceCondition ProposeForceCondition( const G4ForceCondition wrappedProcessCondition )
123  {
124  return wrappedProcessCondition;
125  }
126  // ----
127  // ---- I.2) Methods called in at the AlongStepGetPhysicalInteractionLength(...) level :
128  // ----
129  // ------ o Operation can optionnally limit GPIL Along Step:
130  virtual G4double ProposeAlongStepLimit( const G4BiasingProcessInterface* /* callingProcess */ ) { return DBL_MAX; }
131  // ------ o Operation can propose a GPILSelection in the AlongStepGPIL
132  // ------ this selection superseeded the wrapped process selection
133  // ------ if the wrapped process exists, and if has along methods:
134  virtual G4GPILSelection ProposeGPILSelection( const G4GPILSelection wrappedProcessSelection )
135  {return wrappedProcessSelection;}
136 
137  // ----
138  // ---- I.3) Methods called in at the AlongStepDoIt(...) level :
139  // ----
140  // ------ o Helper method to inform the operation of the move made in the along, and related non-interaction weight
141  // ------ applied to the primary track for this move:
142  virtual void AlongMoveBy( const G4BiasingProcessInterface* /* callingProcess */,
143  const G4Step* /* step */,
144  G4double /* weightForNonInteraction */ ) {}
145 
146  // ----
147  // ---- I.4) Methods called in at the PostStepDoIt(...) level :
148  // ----
149  // ------ o Method allowing the operation to prevent wrapped process interaction to happen:
150  // ------ - has to return true if interaction is denied by the operation
151  // ------ - the primary is left unchanged in this case
152  // ------ - but its weight becomes the one given by proposedTrackWeight.
153  virtual G4bool DenyProcessPostStepDoIt( const G4BiasingProcessInterface* /* callingProcess */,
154  const G4Track* /* track */,
155  const G4Step* /* step */,
156  G4double& /* proposedTrackWeight */ )
157  {return false;}
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 wrappedProcess can be accessed through the G4BiasingProcessInterface if needed.
167  // ------ This can be used in conjonction with an occurence biasing, provided this final
168  // ------ state biasing is uncorrelated with the occurence biasing (as single multiplication
169  // ------ of weights occur between these two biasings).
170  virtual G4VParticleChange* ApplyFinalStateBiasing( const G4BiasingProcessInterface* /* callingProcess */,
171  const G4Track* /* track */,
172  const G4Step* /* step */) = 0;
173 
174 
175  // ---- III. Biasing of the process along step final state:
176  // --------------------------------------------------------
177  // ---- Unprovided for now : requires significant developments.
178 
179 
180 
181  // ***************************************************
182  // -- Methods for non-physics-based biasing operation:
183  // ***************************************************
184  // ----
185  // ---- If the operation is returned to the G4BiasingProcessInterface process by the
186  // ---- ProposeNonPhysicsBiasingOperation(...)/GetProposedNonPhysicsBiasingOperation(...) method
187  // ---- of the biasing operator, all methods below will be called for this operation.
188  // -----
189  // ---- 1) Method called in at the PostStepGetPhysicalInteractionLength(...) level :
190  // ----
191  // ---- o Return to the distance at which the operation should be applied, or may
192  // ---- play with the force condition flags.
193  virtual G4double DistanceToApplyOperation( const G4Track* /* track */,
194  G4double /* previousStepSize */,
195  G4ForceCondition* /* condition */) = 0;
196  // ----
197  // ---- 2) Method called in at the PostStepDoIt(...) level :
198  // ----
199  // ---- o Generate the final state for biasing (eg: splitting, killing, etc.)
200  virtual G4VParticleChange* GenerateBiasingFinalState( const G4Track* /* track */,
201  const G4Step* /* step */) = 0;
202 
203 
204  // ----------------------------------------
205  // -- public interface and utility methods:
206  // ----------------------------------------
207 public:
208  const G4String& GetName() const {return fName;}
209  std::size_t GetUniqueID() const {return fUniqueID;}
210 
211 
212 private:
213  const G4String fName;
214  // -- better would be to have fUniqueID const, but pb on windows with constructor.
215  std::size_t fUniqueID;
216 };
217 
218 #endif
virtual G4bool DenyProcessPostStepDoIt(const G4BiasingProcessInterface *, const G4Track *, const G4Step *, G4double &)
const XML_Char * name
Definition: expat.h:151
virtual G4VParticleChange * ApplyFinalStateBiasing(const G4BiasingProcessInterface *, const G4Track *, const G4Step *)=0
virtual G4double ProposeAlongStepLimit(const G4BiasingProcessInterface *)
bool G4bool
Definition: G4Types.hh:79
virtual G4ForceCondition ProposeForceCondition(const G4ForceCondition wrappedProcessCondition)
Definition: G4Step.hh:76
virtual G4VParticleChange * GenerateBiasingFinalState(const G4Track *, const G4Step *)=0
const G4String & GetName() const
virtual G4GPILSelection ProposeGPILSelection(const G4GPILSelection wrappedProcessSelection)
virtual void AlongMoveBy(const G4BiasingProcessInterface *, const G4Step *, G4double)
virtual const G4VBiasingInteractionLaw * ProvideOccurenceBiasingInteractionLaw(const G4BiasingProcessInterface *)=0
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