Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4VBiasingOperator.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 pilot the biasing in a logical volume. This
35 // class is for *making decisions* on biasing operations to be applied.
36 // These ones are represented by the G4VBiasingOperation class.
37 // The volume in which biasing is applied is specified by the
38 // AttachTo(const G4LogicalVolume *) method. This has to be specified
39 // at detector construction time in the method ConstructSDandField() of
40 // G4VUsedDetectorConstruction.
41 //
42 // At tracking time the biasing operator is messaged by each
43 // G4BiasingProcessInterface object attached to the current track. For
44 // example, if three physics processes are under biasing, and if an
45 // additional G4BiasingProcessInterface is present to handle non-physics
46 // based biasing (splitting, killing), the operator will be messaged by
47 // these four G4BiasingProcessInterface objects.
48 // The idendity of the calling G4BiasingProcessInterface is known
49 // to the G4VBiasingOperator by passing this process pointer to the
50 // operator.
51 //
52 // ** Mandatory methods: **
53 //
54 // Three types of biasing are to be decided by the G4VBiasingOperator:
55 //
56 // 1) non-physics-based biasing:
57 // -----------------------------
58 // Meant for pure killing/splitting/etc. biasing operations, not
59 // associated to a physics process:
60 //
61 // virtual G4VBiasingOperation* ProposeNonPhysicsBiasingOperation( const G4Track* track,
62 // const G4BiasingProcessInterface* callingProcess ) = 0;
63 //
64 // Arguments are the current track, and the G4BiasingProcessInterface
65 // pointer making the call to the operator. In this case, this process
66 // does not wrap a physics process and
67 // callingProcess->GetWrappedProcess() == 0.
68 //
69 // The G4VBiasingOperation pointer returned is the operation to be
70 // applied. Zero can be returned. This operation will limit the
71 // step and propose a final state.
72 //
73 // This method is the first operator method called, it is called at the
74 // by the PostStepGetPhysicalInterationLenght(...) method of the
75 // G4BiasingProcessInterface.
76 //
77 // 2) physics-based biasing:
78 // -------------------------
79 // Physics-based biasing operations are of two types:
80 // - biasing of the physics process occurence interaction law
81 // - biasing of the physics process final state production
82 //
83 // a) The biasing of the occurence interaction law is proposed by:
84 //
85 // virtual G4VBiasingOperation* ProposeOccurenceBiasingOperation( const G4Track* track,
86 // const G4BiasingProcessInterface* callingProcess ) = 0;
87 // The current G4Track pointer and the G4BiasingProcessInterface
88 // pointer of the process calling the operator are passed. The
89 // G4BiasingProcessInterface process wraps an actual physics process
90 // which pointer can be obtained with
91 // callingProcess->GetWrappedProcess() .
92 //
93 // The biasing operation returned will be asked for its biasing
94 // interaction by the calling process, which will be a const object
95 // for the process. All setup and sampling regarding this law should be done
96 // in the operator before returning the related operation to the process.
97 //
98 // This method is the second operator one called in a step, it is called by
99 // the PostStepGetPhysicalInterationLenght(...) method of the
100 // G4BiasingProcessInterface.
101 //
102 // b) The biasing of the physics process final state is proposed by:
103 //
104 // virtual G4VBiasingOperation* ProposeFinalStateBiasingOperation( const G4Track* track,
105 // const G4BiasingProcessInterface* callingProcess ) = 0;
106 //
107 // The operator can propose a biasing operation that will handle the
108 // physic process final state biasing. As in previous case a) the
109 // G4BiasingProcessInterface process wraps an actual physics process
110 // which pointer can be obtained with:
111 // callingProcess->GetWrappedProcess() .
112 //
113 // Cases a) and b) are handled independently, and one or two of these
114 // biasing types can be provided in the same step.
115 //
116 // This method is the last operator one called in a step, it is called
117 // by the PostStepDoIt(...) method of the G4BiasingProcessInterface.
118 //
119 //
120 // ** Optional methods: **
121 //
122 // At the end of the step, the operator is messaged by the G4BiasingProcessInterface
123 // for operation(s) which have been applied during the step. One of the two following
124 // methods is called:
125 //
126 // virtual void OperationApplied( const G4BiasingProcessInterface* callingProcess,
127 // G4BiasingAppliedCase biasingCase,
128 // G4VBiasingOperation* operationApplied,
129 // const G4VParticleChange* particleChangeProduced );
130 // At most a single biasing operation was applied by the process:
131 // - a non-physics biasing operation was applied, biasingCase == BAC_NonPhysics ;
132 // - physics-based biasing:
133 // - the operator requested no biasing operations, and did let the physics
134 // process go : biasingCase == BAC_None;
135 // - a single final state biasing was proposed, with no concomittant occurence:
136 // biasingCase == BAC_FinalState;
137 // The operation applied and final state passed to the tracking (particleChangeProduced) are
138 // passed as information to the operator.
139 //
140 // virtual void OperationApplied( const G4BiasingProcessInterface* callingProcess,
141 // G4BiasingAppliedCase biasingCase,
142 // G4VBiasingOperation* occurenceOperationApplied,
143 // G4double weightForOccurenceInteraction,
144 // G4VBiasingOperation* finalStateOperationApplied,
145 // const G4VParticleChange* particleChangeProduced );
146 // This method is called in case an occurence biasing operation has been applied during the step.
147 // The biasingCase value is then the one of the final state biasing, if any : depending on if the
148 // occurence operation was applied alone and together with a final state operation, the
149 // biasingCase will take values:
150 // - occurence biasing alone : biasingCase == BAC_None ;
151 // in which case finalStateOperationApplied == 0;
152 // - occurence biasing + final state biasing : biasingCase == BAC_FinalState;
153 // The particleChangeProduced is the one *before* application of the weight for occurence : hence
154 // either the particle change of the (analog) physics process, or the biased final state, resulting
155 // from the biasing by the finalStateOperationApplied operation.
156 //
157 //
158 // ----------------G4VBiasingOperation ----------------
159 //
160 // Author: M.Verderi (LLR), November 2013
161 //
162 // --------------------------------------------------------------------
163 
164 #ifndef G4VBiasingOperator_hh
165 #define G4VBiasingOperator_hh 1
166 
167 #include "globals.hh"
168 
169 class G4VBiasingOperation;
170 class G4Track;
172 class G4LogicalVolume;
173 class G4VParticleChange;
175 #include <map>
176 #include <vector>
177 #include "G4BiasingAppliedCase.hh"
178 #include "G4Cache.hh"
179 
180 
182 
183  // -- State machine used to inform operators
184  // -- about run starting.
185  // -- Defined at the end of this file.
187 
188 public:
189  // ---------------
190  // -- Constructor:
191  // ---------------
193  virtual ~G4VBiasingOperator();
194 
195  // ----------------------------------------------
196  // -- abstract and user interface to sub-classes:
197  // ----------------------------------------------
198 protected:
199  // -- mandatory methods to let the operator tell about biasing operations to be applied:
200  // -------------------------------------------------------------------------------------
201  // -- These three methods have the same arguments passed : the current G4Track pointer, and the pointer of the
202  // -- G4BiasingProcessInterface instance calling this biasing operator. This same biasing operator will be called by each
203  // -- of the G4BiasingProcessInterface instances, meaning for example that:
204  // -- - if one G4BiasingProcessInterface with no wrapped physics process exits, ProposeNonPhysicsBiasingOperation(...)
205  // -- will be called one time at the beginning of the step,
206  // -- - if three G4BiasingProcessInterface instances exist, each of these one wrapping a physics process (eg
207  // -- conversion, Compton, photo-electric), ProposeOccurenceBiasingOperation(...) will be called three times,
208  // -- by each of these instances, at the beginning of the step and ProposeFinalStateBiasingOperation(...) will
209  // -- also be called by each of these instances, at the PostStepDoIt level.
210  // -- If a null pointer is returned, the analog -unbiased- behavior is adopted.
211  // -- non-physics-based biasing:
212  // -----------------------------
213  // -- [ First operator method called, at the PostStepGetPhysicalInterationLenght(...) level. ]
214  virtual G4VBiasingOperation* ProposeNonPhysicsBiasingOperation( const G4Track* track, const G4BiasingProcessInterface* callingProcess ) = 0;
215  // -- physics-based biasing:
216  // -------------------------
217  // -- Method to propose an occurence biasing operation : ie a change of the interaction length distribution. The proposed
218  // -- biasing operation will then be asked for its interaction law.
219  // -- Note that *** all sanity checks regarding the operation and its interaction law will have to have been performed
220  // -- before returning the biasing operation pointer *** as no corrective/aborting actions will be possible beyond this point.
221  // -- The informations provided by the G4BiasingProcessInterface calling process (previous occurence operation, previous step length,
222  // -- etc.) might be useful for doing this. They will be useful also to decide with continuing with a same operation proposed
223  // -- in the previous step, updating the interaction law taking into account the new G4Track state and the previous step size.
224  // -- [ Second operator method called, at the PostStepGetPhysicalInterationLenght(...) level. ]
225  virtual G4VBiasingOperation* ProposeOccurenceBiasingOperation( const G4Track* track, const G4BiasingProcessInterface* callingProcess ) = 0;
226  // -- [ Third operator method called, at the PostStepDoIt(...) level. ]
227  virtual G4VBiasingOperation* ProposeFinalStateBiasingOperation( const G4Track* track, const G4BiasingProcessInterface* callingProcess ) = 0;
228 
229 protected:
230  // -- optional methods for further information passed to the operator:
231  // -------------------------------------------------------------------
232  // ---- report to operator about the operation applied, the biasingCase value provides the case of biasing applied:
233  virtual void OperationApplied( const G4BiasingProcessInterface* callingProcess, G4BiasingAppliedCase biasingCase,
234  G4VBiasingOperation* operationApplied, const G4VParticleChange* particleChangeProduced );
235  // ---- same as above, report about the operation applied, for the case an occurence biasing was applied, together or not with a final state biasing.
236  // ---- The variable biasingCase tells if the final state is a biased one or not. **But in all cases**, this call happens only
237  // ---- for an occurence biaising : ie the occurence weight is applied on top of the particleChangeProduced, which is the particle
238  // ---- *before* the weight application for occurence biasing.
239  virtual void OperationApplied( const G4BiasingProcessInterface* callingProcess, G4BiasingAppliedCase biasingCase,
240  G4VBiasingOperation* occurenceOperationApplied, G4double weightForOccurenceInteraction,
241  G4VBiasingOperation* finalStateOperationApplied, const G4VParticleChange* particleChangeProduced );
242 protected:
243  // ---- method to inform operator that its biasing control is over (exit volume, or end of tracking):
244  // ---- [Called at the beginning of next step, or at the end of tracking.]
245  virtual void ExitBiasing( const G4Track* track, const G4BiasingProcessInterface* callingProcess );
246 
247 
248 protected:
249  // -----------------------------------
250  // -- Delegation to an other operator:
251  // -----------------------------------
252  // -- An operator may wish to select a sequence of operations already implemented in an
253  // -- existing biasing operator. In this case, this operator can delegate its work to
254  // -- the "delegated" one by calling DelegateTo( G4VBiasingOperation* delegated );
255  // -- §§ Should we have:
256  // -- §§ - a "step delegation" -where the delegation is made for the current step only-
257  // -- §§ - a long delegation where the delegation can hold over several steps, as long as
258  // -- §§ the scheme is not completed. [let's call it "scheme delegation"]
259  // -- §§ In this case the "execution/delegated" operator might switch off back the
260  // -- §§ delegation from the "delegator" when it knows it has done its work.
261  // -- §§ Add a private SetDelegator( G4VBiasingOperator* ) method, call on the delegated
262  // -- §§ operator.
263  // -- §§ For a step long delegation, the ReportOperationApplied should be used to "unset"
264  // -- §§ the delegation. For a scheme long delegation, the delegater operator will unset
265  // -- §§ itself has delegation. Likely to happen in the ReportOperationApplied as well,
266  // -- §§ but not sure it is mandatory though.
267 
268 
269 public:
270  // ---- Configure() is called in sequential mode or for master thread in MT mode.
271  // ---- It is in particular aimed at registering ID's to physics model at run initialization.
272  virtual void Configure() {}
273  // ---- ConfigureForWorker() is called in MT mode only, and only for worker threads.
274  // ---- It is not not to be used to register ID's to physics model catalog.
275  virtual void ConfigureForWorker() {}
276  // ---- inform the operator of the start of the run:
277  virtual void StartRun() {}
278  // ---- inform the operator of the start (end) of the tracking of a new track:
279  virtual void StartTracking( const G4Track* /* track */ ) {}
280  virtual void EndTracking() {}
281 
282 
283 
284  // --------------------
285  // -- public interface:
286  // --------------------
287  // -- needed by user:
288 public:
289  const G4String GetName() const {return fName;}
290  void AttachTo( const G4LogicalVolume* ); // -- attach to single volume
291 
292  G4BiasingAppliedCase GetPreviousBiasingAppliedCase() const {return fPreviousBiasingAppliedCase;}
293  // -- all operators (might got to a manager):
294  static const std::vector < G4VBiasingOperator* >& GetBiasingOperators() {return fOperators.Get();}
295  // -- get operator associated to a logical volume:
296  static G4VBiasingOperator* GetBiasingOperator( const G4LogicalVolume* ); // -- might go to a manager ; or moved to volume
297 
298 
299 
300  // -- used by biasing process interface, or used by an other operator (not expected to be invoked differently than with these two cases):
301 public:
305  void ExitingBiasing( const G4Track* track, const G4BiasingProcessInterface* callingProcess );
306 
307 public:
308  void ReportOperationApplied( const G4BiasingProcessInterface* callingProcess, G4BiasingAppliedCase biasingCase,
309  G4VBiasingOperation* operationApplied, const G4VParticleChange* particleChangeProduced );
310  void ReportOperationApplied( const G4BiasingProcessInterface* callingProcess, G4BiasingAppliedCase biasingCase,
311  G4VBiasingOperation* occurenceOperationApplied, G4double weightForOccurenceInteraction,
312  G4VBiasingOperation* finalStateOperationApplied, const G4VParticleChange* particleChangeProduced );
313 
314 
315 public:
316  const G4VBiasingOperation* GetPreviousNonPhysicsAppliedOperation() {return fPreviousAppliedNonPhysicsBiasingOperation;}
317 
318 
319 private:
320  const G4String fName;
321  // -- thread local:
322  // static std::map< const G4LogicalVolume*, G4VBiasingOperator* > fLogicalToSetupMap;
324  // -- thread local:
325  static G4VectorCache<G4VBiasingOperator* > fOperators;
326  // static std::vector < G4VBiasingOperator* > fOperators;
327 
328  // -- thread local:
329  static G4Cache< G4BiasingOperatorStateNotifier* > fStateNotifier;
330 
331 
332  // -- For this operator:
333  std::vector< const G4LogicalVolume* > fRootVolumes;
334  std::map < const G4LogicalVolume*, G4int > fDepthInTree;
335 
336  // -- current operation:
337  G4VBiasingOperation* fOccurenceBiasingOperation;
338  G4VBiasingOperation* fFinalStateBiasingOperation;
339  G4VBiasingOperation* fNonPhysicsBiasingOperation;
340 
341  // -- previous operations:
342  const G4VBiasingOperation* fPreviousProposedOccurenceBiasingOperation;
343  const G4VBiasingOperation* fPreviousProposedFinalStateBiasingOperation;
344  const G4VBiasingOperation* fPreviousProposedNonPhysicsBiasingOperation;
345  const G4VBiasingOperation* fPreviousAppliedOccurenceBiasingOperation;
346  const G4VBiasingOperation* fPreviousAppliedFinalStateBiasingOperation;
347  const G4VBiasingOperation* fPreviousAppliedNonPhysicsBiasingOperation;
348  G4BiasingAppliedCase fPreviousBiasingAppliedCase;
349 
350 };
351 
352 // -- state machine to get biasing operators
353 // -- messaged at the beginning of runs:
354 #include "G4VStateDependent.hh"
356 public:
359 public:
360  G4bool Notify(G4ApplicationState requestedState);
361 private:
362  G4ApplicationState fPreviousState;
363 };
364 
365 #endif
const XML_Char * name
Definition: expat.h:151
G4VBiasingOperation * GetProposedFinalStateBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)
virtual void OperationApplied(const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced)
G4String fName
Definition: G4AttUtils.hh:55
G4BiasingAppliedCase GetPreviousBiasingAppliedCase() const
value_type & Get() const
Definition: G4Cache.hh:282
G4VBiasingOperation * GetProposedOccurenceBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)
virtual G4VBiasingOperation * ProposeFinalStateBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0
G4bool Notify(G4ApplicationState requestedState)
virtual void ExitBiasing(const G4Track *track, const G4BiasingProcessInterface *callingProcess)
G4BiasingAppliedCase
virtual void Configure()
bool G4bool
Definition: G4Types.hh:79
virtual G4VBiasingOperation * ProposeOccurenceBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0
void AttachTo(const G4LogicalVolume *)
virtual G4VBiasingOperation * ProposeNonPhysicsBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0
static G4VBiasingOperator * GetBiasingOperator(const G4LogicalVolume *)
G4VBiasingOperator(G4String name)
virtual void StartTracking(const G4Track *)
static const std::vector< G4VBiasingOperator * > & GetBiasingOperators()
double G4double
Definition: G4Types.hh:76
void ReportOperationApplied(const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced)
G4VBiasingOperation * GetProposedNonPhysicsBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)
void ExitingBiasing(const G4Track *track, const G4BiasingProcessInterface *callingProcess)
virtual void ConfigureForWorker()
const G4VBiasingOperation * GetPreviousNonPhysicsAppliedOperation()
virtual void EndTracking()
G4ApplicationState
const G4String GetName() const