Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-10-21 08:57:18

0001 //
0002 // ********************************************************************
0003 // * License and Disclaimer                                           *
0004 // *                                                                  *
0005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
0006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
0007 // * conditions of the Geant4 Software License,  included in the file *
0008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
0009 // * include a list of copyright holders.                             *
0010 // *                                                                  *
0011 // * Neither the authors of this software system, nor their employing *
0012 // * institutes,nor the agencies providing financial support for this *
0013 // * work  make  any representation or  warranty, express or implied, *
0014 // * regarding  this  software system or assume any liability for its *
0015 // * use.  Please see the license in the file  LICENSE  and URL above *
0016 // * for the full disclaimer and the limitation of liability.         *
0017 // *                                                                  *
0018 // * This  code  implementation is the result of  the  scientific and *
0019 // * technical work of the GEANT4 collaboration.                      *
0020 // * By using,  copying,  modifying or  distributing the software (or *
0021 // * any work based  on the software)  you  agree  to acknowledge its *
0022 // * use  in  resulting  scientific  publications,  and indicate your *
0023 // * acceptance of all terms of the Geant4 Software license.          *
0024 // ********************************************************************
0025 //
0026 // G4VBiasingOperator, G4BiasingOperatorStateNotifier
0027 //
0028 // Class Description:
0029 //
0030 // An abstract class to pilot the biasing in a logical volume. This
0031 // class is for *making decisions* on biasing operations to be applied.
0032 // These ones are represented by the G4VBiasingOperation class.
0033 // The volume in which biasing is applied is specified by the
0034 // AttachTo(const G4LogicalVolume *) method. This has to be specified
0035 // at detector construction time in the method ConstructSDandField() of
0036 // G4VUsedDetectorConstruction.
0037 //
0038 // At tracking time the biasing operator is messaged by each
0039 // G4BiasingProcessInterface object attached to the current track. For
0040 // example, if three physics processes are under biasing, and if an
0041 // additional G4BiasingProcessInterface is present to handle non-physics
0042 // based biasing (splitting, killing), the operator will be messaged by
0043 // these four G4BiasingProcessInterface objects.
0044 // The idendity of the calling G4BiasingProcessInterface is known
0045 // to the G4VBiasingOperator by passing this process pointer to the
0046 // operator.
0047 //
0048 // ** Mandatory methods: **
0049 //
0050 // Three types of biasing are to be decided by the G4VBiasingOperator:
0051 //
0052 //   1) non-physics-based biasing:
0053 //   -----------------------------
0054 //   Meant for pure killing/splitting/etc. biasing operations, not 
0055 //   associated to a physics process:
0056 //
0057 //   virtual G4VBiasingOperation* ProposeNonPhysicsBiasingOperation( const G4Track* track,
0058 //                                                                   const G4BiasingProcessInterface* callingProcess ) = 0;
0059 //
0060 //   Arguments are the current track, and the G4BiasingProcessInterface
0061 //   pointer making the call to the operator. In this case, this process
0062 //   does not wrap a physics process and 
0063 //                callingProcess->GetWrappedProcess() == 0.
0064 //
0065 //   The G4VBiasingOperation pointer returned is the operation to be
0066 //   applied. Zero can be returned. This operation will limit the
0067 //   step and propose a final state.
0068 //   
0069 //   This method is the first operator method called, it is called at the
0070 //   by the PostStepGetPhysicalInterationLength(...) method of the
0071 //   G4BiasingProcessInterface.
0072 //
0073 //   2) physics-based biasing:
0074 //   -------------------------
0075 //   Physics-based biasing operations are of two types:
0076 //     - biasing of the physics process occurrence interaction law
0077 //     - biasing of the physics process final state production
0078 //
0079 //   a) The biasing of the occurrence interaction law is proposed by:
0080 //
0081 //   virtual G4VBiasingOperation*  ProposeOccurenceBiasingOperation( const G4Track* track, 
0082 //                                                                   const G4BiasingProcessInterface* callingProcess ) = 0;
0083 //   The current G4Track pointer and the G4BiasingProcessInterface
0084 //   pointer of the process calling the operator are passed. The
0085 //   G4BiasingProcessInterface process wraps an actual physics process
0086 //   which pointer can be obtained with
0087 //                callingProcess->GetWrappedProcess() .
0088 //
0089 //   The biasing operation returned will be asked for its biasing
0090 //   interaction by the calling process, which will be a const object
0091 //   for the process. All setup and sampling regarding this law should be done
0092 //   in the operator before returning the related operation to the process.
0093 //
0094 //   This method is the second operator one called in a step, it is called by
0095 //   the PostStepGetPhysicalInterationLength(...) method of the
0096 //   G4BiasingProcessInterface.
0097 //
0098 //   b) The biasing of the physics process final state is proposed by:
0099 //
0100 //   virtual G4VBiasingOperation* ProposeFinalStateBiasingOperation( const G4Track* track,
0101 //                                                                   const G4BiasingProcessInterface* callingProcess ) = 0;
0102 //
0103 //   The operator can propose a biasing operation that will handle the
0104 //   physic process final state biasing. As in previous case a) the
0105 //   G4BiasingProcessInterface process wraps an actual physics process
0106 //   which pointer can be obtained with:
0107 //                callingProcess->GetWrappedProcess() .
0108 //
0109 //   Cases a) and b) are handled independently, and one or two of these
0110 //   biasing types can be provided in the same step.
0111 //
0112 //   This method is the last operator one called in a step, it is called
0113 //   by the PostStepDoIt(...) method of the G4BiasingProcessInterface.
0114 //
0115 //
0116 // ** Optional methods: **
0117 //
0118 //     At the end of the step, the operator is messaged by the G4BiasingProcessInterface
0119 // for operation(s) which have been applied during the step. One of the two following
0120 // methods is called:
0121 //
0122 // virtual void OperationApplied( const G4BiasingProcessInterface* callingProcess,
0123 //                                G4BiasingAppliedCase biasingCase,
0124 //                G4VBiasingOperation* operationApplied,
0125 //                                const G4VParticleChange* particleChangeProduced );
0126 // At most a single biasing operation was applied by the process:
0127 //    - a non-physics biasing operation was applied, biasingCase == BAC_NonPhysics ;
0128 //    - physics-based biasing:
0129 //      - the operator requested no biasing operations, and did let the physics
0130 //        process go : biasingCase ==  BAC_None;
0131 //      - a single final state biasing was proposed, with no concomittant occurrence:
0132 //        biasingCase ==  BAC_FinalState;
0133 // The operation applied and final state passed to the tracking (particleChangeProduced) are
0134 // passed as information to the operator.
0135 // 
0136 // virtual void OperationApplied( const G4BiasingProcessInterface* callingProcess,
0137 //                                G4BiasingAppliedCase biasingCase,
0138 //                G4VBiasingOperation* occurenceOperationApplied,
0139 //                                G4double weightForOccurenceInteraction,
0140 //                G4VBiasingOperation* finalStateOperationApplied,
0141 //                                const G4VParticleChange* particleChangeProduced );
0142 // This method is called in case an occurrence biasing operation has been applied during the step.
0143 // The biasingCase value is then the one of the final state biasing, if any : depending on if the
0144 // occurrence operation was applied alone and together with a final state operation, the
0145 // biasingCase will take values:
0146 //     - occurrence biasing alone : biasingCase == BAC_None ;
0147 //       in which case finalStateOperationApplied == 0;
0148 //     - occurrence biasing + final state biasing : biasingCase ==  BAC_FinalState;
0149 // The particleChangeProduced is the one *before* application of the weight for occurrence : hence
0150 // either the particle change of the (analog) physics process, or the biased final state, resulting
0151 // from the biasing by the finalStateOperationApplied operation.
0152 //
0153 // Author: M.Verderi (LLR), November 2013
0154 // --------------------------------------------------------------------
0155 #ifndef G4VBiasingOperator_hh
0156 #define G4VBiasingOperator_hh 1
0157 
0158 #include "globals.hh"
0159 #include "G4BiasingAppliedCase.hh"
0160 #include "G4Cache.hh"
0161 #include "G4VStateDependent.hh"
0162 
0163 #include <map>
0164 #include <vector>
0165 
0166 class G4VBiasingOperation;
0167 class G4Track;
0168 class G4BiasingProcessInterface;
0169 class G4LogicalVolume;
0170 class G4VParticleChange;
0171 class G4BiasingOperatorStateNotifier;
0172 
0173 class G4VBiasingOperator
0174 {
0175   // -- State machine used to inform operators
0176   // -- about run starting.
0177   // -- Defined at the end of this file.
0178   friend class G4BiasingOperatorStateNotifier;
0179   
0180   public:
0181 
0182     // ---------------
0183     // -- Constructor:
0184     // ---------------
0185     G4VBiasingOperator(const G4String& name);
0186     virtual ~G4VBiasingOperator() = default;
0187   
0188     // ---- Configure() is called in sequential mode or for master thread in MT mode.
0189     // ---- It is in particular aimed at registering ID's to physics model at run initialization.
0190     virtual void Configure() {}
0191     // ---- ConfigureForWorker() is called in MT mode only, and only for worker threads.
0192     // ---- It is not not to be used to register ID's to physics model catalog.
0193     virtual void ConfigureForWorker() {}
0194     // ---- inform the operator of the start of the run:
0195     virtual void StartRun() {}
0196     // ---- inform the operator of the start (end) of the tracking of a new track:
0197     virtual void StartTracking( const G4Track* /* track */ ) {}
0198     virtual void EndTracking() {}
0199 
0200     // --------------------
0201     // -- public interface:
0202     // --------------------
0203     // -- needed by user:
0204     const G4String& GetName() const { return fName; }
0205     void AttachTo( const G4LogicalVolume* ); // -- attach to single volume 
0206 
0207     G4BiasingAppliedCase GetPreviousBiasingAppliedCase() const
0208       { return fPreviousBiasingAppliedCase; }
0209     // -- all operators (might got to a manager):
0210     static const std::vector < G4VBiasingOperator* >& GetBiasingOperators();
0211       // -- get operator associated to a logical volume:
0212     static G4VBiasingOperator* GetBiasingOperator( const G4LogicalVolume* );
0213       // -- might go to a manager ; or moved to volume
0214 
0215     // -- used by biasing process interface, or used by another operator (not expected to be invoked differently than with these two cases):
0216     G4VBiasingOperation* GetProposedOccurenceBiasingOperation( const G4Track* track,
0217                          const G4BiasingProcessInterface* callingProcess );
0218     G4VBiasingOperation* GetProposedFinalStateBiasingOperation( const G4Track* track,
0219                          const G4BiasingProcessInterface* callingProcess );
0220     G4VBiasingOperation* GetProposedNonPhysicsBiasingOperation( const G4Track* track,
0221                          const G4BiasingProcessInterface* callingProcess );
0222     void ExitingBiasing( const G4Track* track,
0223                          const G4BiasingProcessInterface* callingProcess );
0224 
0225     void ReportOperationApplied( const G4BiasingProcessInterface* callingProcess,
0226                          G4BiasingAppliedCase biasingCase,
0227                          G4VBiasingOperation* operationApplied,
0228                          const G4VParticleChange* particleChangeProduced );
0229     void ReportOperationApplied( const G4BiasingProcessInterface* callingProcess,
0230                          G4BiasingAppliedCase biasingCase,
0231                          G4VBiasingOperation* occurenceOperationApplied,
0232                          G4double weightForOccurenceInteraction,
0233                          G4VBiasingOperation* finalStateOperationApplied,
0234                          const G4VParticleChange* particleChangeProduced );
0235   
0236     const G4VBiasingOperation* GetPreviousNonPhysicsAppliedOperation()
0237       { return  fPreviousAppliedNonPhysicsBiasingOperation; }
0238 
0239   protected:
0240 
0241     // -- mandatory methods to let the operator tell about biasing operations to be applied:
0242     // -------------------------------------------------------------------------------------
0243     // -- These three methods have the same arguments passed : the current G4Track pointer, and the pointer of the
0244     // -- G4BiasingProcessInterface instance calling this biasing operator. This same biasing operator will be called by each
0245     // -- of the G4BiasingProcessInterface instances, meaning for example that:
0246     // --     - if one G4BiasingProcessInterface with no wrapped physics process exits, ProposeNonPhysicsBiasingOperation(...)
0247     // --       will be called one time at the beginning of the step,
0248     // --     - if three G4BiasingProcessInterface instances exist, each of these one wrapping a physics process (eg
0249     // --       conversion, Compton, photo-electric), ProposeOccurenceBiasingOperation(...) will be called three times,
0250     // --       by each of these instances, at the beginning of the step and ProposeFinalStateBiasingOperation(...) will
0251     // --       also be called by each of these instances, at the PostStepDoIt level.
0252     // -- If a null pointer is returned, the analog -unbiased- behavior is adopted.
0253     // -- non-physics-based biasing:
0254     // -----------------------------
0255     // -- [ First operator method called, at the PostStepGetPhysicalInterationLength(...) level. ]
0256     virtual G4VBiasingOperation* ProposeNonPhysicsBiasingOperation( const G4Track* track, const G4BiasingProcessInterface* callingProcess ) = 0;
0257 
0258     // -- physics-based biasing:
0259     // -------------------------
0260     // -- Method to propose an occurrence biasing operation : ie a change of the interaction length distribution. The proposed
0261     // -- biasing operation will then be asked for its interaction law.
0262     // -- Note that *** all sanity checks regarding the operation and its interaction law will have to have been performed
0263     // -- before returning the biasing operation pointer *** as no corrective/aborting actions will be possible beyond this point.
0264     // -- The informations provided by the G4BiasingProcessInterface calling process (previous occurrence operation, previous step length,
0265     // -- etc.) might be useful for doing this. They will be useful also to decide with continuing with a same operation proposed
0266     // -- in the previous step, updating the interaction law taking into account the new G4Track state and the previous step size.
0267     // -- [ Second operator method called, at the PostStepGetPhysicalInterationLength(...) level. ]
0268     virtual G4VBiasingOperation* ProposeOccurenceBiasingOperation( const G4Track* track, const G4BiasingProcessInterface* callingProcess ) = 0;
0269 
0270     // -- [ Third operator method called, at the PostStepDoIt(...) level. ]
0271     virtual G4VBiasingOperation* ProposeFinalStateBiasingOperation( const G4Track* track, const G4BiasingProcessInterface* callingProcess ) = 0;
0272 
0273     // -- optional methods for further information passed to the operator:
0274     // -------------------------------------------------------------------
0275     // ---- report to operator about the operation applied, the biasingCase value provides the case of biasing applied:
0276     virtual void OperationApplied( const G4BiasingProcessInterface* callingProcess, G4BiasingAppliedCase biasingCase,
0277                    G4VBiasingOperation* operationApplied, const G4VParticleChange* particleChangeProduced );
0278     // ---- same as above, report about the operation applied, for the case an occurrence biasing was applied, together or not with a final state biasing.
0279     // ---- The variable biasingCase tells if the final state is a biased one or not. **But in all cases**, this call happens only
0280     // ---- for an occurrence biasing : ie the occurrence weight is applied on top of the particleChangeProduced, which is the particle
0281     // ---- *before* the weight application for occurence biasing.
0282     virtual void OperationApplied( const G4BiasingProcessInterface* callingProcess, G4BiasingAppliedCase biasingCase,
0283                    G4VBiasingOperation* occurenceOperationApplied, G4double weightForOccurenceInteraction,
0284                    G4VBiasingOperation* finalStateOperationApplied, const G4VParticleChange* particleChangeProduced );
0285 
0286     // ---- method to inform operator that its biasing control is over (exit volume, or end of tracking):
0287     // ---- [Called at the beginning of next step, or at the end of tracking.]
0288     virtual void ExitBiasing( const G4Track* track, const G4BiasingProcessInterface* callingProcess );
0289 
0290     // ----------------------------------
0291     // -- Delegation to another operator:
0292     // ----------------------------------
0293     // -- An operator may wish to select a sequence of operations already implemented in an
0294     // -- existing biasing operator. In this case, this operator can delegate its work to
0295     // -- the "delegated" one by calling DelegateTo( G4VBiasingOperation* delegated );
0296     // -- Should we have:
0297     // --    - a "step delegation" -where the delegation is made for the current step only-
0298     // --    - a long delegation where the delegation can hold over several steps, as long as
0299     // --      the scheme is not completed. [let's call it "scheme delegation"]
0300     // --      In this case the "execution/delegated" operator might switch off back the
0301     // --      delegation from the "delegator" when it knows it has done its work.
0302     // -- Add a private SetDelegator( G4VBiasingOperator* ) method, call on the delegated
0303     // -- operator.
0304     // -- For a step long delegation, the ReportOperationApplied should be used to "unset"
0305     // -- the delegation. For a scheme long delegation, the delegater operator will unset
0306     // -- itself has delegation. Likely to happen in the ReportOperationApplied as well,
0307     // -- but not sure it is mandatory though.
0308   
0309   private:
0310 
0311     const G4String fName;
0312     // -- thread local:
0313     static G4MapCache< const G4LogicalVolume*, G4VBiasingOperator* > fLogicalToSetupMap;
0314     // -- thread local:
0315     static G4VectorCache<G4VBiasingOperator* > fOperators;
0316 
0317     // -- thread local:
0318     static G4Cache< G4BiasingOperatorStateNotifier* > fStateNotifier;
0319 
0320     // -- For this operator:
0321     std::vector< const G4LogicalVolume* > fRootVolumes;
0322     std::map< const G4LogicalVolume*, G4int > fDepthInTree;
0323   
0324     // -- current operation:
0325     G4VBiasingOperation* fOccurenceBiasingOperation = nullptr;
0326     G4VBiasingOperation* fFinalStateBiasingOperation = nullptr;
0327     G4VBiasingOperation* fNonPhysicsBiasingOperation = nullptr;
0328   
0329     // -- previous operations:
0330     const G4VBiasingOperation* fPreviousProposedOccurenceBiasingOperation = nullptr;
0331     const G4VBiasingOperation* fPreviousProposedFinalStateBiasingOperation = nullptr;
0332     const G4VBiasingOperation* fPreviousProposedNonPhysicsBiasingOperation = nullptr;
0333     const G4VBiasingOperation* fPreviousAppliedOccurenceBiasingOperation = nullptr;
0334     const G4VBiasingOperation* fPreviousAppliedFinalStateBiasingOperation = nullptr;
0335     const G4VBiasingOperation* fPreviousAppliedNonPhysicsBiasingOperation = nullptr;
0336     G4BiasingAppliedCase fPreviousBiasingAppliedCase = BAC_None;
0337 };
0338 
0339 // -- state machine to get biasing operators
0340 // -- messaged at the beginning of runs:
0341 
0342 class G4BiasingOperatorStateNotifier : public G4VStateDependent
0343 {
0344   public:
0345 
0346     G4BiasingOperatorStateNotifier();
0347     ~G4BiasingOperatorStateNotifier() = default;
0348 
0349     G4bool Notify(G4ApplicationState requestedState);
0350 
0351   private:
0352 
0353     G4ApplicationState fPreviousState;
0354 };
0355 
0356 #endif