Geant4.10
 All Data Structures 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 
53 #ifndef G4BiasingProcessInterface_h
54 #define G4BiasingProcessInterface_h
55 
56 #include "globals.hh"
57 #include "G4VProcess.hh"
58 #include "G4Cache.hh"
59 
62 class G4VBiasingOperator;
65 class G4ParticleChange;
67 
69 public:
70  // --------------------------------------------------------------------------------
71  // -- constructor for dealing with biasing options not affecting physics processes:
72  // --------------------------------------------------------------------------------
73  G4BiasingProcessInterface( G4String name = "biasWrapper(0)" );
74  // ---------------------------------------------------------------
75  // -- constructor to transform the behaviour of a physics process:
76  // ---------------------------------------------------------------
77  // -- wrappedProcess pointer MUST NOT be null.
78  G4BiasingProcessInterface(G4VProcess* wrappedProcess,
79  G4bool wrappedIsAtRest, G4bool wrappedIsAlongStep, G4bool wrappedIsPostStep,
80  G4String useThisName = "");
82  // -- pointer of wrapped physics process:
83  G4VProcess* GetWrappedProcess() const {return fWrappedProcess;}
84 
85  // ---------------------
86  // -- Biasing interface:
87  // ---------------------
88  // -- Helper methods:
89  // ------------------
90  // -- Current step and previous step biasing operator, if any:
91  G4VBiasingOperator* GetCurrentBiasingOperator() const {return fCurrentBiasingOperator; }
92  G4VBiasingOperator* GetPreviousBiasingOperator() const {return fPreviousBiasingOperator; }
93  // -- current and previous operation:
94  G4VBiasingOperation* GetCurrentOccurenceBiasingOperation() const { return fOccurenceBiasingOperation; }
95  G4VBiasingOperation* GetPreviousOccurenceBiasingOperation() const { return fPreviousOccurenceBiasingOperation; }
96  G4VBiasingOperation* GetCurrentFinalStateBiasingOperation() const { return fFinalStateBiasingOperation; }
97  G4VBiasingOperation* GetPreviousFinalStateBiasingOperation() const { return fPreviousFinalStateBiasingOperation; }
98  G4VBiasingOperation* GetCurrentNonPhysicsBiasingOperation() const { return fNonPhysicsBiasingOperation; }
99  G4VBiasingOperation* GetPreviousNonPhysicsBiasingOperation() const { return fPreviousNonPhysicsBiasingOperation; }
100 
101 
102 
103  // ------------------
104  // -- Helper methods:
105  // ------------------
106  // ---- Tell is this process is first/last in the PostStep GPIL and DoIt lists.
107  // ---- If physOnly is true, only wrapper for physics processes are considered,
108  // ---- otherwise all G4BiasingProcessInterface processes of this particle are
109  // ---- considered.
110  // ---- These methods just return the corresponding flag values setup at
111  // ---- initialization phase by the next four ones.
112  // ---- Will not be updated if processes are activate/unactivated on the fly.
113  // ---- Use next methods (less fast) instead in this case.
114  G4bool GetIsFirstPostStepGPILInterface(G4bool physOnly = true) const;
115  G4bool GetIsLastPostStepGPILInterface(G4bool physOnly = true) const;
116  G4bool GetIsFirstPostStepDoItInterface(G4bool physOnly = true) const;
117  G4bool GetIsLastPostStepDoItInterface(G4bool physOnly = true) const;
118  // ---- Determine if the process is first/last in the PostStep GPIL and DoIt lists.
119  G4bool IsFirstPostStepGPILInterface(G4bool physOnly = true) const;
120  G4bool IsLastPostStepGPILInterface(G4bool physOnly = true) const;
121  G4bool IsFirstPostStepDoItInterface(G4bool physOnly = true) const;
122  G4bool IsLastPostStepDoItInterface(G4bool physOnly = true) const;
123  // -- Information about wrapped process:
124  G4bool GetWrappedProcessIsAtRest() const { return fWrappedProcessIsAtRest; }
125  G4bool GetWrappedProcessIsAlong() const { return fWrappedProcessIsAlong; }
126  G4bool GetWrappedProcessIsPost() const { return fWrappedProcessIsPost; }
127 
128 
129  // -- Information methods:
130  G4double GetPreviousStepSize() const { return fPreviousStepSize;}
131  G4double GetCurrentMinimumStep() const { return fCurrentMinimumStep;}
132  G4double GetProposedSafety() const { return fProposedSafety;}
133  void SetProposedSafety(G4double sft) { fProposedSafety = sft;}
134  // -- return the actual PostStep and AlongStep limits returned by the process to the tracking :
135  G4double GetPostStepGPIL() const { return fBiasingPostStepGPIL; }
136  G4double GetAlongStepGPIL() const { return fWrappedProcessAlongStepGPIL; }
137 
138 
139  // --------------------------------------------------------------
140  // -- G4VProcess interface --
141  // --------------------------------------------------------------
142 public:
143  // -- Start/End tracking:
144  void StartTracking(G4Track* track);
145  void EndTracking();
146 
147  // -- PostStep methods:
149  G4double previousStepSize,
151  virtual G4VParticleChange* PostStepDoIt(const G4Track& track,
152  const G4Step& step);
153  // -- AlongStep methods:
155  G4double previousStepSize,
156  G4double currentMinimumStep,
157  G4double& proposedSafety,
158  G4GPILSelection* selection);
159  virtual G4VParticleChange* AlongStepDoIt(const G4Track& track,
160  const G4Step& step);
161  // -- AtRest methods
164  virtual G4VParticleChange* AtRestDoIt(const G4Track&,
165  const G4Step&);
166 
167  virtual G4bool IsApplicable(const G4ParticleDefinition& pd);
168  virtual void BuildPhysicsTable(const G4ParticleDefinition& pd);
169  virtual void PreparePhysicsTable(const G4ParticleDefinition& pd);
171  const G4String& s, G4bool f);
173  const G4String& s, G4bool f);
174  // --
175  virtual void SetProcessManager(const G4ProcessManager*);
176  virtual const G4ProcessManager* GetProcessManager();
177  // --
178  virtual void ResetNumberOfInteractionLengthLeft();
179  // virtual void ClearNumberOfInteractionLengthLeft();
180  // --
181  // virtual void DumpInfo() const;
182 
183  virtual void SetMasterProcess(G4VProcess* masterP);
184  virtual void BuildWorkerPhysicsTable(const G4ParticleDefinition& pd);
185  virtual void PrepareWorkerPhysicsTable(const G4ParticleDefinition& pd);
186 
187 
188 private:
189  // ---- Setup first/last flags:
190  void SetUpFirstLastFlags();
191  void ResetForUnbiasedTracking();
192 
193  G4Track* fCurrentTrack;
194  G4double fPreviousStepSize;
195  G4double fCurrentMinimumStep;
196  G4double fProposedSafety;
197 
198  G4VBiasingOperator* fCurrentBiasingOperator;
199  G4VBiasingOperator* fPreviousBiasingOperator;
200  G4VBiasingOperation* fOccurenceBiasingOperation;
201  G4VBiasingOperation* fFinalStateBiasingOperation;
202  G4VBiasingOperation* fNonPhysicsBiasingOperation;
203  G4VBiasingOperation* fPreviousOccurenceBiasingOperation;
204  G4VBiasingOperation* fPreviousFinalStateBiasingOperation;
205  G4VBiasingOperation* fPreviousNonPhysicsBiasingOperation;
206 
207  G4bool fResetWrappedProcessInteractionLength;
208 
209  G4VProcess* fWrappedProcess;
210  const G4bool fIsPhysicsBasedBiasing;
211  const G4bool fWrappedProcessIsAtRest;
212  const G4bool fWrappedProcessIsAlong;
213  const G4bool fWrappedProcessIsPost;
214 
215 
216  G4double fWrappedProcessPostStepGPIL;
217  G4double fBiasingPostStepGPIL;
218  G4double fWrappedProcessInteractionLength; // -- inverse of analog cross-section
219  G4ForceCondition fWrappedProcessForceCondition;
220  G4ForceCondition fBiasingForceCondition;
221  G4double fWrappedProcessAlongStepGPIL;
222  G4double fBiasingAlongStepGPIL;
223  G4GPILSelection fWrappedProcessGPILSelection;
224  G4GPILSelection fBiasingGPILSelection;
225 
226  const G4VBiasingInteractionLaw* fBiasingInteractionLaw;
227  const G4VBiasingInteractionLaw* fPreviousBiasingInteractionLaw;
228  G4InteractionLawPhysical* fPhysicalInteractionLaw;
229  G4ParticleChangeForOccurenceBiasing* fOccurenceBiasingParticleChange;
230  G4ParticleChange* fParticleChange; // -- to be changed with a light version for weight only
231  G4ParticleChangeForNothing* fDummyParticleChange;
232  G4bool fFirstLastFlags[8];
233  G4int IdxFirstLast(G4int firstLast, G4int GPILDoIt, G4int physAll) const
234  {
235  // -- be careful : all arguments are *assumed* to be 0 or 1. No check
236  // -- for that is provided. Should be of pure internal usage.
237  return 4*firstLast + 2*GPILDoIt + physAll;
238  }
239  // -- TO DO : let some common operation be done only once, by the first
240  // -- process in GPIL. Nothing done yet, may help to reduce overhead.
241  // -- For example each process instance fetch for the current operator,
242  // -- which could be done once, and with pointer shared.
243  // -- Idea is to have a class containing the data members shared among
244  // -- the instances.
245  // -- Could make these shared data members shared only per particle type
246  // -- (at present, with a single static, any instance of a process shares
247  // -- its members with any other instance, regardless of particle type.)
248  G4bool fIamFirstGPIL;
249 
250 
251  // -- MUST be **thread local**:
252  static G4Cache<G4bool> fResetInteractionLaws;
253  static G4Cache<G4bool> fCommonStart;
254  static G4Cache<G4bool> fCommonEnd;
255 
256  // -- Maintain lists of interfaces attached to a same particle (ie
257  // -- to a same G4ProcessManager ).
258  // -- This is used when biasing operations in different processes need to
259  // -- cooperate (like for knowing the total cross-section for example).
260  std::vector< G4BiasingProcessInterface* >* fCoInterfaces;
261  // -- thread local:
262  static G4MapCache< const G4ProcessManager*,
263  std::vector< G4BiasingProcessInterface* > > fManagerInterfaceMap;
264  const G4ProcessManager* fProcessManager;
265 
266 };
267 
268 #endif
G4double condition(const G4ErrorSymMatrix &m)
G4VBiasingOperator * GetCurrentBiasingOperator() const
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)
const XML_Char * s
virtual G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &step)
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &step)
const XML_Char * name
G4VBiasingOperation * GetPreviousNonPhysicsBiasingOperation() const
virtual void SetProcessManager(const G4ProcessManager *)
int G4int
Definition: G4Types.hh:78
G4VProcess * GetWrappedProcess() const
G4bool GetIsFirstPostStepGPILInterface(G4bool physOnly=true) const
G4bool GetIsLastPostStepDoItInterface(G4bool physOnly=true) const
G4bool GetIsLastPostStepGPILInterface(G4bool physOnly=true) const
bool G4bool
Definition: G4Types.hh:79
G4bool GetIsFirstPostStepDoItInterface(G4bool physOnly=true) 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()
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
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)