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