Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/biasing/management/include/G4VBiasingOperator.hh

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

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