Geant4-11
G4BiasingProcessInterface.hh
Go to the documentation of this file.
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//
27//
28//--------------------------------------------------------------------
29//
30// G4BiasingProcessInterface
31//
32// Class Description:
33// A wrapper process making the interface between the tracking
34// and the G4VBiasingOperator objects attached to volumes.
35// If this process holds a physics process, it forwards
36// tracking calls to this process in volume where not biasing
37// occurs. In volumes with biasing (with a G4VBiasingOperator
38// attached) the process gets what to do messaging the biasing
39// operator :
40// - at the PostStepGPIL level, for getting an occurrence biasing
41// operation. If such an operation is returned to the process
42// this operation will be messaged at several places.
43// - at the PostStepDoIt level, to get a possible final state
44// biasing operation
45// If the process does not hold a physics process, it is meant
46// as handling "non physics" biasing operations: pure splitting
47// or pure kiling for example (ie not brem splitting).
48//
49//--------------------------------------------------------------------
50// Initial version Sep. 2013 M. Verderi
51// Use of "shared data" class Sep. 2014 M. Verderi
52
53
54#ifndef G4BiasingProcessInterface_h
55#define G4BiasingProcessInterface_h
56
57#include "globals.hh"
58#include "G4VProcess.hh"
59#include "G4Cache.hh"
61
68
70public:
71 // --------------------------------------------------------------------------------
72 // -- constructor for dealing with biasing options not affecting physics processes:
73 // --------------------------------------------------------------------------------
74 G4BiasingProcessInterface( G4String name = "biasWrapper(0)" );
75 // ---------------------------------------------------------------
76 // -- constructor to transform the behaviour of a physics process:
77 // ---------------------------------------------------------------
78 // -- wrappedProcess pointer MUST NOT be null.
80 G4bool wrappedIsAtRest, G4bool wrappedIsAlongStep, G4bool wrappedIsPostStep,
81 G4String useThisName = "");
83 // -- pointer of wrapped physics process:
85
86 // ---------------------
87 // -- Biasing interface:
88 // ---------------------
89 // -- Helper methods:
90 // ------------------
91 // -- Current step and previous step biasing operator, if any:
92 G4VBiasingOperator* GetCurrentBiasingOperator() const {return fSharedData-> fCurrentBiasingOperator; }
94 // -- current and previous operation:
101
102 // -- Lists of processes cooperating under a same particle type/G4ProcessManager.
103 // -- The vector ordering is:
104 // -- - random, before first "run/beamOn"
105 // -- - that of the PostStepGetPhysicalInteractionLength() once "/run/beamOn" has been issued
106 const std::vector< const G4BiasingProcessInterface* >& GetBiasingProcessInterfaces() const
107 {return fSharedData-> fPublicBiasingProcessInterfaces;}
108 const std::vector< const G4BiasingProcessInterface* >& GetPhysicsBiasingProcessInterfaces() const
109 {return fSharedData-> fPublicPhysicsBiasingProcessInterfaces;}
110 const std::vector< const G4BiasingProcessInterface* >& GetNonPhysicsBiasingProcessInterfaces() const
112
113 // -- Get shared data for this process:
115 // -- Get shared data associated to a G4ProcessManager:
117
118
119
120 // ------------------
121 // -- Helper methods:
122 // ------------------
123 // ---- Tell is this process is first/last in the PostStep GPIL and DoIt lists.
124 // ---- If physOnly is true, only wrapper for physics processes are considered,
125 // ---- otherwise all G4BiasingProcessInterface processes of this particle are
126 // ---- considered.
127 // ---- These methods just return the corresponding flag values setup at
128 // ---- initialization phase by the next four ones.
129 // ---- Will not be updated if processes are activate/unactivated on the fly.
130 // ---- Use next methods (less fast) instead in this case.
131 G4bool GetIsFirstPostStepGPILInterface(G4bool physOnly = true) const;
132 G4bool GetIsLastPostStepGPILInterface(G4bool physOnly = true) const;
133 G4bool GetIsFirstPostStepDoItInterface(G4bool physOnly = true) const;
134 G4bool GetIsLastPostStepDoItInterface(G4bool physOnly = true) const;
135 // ---- Determine if the process is first/last in the PostStep GPIL and DoIt lists.
136 G4bool IsFirstPostStepGPILInterface(G4bool physOnly = true) const;
137 G4bool IsLastPostStepGPILInterface(G4bool physOnly = true) const;
138 G4bool IsFirstPostStepDoItInterface(G4bool physOnly = true) const;
139 G4bool IsLastPostStepDoItInterface(G4bool physOnly = true) const;
140 // -- Information about wrapped process:
144
145
146 // -- Information methods:
151 // -- return the actual PostStep and AlongStep limits returned by the process to the tracking :
154
155
156 // --------------------------------------------------------------
157 // -- G4VProcess interface --
158 // --------------------------------------------------------------
159public:
160 // -- Start/End tracking:
161 void StartTracking(G4Track* track);
162 void EndTracking();
163
164 // -- PostStep methods:
166 G4double previousStepSize,
168 virtual G4VParticleChange* PostStepDoIt(const G4Track& track,
169 const G4Step& step);
170 // -- AlongStep methods:
172 G4double previousStepSize,
173 G4double currentMinimumStep,
174 G4double& proposedSafety,
175 G4GPILSelection* selection);
176 virtual G4VParticleChange* AlongStepDoIt(const G4Track& track,
177 const G4Step& step);
178 // -- AtRest methods
181 virtual G4VParticleChange* AtRestDoIt(const G4Track&,
182 const G4Step&);
183
184 virtual G4bool IsApplicable(const G4ParticleDefinition& pd);
185 virtual void BuildPhysicsTable(const G4ParticleDefinition& pd);
186 virtual void PreparePhysicsTable(const G4ParticleDefinition& pd);
188 const G4String& s, G4bool f);
190 const G4String& s, G4bool f);
191 // --
192 virtual void SetProcessManager(const G4ProcessManager*);
193 virtual const G4ProcessManager* GetProcessManager();
194 // --
196 // virtual void ClearNumberOfInteractionLengthLeft();
197 // --
198 // virtual void DumpInfo() const;
199
200 virtual void SetMasterProcess(G4VProcess* masterP);
201 virtual void BuildWorkerPhysicsTable(const G4ParticleDefinition& pd);
202 virtual void PrepareWorkerPhysicsTable(const G4ParticleDefinition& pd);
203
204
205private:
206 // ---- Internal utility methods:
207 void SetUpFirstLastFlags();
210
215
222
224
230
231
234 G4double fWrappedProcessInteractionLength; // -- inverse of analog cross-section
241
248 G4int IdxFirstLast(G4int firstLast, G4int GPILDoIt, G4int physAll) const
249 {
250 // -- be careful : all arguments are *assumed* to be 0 or 1. No check
251 // -- for that is provided. Should be of pure internal usage.
252 return 4*firstLast + 2*GPILDoIt + physAll;
253 }
254 // -- method used to anticipate stepping manager calls to PostStepGPIL
255 // -- of wrapped processes : this method calls wrapped process PostStepGPIL
256 // -- and caches results for PostStepGPIL and condition.
257 void InvokeWrappedProcessPostStepGPIL( const G4Track& track,
258 G4double previousStepSize,
260 // -- the instance being "firstGPIL" does work shared by other instances:
262
263
264 // -- MUST be **thread local**:
269
271
272
273 // -- the data shared among processes attached to a same process manager:
275
276};
277
278#endif
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
G4GPILSelection
static constexpr double s
Definition: G4SIunits.hh:154
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const std::vector< const G4BiasingProcessInterface * > & GetNonPhysicsBiasingProcessInterfaces() const
static G4Cache< G4bool > fCommonStart
G4bool GetIsFirstPostStepGPILInterface(G4bool physOnly=true) const
const std::vector< const G4BiasingProcessInterface * > & GetPhysicsBiasingProcessInterfaces() const
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
G4bool IsLastPostStepGPILInterface(G4bool physOnly=true) const
static G4Cache< G4bool > fResetInteractionLaws
G4BiasingProcessSharedData * fSharedData
const std::vector< const G4BiasingProcessInterface * > & GetBiasingProcessInterfaces() const
virtual void SetProcessManager(const G4ProcessManager *)
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
virtual void SetMasterProcess(G4VProcess *masterP)
const G4VBiasingInteractionLaw * fPreviousBiasingInteractionLaw
static G4Cache< G4bool > fDoCommonConfigure
void InvokeWrappedProcessPostStepGPIL(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
virtual G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &step)
G4VBiasingOperation * GetPreviousFinalStateBiasingOperation() const
const G4ProcessManager * fProcessManager
G4bool GetIsLastPostStepDoItInterface(G4bool physOnly=true) const
virtual void PrepareWorkerPhysicsTable(const G4ParticleDefinition &pd)
G4int IdxFirstLast(G4int firstLast, G4int GPILDoIt, G4int physAll) const
G4bool GetIsLastPostStepGPILInterface(G4bool physOnly=true) const
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &step)
G4VBiasingOperator * GetPreviousBiasingOperator() const
G4ParticleChangeForNothing * fDummyParticleChange
G4VBiasingOperation * GetCurrentFinalStateBiasingOperation() const
G4VBiasingOperation * fPreviousOccurenceBiasingOperation
G4bool GetIsFirstPostStepDoItInterface(G4bool physOnly=true) const
G4bool IsFirstPostStepDoItInterface(G4bool physOnly=true) const
virtual G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
static G4Cache< G4bool > fCommonEnd
G4bool IsLastPostStepDoItInterface(G4bool physOnly=true) const
G4VBiasingOperator * GetCurrentBiasingOperator() const
G4VBiasingOperation * fFinalStateBiasingOperation
G4BiasingProcessInterface(G4String name="biasWrapper(0)")
G4VBiasingOperation * GetCurrentNonPhysicsBiasingOperation() const
const G4BiasingProcessSharedData * GetSharedData() const
virtual G4bool RetrievePhysicsTable(const G4ParticleDefinition *pd, const G4String &s, G4bool f)
virtual const G4ProcessManager * GetProcessManager()
virtual void PreparePhysicsTable(const G4ParticleDefinition &pd)
G4VBiasingOperation * fNonPhysicsBiasingOperation
G4VBiasingOperation * fPreviousNonPhysicsBiasingOperation
G4ParticleChangeForOccurenceBiasing * fOccurenceBiasingParticleChange
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
virtual void BuildPhysicsTable(const G4ParticleDefinition &pd)
G4VBiasingOperation * GetPreviousOccurenceBiasingOperation() const
virtual G4bool StorePhysicsTable(const G4ParticleDefinition *pd, const G4String &s, G4bool f)
G4InteractionLawPhysical * fPhysicalInteractionLaw
const G4VBiasingInteractionLaw * fBiasingInteractionLaw
virtual G4bool IsApplicable(const G4ParticleDefinition &pd)
G4bool IsFirstPostStepGPILInterface(G4bool physOnly=true) const
G4VBiasingOperation * fOccurenceBiasingOperation
G4VBiasingOperation * GetCurrentOccurenceBiasingOperation() const
G4VProcess * GetWrappedProcess() const
G4VBiasingOperation * GetPreviousNonPhysicsBiasingOperation() const
virtual void BuildWorkerPhysicsTable(const G4ParticleDefinition &pd)
G4VBiasingOperation * fPreviousFinalStateBiasingOperation
G4VBiasingOperator * fPreviousBiasingOperator
std::vector< const G4BiasingProcessInterface * > fPublicNonPhysicsBiasingProcessInterfaces
Definition: G4Step.hh:62
const char * name(G4int ptype)