Geant4 Cross Reference |
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 // G4VBiasingOperation 27 // 28 // Class Description: 29 // 30 // An abstract class to model the behavior of any type of biasing : 31 // physics-based biasing (change of physics process behavior) or non- 32 // physics-based one, like splitting, killing. 33 // 34 // o The change of behavior of a physics process can be: 35 // - a change of the PostStep interaction probabilty, so-called 36 // occurrence biasing 37 // - a change in final state production 38 // - both, provided above two are uncorrelated. 39 // o The change of occurrence is driven by providing a biasing interaction 40 // law (G4VBiasingInteractionLaw) that is used in place of the analog 41 // exponential law. 42 // This change of occurrence is controlled through many handles. 43 // o The change in final state production is made through one single 44 // method the user is fully responsible of. 45 // 46 // o Non-physics-based biasing is controlled by two methods : one to 47 // specify where this biasing should happen, and one for generating 48 // the related final-state. 49 // 50 // Author: Marc Verderi (LLR), November 2013 51 // -------------------------------------------------------------------- 52 #ifndef G4VBiasingOperation_hh 53 #define G4VBiasingOperation_hh 1 54 55 #include "globals.hh" 56 #include "G4ForceCondition.hh" 57 #include "G4GPILSelection.hh" 58 59 class G4VParticleChange; 60 class G4Track; 61 class G4Step; 62 class G4VBiasingInteractionLaw; 63 class G4VProcess; 64 class G4BiasingProcessInterface; 65 66 class G4VBiasingOperation 67 { 68 public: 69 70 // -- Constructor: 71 G4VBiasingOperation(const G4String& name); 72 73 // -- destructor: 74 virtual ~G4VBiasingOperation() = default; 75 76 // ----------------------------- 77 // -- Interface to sub-classes : 78 // ----------------------------- 79 // -- 80 // ************************************* 81 // ** Methods for physics-based biasing: 82 // ************************************* 83 // -- 84 // ---- I. Biasing of the process occurrence: 85 // ----------------------------------------- 86 // ---- The biasing of the process occurrence regards the occurrence of the PostStepDoIt 87 // ---- behavior. But the weight is manipulated by both AlongStep methods (weight for 88 // ---- non-interaction) and PostStep methods (weight for interaction). For this 89 // ---- reason, occurrence biasing is handled by both AlongStep and PostStep methods. 90 // ---- 91 // ---- If the operation is returned to the G4BiasingProcessInterface process by the 92 // ---- ProposeOccurenceBiasingOperation(...)/GetProposedOccurenceBiasingOperation(...) method 93 // ---- of the biasing operator, all methods below will be called for this operation. 94 // ---- 95 // ---- I.1) Methods called in at the PostStepGetPhysicalInteractionLength(...) level : 96 // ---- 97 // ------ o Main and mandatory method for biasing of the PostStep process biasing occurrence : 98 // ------ - propose an interaction law to be substituted to the process that is biased 99 // ------ - the operation is told which is the G4BiasingProcessInterface calling it with 100 // ------ callingProcess argument. 101 // ------ - the returned law will have to have been sampled prior to be returned as it will be 102 // ------ asked for its GetSampledInteractionLength() by the callingProcess. 103 // ------ - the operation can propose a force condition in the PostStepGPIL (the passed value 104 // ------ to the operation is the one of the wrapped process, if proposeForceCondition is 105 // ------ unchanged, this same value will be used as the biasing foroce condition) 106 virtual const G4VBiasingInteractionLaw* 107 ProvideOccurenceBiasingInteractionLaw( const G4BiasingProcessInterface* /* callingProcess */ , 108 G4ForceCondition& /* proposeForceCondition */ ) = 0; 109 // ---- 110 // ---- I.2) Methods called in at the AlongStepGetPhysicalInteractionLength(...) level : 111 // ---- 112 // ------ o Operation can optionnally limit GPIL Along Step: 113 virtual G4double ProposeAlongStepLimit( const G4BiasingProcessInterface* /* callingProcess */ ) 114 { return DBL_MAX; } 115 116 // ------ o Operation can propose a GPILSelection in the AlongStepGPIL 117 // ------ this selection superseeded the wrapped process selection 118 // ------ if the wrapped process exists, and if has along methods: 119 virtual G4GPILSelection ProposeGPILSelection( const G4GPILSelection wrappedProcessSelection ) 120 { return wrappedProcessSelection; } 121 122 // ---- 123 // ---- I.3) Methods called in at the AlongStepDoIt(...) level : 124 // ---- 125 // ------ o Helper method to inform the operation of the move made in the along, and related non-interaction weight 126 // ------ applied to the primary track for this move: 127 virtual void AlongMoveBy( const G4BiasingProcessInterface* /* callingProcess */, 128 const G4Step* /* step */, 129 G4double /* weightForNonInteraction */ ) {} 130 131 // ---- II. Biasing of the process post step final state: 132 // ------------------------------------------------------ 133 // ------ Mandatory method for biasing of the PostStepDoIt of the wrapped process 134 // ------ holds by the G4BiasingProcessInterface callingProcess. 135 // ------ User has full freedom for the particle change returned, and is reponsible for 136 // ------ the correctness of weights set to tracks. 137 // ------ The forcedBiasedFinalState should be left as is (ie false) in general. In this 138 // ------ way, if an occurrence biasing is also applied in the step, the weight correction 139 // ------ for it will be applied. If returned forceBiasedFinalState is returned true, then 140 // ------ the returned particle change will be returned as is to the stepping. Full 141 // ------ responsibility of the weight correctness is taken by the biasing operation. 142 // ------ The wrappedProcess can be accessed through the G4BiasingProcessInterface if needed. 143 // ------ This can be used in conjunction with an occurrence biasing, provided this final 144 // ------ state biasing is uncorrelated with the occurrence biasing (as single multiplication 145 // ------ of weights occur between these two biasings). 146 virtual G4VParticleChange* ApplyFinalStateBiasing( const G4BiasingProcessInterface* /* callingProcess */, 147 const G4Track* /* track */, 148 const G4Step* /* step */, 149 G4bool& /* forceBiasedFinalState */) = 0; 150 151 // ---- III. Biasing of the process along step final state: 152 // -------------------------------------------------------- 153 // ---- Unprovided for now : requires significant developments. 154 155 // *************************************************** 156 // -- Methods for non-physics-based biasing operation: 157 // *************************************************** 158 // ---- 159 // ---- If the operation is returned to the G4BiasingProcessInterface process by the 160 // ---- ProposeNonPhysicsBiasingOperation(...)/GetProposedNonPhysicsBiasingOperation(...) method 161 // ---- of the biasing operator, all methods below will be called for this operation. 162 // ----- 163 // ---- 1) Method called in at the PostStepGetPhysicalInteractionLength(...) level : 164 // ---- 165 // ---- o Return to the distance at which the operation should be applied, or may 166 // ---- play with the force condition flags. 167 virtual G4double DistanceToApplyOperation( const G4Track* /* track */, 168 G4double /* previousStepSize */, 169 G4ForceCondition* /* condition */) = 0; 170 // ---- 171 // ---- 2) Method called in at the PostStepDoIt(...) level : 172 // ---- 173 // ---- o Generate the final state for biasing (eg: splitting, killing, etc.) 174 virtual G4VParticleChange* GenerateBiasingFinalState( const G4Track* /* track */, 175 const G4Step* /* step */) = 0; 176 177 // ---------------------------------------- 178 // -- public interface and utility methods: 179 // ---------------------------------------- 180 181 const G4String& GetName() const { return fName; } 182 std::size_t GetUniqueID() const { return fUniqueID; } 183 184 private: 185 186 const G4String fName; 187 // -- better would be to have fUniqueID const, but pb on windows with constructor. 188 std::size_t fUniqueID; 189 }; 190 191 #endif 192