Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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;
69 
71 public:
72  // --------------------------------------------------------------------------------
73  // -- constructor for dealing with biasing options not affecting physics processes:
74  // --------------------------------------------------------------------------------
75  G4BiasingProcessInterface( G4String name = "biasWrapper(0)" );
76  // ---------------------------------------------------------------
77  // -- constructor to transform the behaviour of a physics process:
78  // ---------------------------------------------------------------
79  // -- wrappedProcess pointer MUST NOT be null.
80  G4BiasingProcessInterface(G4VProcess* wrappedProcess,
81  G4bool wrappedIsAtRest, G4bool wrappedIsAlongStep, G4bool wrappedIsPostStep,
82  G4String useThisName = "");
84  // -- pointer of wrapped physics process:
85  G4VProcess* GetWrappedProcess() const {return fWrappedProcess;}
86 
87  // ---------------------
88  // -- Biasing interface:
89  // ---------------------
90  // -- Helper methods:
91  // ------------------
92  // -- Current step and previous step biasing operator, if any:
93  G4VBiasingOperator* GetCurrentBiasingOperator() const {return fSharedData-> fCurrentBiasingOperator; }
94  G4VBiasingOperator* GetPreviousBiasingOperator() const {return fSharedData->fPreviousBiasingOperator; }
95  // -- current and previous operation:
96  G4VBiasingOperation* GetCurrentNonPhysicsBiasingOperation() const { return fNonPhysicsBiasingOperation; }
97  G4VBiasingOperation* GetPreviousNonPhysicsBiasingOperation() const { return fPreviousNonPhysicsBiasingOperation; }
98  G4VBiasingOperation* GetCurrentOccurenceBiasingOperation() const { return fOccurenceBiasingOperation; }
99  G4VBiasingOperation* GetPreviousOccurenceBiasingOperation() const { return fPreviousOccurenceBiasingOperation; }
100  G4VBiasingOperation* GetCurrentFinalStateBiasingOperation() const { return fFinalStateBiasingOperation; }
101  G4VBiasingOperation* GetPreviousFinalStateBiasingOperation() const { return fPreviousFinalStateBiasingOperation; }
102 
103  // -- Lists of processes cooperating under a same particle type/G4ProcessManager.
104  // -- The vector ordering is:
105  // -- - random, before first "run/beamOn"
106  // -- - that of the PostStepGetPhysicalInteractionLength() once "/run/beamOn" has been issued
107  const std::vector< const G4BiasingProcessInterface* >& GetBiasingProcessInterfaces() const
108  {return fSharedData-> fPublicBiasingProcessInterfaces;}
109  const std::vector< const G4BiasingProcessInterface* >& GetPhysicsBiasingProcessInterfaces() const
110  {return fSharedData-> fPublicPhysicsBiasingProcessInterfaces;}
111  const std::vector< const G4BiasingProcessInterface* >& GetNonPhysicsBiasingProcessInterfaces() const
112  {return fSharedData->fPublicNonPhysicsBiasingProcessInterfaces;}
113 
114  // -- Get shared data for this process:
115  const G4BiasingProcessSharedData* GetSharedData() const { return fSharedData; }
116  // -- Get shared data associated to a G4ProcessManager:
118 
119 
120 
121  // ------------------
122  // -- Helper methods:
123  // ------------------
124  // ---- Tell is this process is first/last in the PostStep GPIL and DoIt lists.
125  // ---- If physOnly is true, only wrapper for physics processes are considered,
126  // ---- otherwise all G4BiasingProcessInterface processes of this particle are
127  // ---- considered.
128  // ---- These methods just return the corresponding flag values setup at
129  // ---- initialization phase by the next four ones.
130  // ---- Will not be updated if processes are activate/unactivated on the fly.
131  // ---- Use next methods (less fast) instead in this case.
132  G4bool GetIsFirstPostStepGPILInterface(G4bool physOnly = true) const;
133  G4bool GetIsLastPostStepGPILInterface(G4bool physOnly = true) const;
134  G4bool GetIsFirstPostStepDoItInterface(G4bool physOnly = true) const;
135  G4bool GetIsLastPostStepDoItInterface(G4bool physOnly = true) const;
136  // ---- Determine if the process is first/last in the PostStep GPIL and DoIt lists.
137  G4bool IsFirstPostStepGPILInterface(G4bool physOnly = true) const;
138  G4bool IsLastPostStepGPILInterface(G4bool physOnly = true) const;
139  G4bool IsFirstPostStepDoItInterface(G4bool physOnly = true) const;
140  G4bool IsLastPostStepDoItInterface(G4bool physOnly = true) const;
141  // -- Information about wrapped process:
142  G4bool GetWrappedProcessIsAtRest() const { return fWrappedProcessIsAtRest; }
143  G4bool GetWrappedProcessIsAlong() const { return fWrappedProcessIsAlong; }
144  G4bool GetWrappedProcessIsPost() const { return fWrappedProcessIsPost; }
145 
146 
147  // -- Information methods:
148  G4double GetPreviousStepSize() const { return fPreviousStepSize;}
149  G4double GetCurrentMinimumStep() const { return fCurrentMinimumStep;}
150  G4double GetProposedSafety() const { return fProposedSafety;}
151  void SetProposedSafety(G4double sft) { fProposedSafety = sft;}
152  // -- return the actual PostStep and AlongStep limits returned by the process to the tracking :
153  G4double GetPostStepGPIL() const { return fBiasingPostStepGPIL; }
154  G4double GetAlongStepGPIL() const { return fWrappedProcessAlongStepGPIL; }
155 
156 
157  // --------------------------------------------------------------
158  // -- G4VProcess interface --
159  // --------------------------------------------------------------
160 public:
161  // -- Start/End tracking:
162  void StartTracking(G4Track* track);
163  void EndTracking();
164 
165  // -- PostStep methods:
167  G4double previousStepSize,
169  virtual G4VParticleChange* PostStepDoIt(const G4Track& track,
170  const G4Step& step);
171  // -- AlongStep methods:
173  G4double previousStepSize,
174  G4double currentMinimumStep,
175  G4double& proposedSafety,
176  G4GPILSelection* selection);
177  virtual G4VParticleChange* AlongStepDoIt(const G4Track& track,
178  const G4Step& step);
179  // -- AtRest methods
182  virtual G4VParticleChange* AtRestDoIt(const G4Track&,
183  const G4Step&);
184 
185  virtual G4bool IsApplicable(const G4ParticleDefinition& pd);
186  virtual void BuildPhysicsTable(const G4ParticleDefinition& pd);
187  virtual void PreparePhysicsTable(const G4ParticleDefinition& pd);
189  const G4String& s, G4bool f);
191  const G4String& s, G4bool f);
192  // --
193  virtual void SetProcessManager(const G4ProcessManager*);
194  virtual const G4ProcessManager* GetProcessManager();
195  // --
196  virtual void ResetNumberOfInteractionLengthLeft();
197  // virtual void ClearNumberOfInteractionLengthLeft();
198  // --
199  // virtual void DumpInfo() const;
200 
201  virtual void SetMasterProcess(G4VProcess* masterP);
202  virtual void BuildWorkerPhysicsTable(const G4ParticleDefinition& pd);
203  virtual void PrepareWorkerPhysicsTable(const G4ParticleDefinition& pd);
204 
205 
206 private:
207  // ---- Internal utility methods:
208  void SetUpFirstLastFlags();
209  void ResetForUnbiasedTracking();
210  void ReorderBiasingVectorAsGPIL();
211 
212  G4Track* fCurrentTrack;
213  G4double fPreviousStepSize;
214  G4double fCurrentMinimumStep;
215  G4double fProposedSafety;
216 
217  G4VBiasingOperation* fOccurenceBiasingOperation;
218  G4VBiasingOperation* fFinalStateBiasingOperation;
219  G4VBiasingOperation* fNonPhysicsBiasingOperation;
220  G4VBiasingOperation* fPreviousOccurenceBiasingOperation;
221  G4VBiasingOperation* fPreviousFinalStateBiasingOperation;
222  G4VBiasingOperation* fPreviousNonPhysicsBiasingOperation;
223 
224  G4bool fResetWrappedProcessInteractionLength;
225 
226  G4VProcess* fWrappedProcess;
227  const G4bool fIsPhysicsBasedBiasing;
228  const G4bool fWrappedProcessIsAtRest;
229  const G4bool fWrappedProcessIsAlong;
230  const G4bool fWrappedProcessIsPost;
231 
232 
233  G4double fWrappedProcessPostStepGPIL;
234  G4double fBiasingPostStepGPIL;
235  G4double fWrappedProcessInteractionLength; // -- inverse of analog cross-section
236  G4ForceCondition fWrappedProcessForceCondition;
237  G4ForceCondition fBiasingForceCondition;
238  G4double fWrappedProcessAlongStepGPIL;
239  G4double fBiasingAlongStepGPIL;
240  G4GPILSelection fWrappedProcessGPILSelection;
241  G4GPILSelection fBiasingGPILSelection;
242 
243  const G4VBiasingInteractionLaw* fBiasingInteractionLaw;
244  const G4VBiasingInteractionLaw* fPreviousBiasingInteractionLaw;
245  G4InteractionLawPhysical* fPhysicalInteractionLaw;
246  G4ParticleChangeForOccurenceBiasing* fOccurenceBiasingParticleChange;
247  G4ParticleChangeForNothing* fDummyParticleChange;
248  G4bool fFirstLastFlags[8];
249  G4int IdxFirstLast(G4int firstLast, G4int GPILDoIt, G4int physAll) const
250  {
251  // -- be careful : all arguments are *assumed* to be 0 or 1. No check
252  // -- for that is provided. Should be of pure internal usage.
253  return 4*firstLast + 2*GPILDoIt + physAll;
254  }
255  // -- method used to anticipate stepping manager calls to PostStepGPIL
256  // -- of wrapped processes : this method calls wrapped process PostStepGPIL
257  // -- and caches results for PostStepGPIL and condition.
258  void InvokeWrappedProcessPostStepGPIL( const G4Track& track,
259  G4double previousStepSize,
261  // -- the instance being "firstGPIL" does work shared by other instances:
262  G4bool fIamFirstGPIL;
263 
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;
272 
273 
274  // -- the data shared among processes attached to a same process manager:
275  G4BiasingProcessSharedData* fSharedData;
276 
277 };
278 
279 #endif
G4double condition(const G4ErrorSymMatrix &m)
G4VBiasingOperator * GetCurrentBiasingOperator() const
const XML_Char * name
Definition: expat.h:151
virtual void BuildWorkerPhysicsTable(const G4ParticleDefinition &pd)
G4bool IsLastPostStepGPILInterface(G4bool physOnly=true) const
G4VBiasingOperation * GetCurrentOccurenceBiasingOperation() const
virtual void PreparePhysicsTable(const G4ParticleDefinition &pd)
virtual G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
virtual G4bool StorePhysicsTable(const G4ParticleDefinition *pd, const G4String &s, G4bool f)
virtual void BuildPhysicsTable(const G4ParticleDefinition &pd)
virtual void PrepareWorkerPhysicsTable(const G4ParticleDefinition &pd)
virtual G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &step)
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &step)
G4VBiasingOperation * GetPreviousNonPhysicsBiasingOperation() const
virtual void SetProcessManager(const G4ProcessManager *)
int G4int
Definition: G4Types.hh:78
G4VProcess * GetWrappedProcess() const
G4bool GetIsFirstPostStepGPILInterface(G4bool physOnly=true) const
const XML_Char * s
Definition: expat.h:262
G4bool GetIsLastPostStepDoItInterface(G4bool physOnly=true) const
G4bool GetIsLastPostStepGPILInterface(G4bool physOnly=true) const
bool G4bool
Definition: G4Types.hh:79
G4bool GetIsFirstPostStepDoItInterface(G4bool physOnly=true) const
const std::vector< const G4BiasingProcessInterface * > & GetBiasingProcessInterfaces() const
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
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
double G4double
Definition: G4Types.hh:76
G4ForceCondition
G4VBiasingOperation * GetCurrentFinalStateBiasingOperation() const
G4bool IsLastPostStepDoItInterface(G4bool physOnly=true) const
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