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 // G4BiasingProcessInterface 27 // 28 // Class Description: 29 // 30 // A wrapper process making the interface between the tracking 31 // and the G4VBiasingOperator objects attached to volumes. 32 // If this process holds a physics process, it forwards 33 // tracking calls to this process in volume where not biasing 34 // occurs. In volumes with biasing (with a G4VBiasingOperator 35 // attached) the process gets what to do messaging the biasing 36 // operator : 37 // - at the PostStepGPIL level, for getting an occurrence biasing 38 // operation. If such an operation is returned to the process 39 // this operation will be messaged at several places. 40 // - at the PostStepDoIt level, to get a possible final state 41 // biasing operation 42 // If the process does not hold a physics process, it is meant 43 // as handling "non physics" biasing operations: pure splitting 44 // or pure kiling for example (ie not brem splitting). 45 // 46 // Author: Marc Verderi, September 2013. 47 // -------------------------------------------------------------------- 48 #ifndef G4BiasingProcessInterface_h 49 #define G4BiasingProcessInterface_h 1 50 51 #include "globals.hh" 52 #include "G4VProcess.hh" 53 #include "G4Cache.hh" 54 #include "G4BiasingProcessSharedData.hh" 55 56 class G4VBiasingInteractionLaw; 57 class G4InteractionLawPhysical; 58 class G4VBiasingOperator; 59 class G4VBiasingOperation; 60 class G4ParticleChangeForOccurenceBiasing; 61 class G4ParticleChangeForNothing; 62 63 class G4BiasingProcessInterface : public G4VProcess 64 { 65 public: 66 67 // -------------------------------------------------------------------------------- 68 // -- constructor for dealing with biasing options not affecting physics processes: 69 // -------------------------------------------------------------------------------- 70 G4BiasingProcessInterface( const G4String& name = "biasWrapper(0)" ); 71 // --------------------------------------------------------------- 72 // -- constructor to transform the behaviour of a physics process: 73 // --------------------------------------------------------------- 74 // -- wrappedProcess pointer MUST NOT be null. 75 G4BiasingProcessInterface(G4VProcess* wrappedProcess, 76 G4bool wrappedIsAtRest, 77 G4bool wrappedIsAlongStep, 78 G4bool wrappedIsPostStep, 79 G4String useThisName = ""); 80 ~G4BiasingProcessInterface(); 81 // -- pointer of wrapped physics process: 82 G4VProcess* GetWrappedProcess() const { return fWrappedProcess; } 83 84 // --------------------- 85 // -- Biasing interface: 86 // --------------------- 87 // -- Helper methods: 88 // ------------------ 89 // -- Current step and previous step biasing operator, if any: 90 G4VBiasingOperator* GetCurrentBiasingOperator() const 91 { return fSharedData->fCurrentBiasingOperator; } 92 G4VBiasingOperator* GetPreviousBiasingOperator() const 93 { return fSharedData->fPreviousBiasingOperator; } 94 // -- current and previous operation: 95 G4VBiasingOperation* GetCurrentNonPhysicsBiasingOperation() const 96 { return fNonPhysicsBiasingOperation; } 97 G4VBiasingOperation* GetPreviousNonPhysicsBiasingOperation() const 98 { return fPreviousNonPhysicsBiasingOperation; } 99 G4VBiasingOperation* GetCurrentOccurenceBiasingOperation() const 100 { return fOccurenceBiasingOperation; } 101 G4VBiasingOperation* GetPreviousOccurenceBiasingOperation() const 102 { return fPreviousOccurenceBiasingOperation; } 103 G4VBiasingOperation* GetCurrentFinalStateBiasingOperation() const 104 { return fFinalStateBiasingOperation; } 105 G4VBiasingOperation* GetPreviousFinalStateBiasingOperation() const 106 { return fPreviousFinalStateBiasingOperation; } 107 108 // -- Lists of processes cooperating under a same particle type/G4ProcessManager. 109 // -- The vector ordering is: 110 // -- - random, before first "run/beamOn" 111 // -- - that of the PostStepGetPhysicalInteractionLength() once "/run/beamOn" has been issued 112 const std::vector< const G4BiasingProcessInterface* >& GetBiasingProcessInterfaces() const 113 { return fSharedData->fPublicBiasingProcessInterfaces; } 114 const std::vector< const G4BiasingProcessInterface* >& GetPhysicsBiasingProcessInterfaces() const 115 { return fSharedData->fPublicPhysicsBiasingProcessInterfaces; } 116 const std::vector< const G4BiasingProcessInterface* >& GetNonPhysicsBiasingProcessInterfaces() const 117 { return fSharedData->fPublicNonPhysicsBiasingProcessInterfaces; } 118 119 // -- Get shared data for this process: 120 const G4BiasingProcessSharedData* GetSharedData() const 121 { return fSharedData; } 122 // -- Get shared data associated to a G4ProcessManager: 123 static const G4BiasingProcessSharedData* GetSharedData( const G4ProcessManager* ); 124 125 // ------------------ 126 // -- Helper methods: 127 // ------------------ 128 // ---- Tell is this process is first/last in the PostStep GPIL and DoIt lists. 129 // ---- If physOnly is true, only wrapper for physics processes are considered, 130 // ---- otherwise all G4BiasingProcessInterface processes of this particle are 131 // ---- considered. 132 // ---- These methods just return the corresponding flag values setup at 133 // ---- initialization phase by the next four ones. 134 // ---- Will not be updated if processes are activate/unactivated on the fly. 135 // ---- Use next methods (less fast) instead in this case. 136 G4bool GetIsFirstPostStepGPILInterface(G4bool physOnly = true) const; 137 G4bool GetIsLastPostStepGPILInterface(G4bool physOnly = true) const; 138 G4bool GetIsFirstPostStepDoItInterface(G4bool physOnly = true) const; 139 G4bool GetIsLastPostStepDoItInterface(G4bool physOnly = true) const; 140 // ---- Determine if the process is first/last in the PostStep GPIL and DoIt lists. 141 G4bool IsFirstPostStepGPILInterface(G4bool physOnly = true) const; 142 G4bool IsLastPostStepGPILInterface(G4bool physOnly = true) const; 143 G4bool IsFirstPostStepDoItInterface(G4bool physOnly = true) const; 144 G4bool IsLastPostStepDoItInterface(G4bool physOnly = true) const; 145 // -- Information about wrapped process: 146 G4bool GetWrappedProcessIsAtRest() const { return fWrappedProcessIsAtRest; } 147 G4bool GetWrappedProcessIsAlong() const { return fWrappedProcessIsAlong; } 148 G4bool GetWrappedProcessIsPost() const { return fWrappedProcessIsPost; } 149 150 // -- Information methods: 151 G4double GetPreviousStepSize() const { return fPreviousStepSize; } 152 G4double GetCurrentMinimumStep() const { return fCurrentMinimumStep; } 153 G4double GetProposedSafety() const { return fProposedSafety; } 154 void SetProposedSafety(G4double sft) { fProposedSafety = sft; } 155 // -- return the actual PostStep and AlongStep limits returned by the process to the tracking : 156 G4double GetPostStepGPIL() const { return fBiasingPostStepGPIL; } 157 G4double GetAlongStepGPIL() const { return fWrappedProcessAlongStepGPIL; } 158 159 // -------------------------------------------------------------- 160 // -- G4VProcess interface -- 161 // -------------------------------------------------------------- 162 163 // -- Start/End tracking: 164 void StartTracking(G4Track* track); 165 void EndTracking(); 166 167 // -- PostStep methods: 168 virtual G4double PostStepGetPhysicalInteractionLength(const G4Track& track, 169 G4double previousStepSize, 170 G4ForceCondition* condition); 171 virtual G4VParticleChange* PostStepDoIt(const G4Track& track, 172 const G4Step& step); 173 // -- AlongStep methods: 174 virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track& track, 175 G4double previousStepSize, 176 G4double currentMinimumStep, 177 G4double& proposedSafety, 178 G4GPILSelection* selection); 179 virtual G4VParticleChange* AlongStepDoIt(const G4Track& track, 180 const G4Step& step); 181 // -- AtRest methods 182 virtual G4double AtRestGetPhysicalInteractionLength(const G4Track&, 183 G4ForceCondition*); 184 virtual G4VParticleChange* AtRestDoIt(const G4Track&, const G4Step&); 185 186 virtual G4bool IsApplicable(const G4ParticleDefinition& pd); 187 virtual void BuildPhysicsTable(const G4ParticleDefinition& pd); 188 virtual void PreparePhysicsTable(const G4ParticleDefinition& pd); 189 virtual G4bool StorePhysicsTable(const G4ParticleDefinition* pd, 190 const G4String& s, G4bool f); 191 virtual G4bool RetrievePhysicsTable(const G4ParticleDefinition* pd, 192 const G4String& s, G4bool f); 193 // -- 194 virtual void SetProcessManager(const G4ProcessManager*); 195 virtual const G4ProcessManager* GetProcessManager(); 196 // -- 197 virtual void ResetNumberOfInteractionLengthLeft(); 198 199 virtual void SetMasterProcess(G4VProcess* masterP); 200 virtual void BuildWorkerPhysicsTable(const G4ParticleDefinition& pd); 201 virtual void PrepareWorkerPhysicsTable(const G4ParticleDefinition& pd); 202 203 private: 204 205 // ---- Internal utility methods: 206 void SetUpFirstLastFlags(); 207 void ResetForUnbiasedTracking(); 208 void ReorderBiasingVectorAsGPIL(); 209 210 G4int IdxFirstLast(G4int firstLast, G4int GPILDoIt, G4int physAll) const 211 { 212 // -- be careful : all arguments are *assumed* to be 0 or 1. No check 213 // -- for that is provided. Should be of pure internal usage. 214 return 4*firstLast + 2*GPILDoIt + physAll; 215 } 216 217 // -- method used to anticipate stepping manager calls to PostStepGPIL 218 // -- of wrapped processes : this method calls wrapped process PostStepGPIL 219 // -- and caches results for PostStepGPIL and condition. 220 void InvokeWrappedProcessPostStepGPIL( const G4Track& track, 221 G4double previousStepSize, 222 G4ForceCondition* condition ); 223 private: 224 225 G4Track* fCurrentTrack = nullptr; 226 G4double fPreviousStepSize = -1.0; 227 G4double fCurrentMinimumStep = -1.0; 228 G4double fProposedSafety = -1.0; 229 230 G4VBiasingOperation* fOccurenceBiasingOperation = nullptr; 231 G4VBiasingOperation* fFinalStateBiasingOperation = nullptr; 232 G4VBiasingOperation* fNonPhysicsBiasingOperation = nullptr; 233 G4VBiasingOperation* fPreviousOccurenceBiasingOperation = nullptr; 234 G4VBiasingOperation* fPreviousFinalStateBiasingOperation = nullptr; 235 G4VBiasingOperation* fPreviousNonPhysicsBiasingOperation = nullptr; 236 237 G4bool fResetWrappedProcessInteractionLength = false; 238 239 G4VProcess* fWrappedProcess = nullptr; 240 const G4bool fIsPhysicsBasedBiasing = false; 241 const G4bool fWrappedProcessIsAtRest = false; 242 const G4bool fWrappedProcessIsAlong = false; 243 const G4bool fWrappedProcessIsPost = false; 244 245 G4double fWrappedProcessPostStepGPIL = -1.0; 246 G4double fBiasingPostStepGPIL = -1.0; 247 G4double fWrappedProcessInteractionLength = -1.0; // -- inverse of analog cross-section 248 G4ForceCondition fWrappedProcessForceCondition = NotForced; 249 G4ForceCondition fBiasingForceCondition = NotForced; 250 G4double fWrappedProcessAlongStepGPIL = -1.0; 251 G4double fBiasingAlongStepGPIL = -1.0; 252 G4GPILSelection fWrappedProcessGPILSelection = NotCandidateForSelection; 253 G4GPILSelection fBiasingGPILSelection = NotCandidateForSelection; 254 255 const G4VBiasingInteractionLaw* fBiasingInteractionLaw = nullptr; 256 const G4VBiasingInteractionLaw* fPreviousBiasingInteractionLaw = nullptr; 257 G4InteractionLawPhysical* fPhysicalInteractionLaw = nullptr; 258 G4ParticleChangeForOccurenceBiasing* fOccurenceBiasingParticleChange = nullptr; 259 G4ParticleChangeForNothing* fDummyParticleChange = nullptr; 260 G4bool fFirstLastFlags[8]; 261 262 // -- the instance being "firstGPIL" does work shared by other instances: 263 G4bool fIamFirstGPIL = false; 264 265 // -- MUST be **thread local**: 266 static G4Cache<G4bool> fResetInteractionLaws; 267 static G4Cache<G4bool> fCommonStart; 268 static G4Cache<G4bool> fCommonEnd; 269 static G4Cache<G4bool> fDoCommonConfigure; 270 271 const G4ProcessManager* fProcessManager = nullptr; 272 273 // -- the data shared among processes attached to a same process manager: 274 G4BiasingProcessSharedData* fSharedData = nullptr; 275 }; 276 277 #endif 278