Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:59:19

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