Geant4  10.02.p01
G4FastSimulationManagerProcess.cc
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: G4FastSimulationManagerProcess.cc 68056 2013-03-13 14:44:48Z gcosmo $
28 //
29 //
30 //---------------------------------------------------------------
31 //
32 // G4FastSimulationProcess.cc
33 //
34 // Description:
35 // The process that triggers the parameterised simulations,
36 // if any.
37 //
38 // History:
39 // August 97: First implementation. Verderi && MoraDeFreitas.
40 // October 06: move to parallel geometry scheme, M. Verderi
41 //---------------------------------------------------------------
42 
43 #include "G4ios.hh"
47 #include "G4PathFinder.hh"
48 #include "G4ParticleChange.hh"
49 #include "G4FieldTrackUpdator.hh"
50 
51 #define PARANOIA
52 
55  G4ProcessType theType) :
56  G4VProcess(processName,theType),
57  fWorldVolume(0),
58  fIsTrackingTime(false),
59  fGhostNavigator(0),
60  fGhostNavigatorIndex(-1),
61  fIsGhostGeometry(false),
62  fGhostSafety(-1.0),
63  fFieldTrack('0'),
64  fFastSimulationManager(0),
65  fFastSimulationTrigger(false)
66 {
67  // -- set Process Sub Type
68  SetProcessSubType(static_cast<int>(FASTSIM_ManagerProcess));
69 
70 
73 
75  if (verboseLevel>0) G4cout << "G4FastSimulationManagerProcess `" << GetProcessName()
76  << "' is created, and will message geometry with world volume `"
77  << fWorldVolume->GetName() << "'." << G4endl;
79 }
80 
81 
84  const G4String& worldVolumeName,
85  G4ProcessType theType) :
86  G4VProcess(processName,theType),
87  fWorldVolume(0),
88  fIsTrackingTime(false),
89  fGhostNavigator(0),
90  fGhostNavigatorIndex(-1),
91  fIsGhostGeometry(false),
92  fGhostSafety(-1.0),
93  fFieldTrack('0'),
94  fFastSimulationManager(0),
95  fFastSimulationTrigger(false)
96 {
97  // -- set Process Sub Type
98  SetProcessSubType(static_cast<int>(FASTSIM_ManagerProcess));
99 
100 
103 
104  SetWorldVolume(worldVolumeName);
105  if (verboseLevel>0) G4cout << "G4FastSimulationManagerProcess `" << GetProcessName()
106  << "' is created, and will message geometry with world volume `"
107  << fWorldVolume->GetName() << "'." << G4endl;
109 }
110 
111 
114  G4VPhysicalVolume* worldVolume,
115  G4ProcessType theType) :
116  G4VProcess(processName,theType),
117  fWorldVolume(0),
118  fIsTrackingTime(false),
119  fGhostNavigator(0),
120  fGhostNavigatorIndex(-1),
121  fIsGhostGeometry(false),
122  fGhostSafety(-1.0),
123  fFieldTrack('0'),
124  fFastSimulationManager(0),
125  fFastSimulationTrigger(false)
126 {
127  // -- set Process Sub Type
128  SetProcessSubType(static_cast<int>(FASTSIM_ManagerProcess));
129 
130 
133 
134  SetWorldVolume(worldVolume);
135  if (verboseLevel>0) G4cout << "G4FastSimulationManagerProcess `" << GetProcessName()
136  << "' is created, and will message geometry with world volume `"
137  << fWorldVolume->GetName() << "'." << G4endl;
139 }
140 
141 
143 {
145 }
146 
147 
148 // -----------------------
149 // User access methods:
150 // -----------------------
152 {
153  if (fIsTrackingTime)
154  {
156  ed << "G4FastSimulationManagerProcess `" << GetProcessName()
157  << "': changing of world volume at tracking time is not allowed." << G4endl;
158  G4Exception("G4FastSimulationManagerProcess::SetWorldVolume(const G4String)",
159  "FastSim002",
160  JustWarning, ed,
161  "Call ignored.");
162  }
163  else
164  {
165  G4VPhysicalVolume* newWorld = fTransportationManager->IsWorldExisting(newWorldName);
166  if (newWorld == 0)
167  {
168  G4ExceptionDescription tellWhatIsWrong;
169  tellWhatIsWrong << "Volume newWorldName = `" << newWorldName
170  << "' is not a parallel world nor the mass world volume."
171  << G4endl;
172  G4Exception("G4FastSimulationManagerProcess::SetWorldVolume(const G4String)",
173  "FastSim003",
175  tellWhatIsWrong);
176  }
177  if (verboseLevel>0)
178  {
179  if (fWorldVolume) G4cout << "G4FastSimulationManagerProcess `" << GetProcessName()
180  << "': changing world volume from '" << fWorldVolume->GetName()
181  << "' to `" << newWorld << "'." << G4endl;
182  else G4cout << "G4FastSimulationManagerProcess `" << GetProcessName()
183  << "': setting world volume from to `"<< newWorld->GetName() << "'." << G4endl;
184  }
185  fWorldVolume = newWorld;
186  }
187 }
188 
189 
191 {
192  if (newWorld) SetWorldVolume(newWorld->GetName());
193  else
194  {
195  G4ExceptionDescription tellWhatIsWrong;
196  tellWhatIsWrong << "Null pointer passed for world volume." << G4endl;
197  G4Exception("G4FastSimulationManagerProcess::SetWorldVolume(const G4VPhysicalVolume* newWorld)",
198  "FastSim004",
200  tellWhatIsWrong);
201  }
202 }
203 
204 
205 // --------------------
206 // Start/End tracking:
207 // --------------------
208 void
211 {
212  fIsTrackingTime = true;
213  fIsFirstStep = true;
214 
215  // -- fetch the navigator (and its index) and activate it:
217  fGhostNavigator = transportationManager->GetNavigator(fWorldVolume);
218  fIsGhostGeometry = (fGhostNavigator != transportationManager->GetNavigatorForTracking());
220  else fGhostNavigatorIndex = -1;
221 
223 }
224 
225 
226 void
229 {
230  fIsTrackingTime = false;
232 }
233 
234 
235 // ------------------------------------------
236 // PostStepGetPhysicalInteractionLength():
237 // ------------------------------------------
238 G4double
241  G4double,
243 {
244  // -- Get current volume, and check for presence of fast simulation manager.
245  // -- For the case of the navigator for tracking (fGhostNavigatorIndex == 0)
246  // -- we use the track volume. This allows the code to be valid for both
247  // -- cases where the PathFinder is used (G4CoupledTranportation) or not
248  // -- (G4Transportation).
249  const G4VPhysicalVolume* currentVolume(0);
251  else currentVolume = track.GetVolume();
252 
253  if ( currentVolume )
254  {
257  {
258  // Ask for trigger:
261  {
262  // Take control over stepping:
263  *condition = ExclusivelyForced;
264  return 0.0;
265  }
266  }
267  }
268 
269  // -- no fast simulation occuring there:
270  *condition = NotForced;
271  return DBL_MAX;
272 }
273 
274 //------------------------------------
275 // PostStepDoIt()
276 //------------------------------------
280  const G4Step&)
281 {
283 
284  // If the particle is still alive, suspend it to force physics re-initialisation:
285  if (finalState->GetTrackStatus() != fStopAndKill) finalState->ProposeTrackStatus(fSuspend);
286 
287  return finalState;
288 }
289 
290 
291 G4double
294  G4double previousStepSize,
295  G4double currentMinimumStep,
296  G4double& proposedSafety,
297  G4GPILSelection* selection)
298 {
299 
300  *selection = NotCandidateForSelection;
301  G4double returnedStep = DBL_MAX;
302 
303  // ---------------------------------------------------
304  // -- Below code valid for ghost geometry, otherwise
305  // -- useless for fast simulation attached to mass
306  // -- geometry. Warn user in case along used for
307  // -- mass geometry ?
308  // --------------------------------------------------
309  if ( fIsGhostGeometry )
310  {
311  static G4ThreadLocal G4FieldTrack *endTrack_G4MT_TLS_ = 0 ; if (!endTrack_G4MT_TLS_) endTrack_G4MT_TLS_ = new G4FieldTrack ('0') ; G4FieldTrack &endTrack = *endTrack_G4MT_TLS_;
312  static G4ThreadLocal ELimited *eLimited_G4MT_TLS_ = 0 ; if (!eLimited_G4MT_TLS_) eLimited_G4MT_TLS_ = new ELimited ; ELimited &eLimited = *eLimited_G4MT_TLS_;
313 
314  if (previousStepSize > 0.) fGhostSafety -= previousStepSize;
315  if (fGhostSafety < 0.) fGhostSafety = 0.0;
316 
317  // ------------------------------------------
318  // Determination of the proposed step length:
319  // ------------------------------------------
320  if (currentMinimumStep <= fGhostSafety && currentMinimumStep > 0.)
321  {
322  // -- No chance to limit the step, as proposed move inside safety
323  returnedStep = currentMinimumStep;
324  proposedSafety = fGhostSafety - currentMinimumStep;
325  }
326  else
327  {
328  // -- Proposed move exceeds safety, need to state
330  returnedStep = fPathFinder->ComputeStep(fFieldTrack,
331  currentMinimumStep,
333  track.GetCurrentStepNumber(),
334  fGhostSafety,
335  eLimited,
336  endTrack,
337  track.GetVolume());
338 
339  if(eLimited == kDoNot) fGhostSafety = fGhostNavigator->ComputeSafety(endTrack.GetPosition()); // -- step no limited by ghost
340  proposedSafety = fGhostSafety;
341  if (eLimited == kUnique || eLimited == kSharedOther) *selection = CandidateForSelection;
342  else if (eLimited == kSharedTransport) returnedStep *= (1.0 + 1.0e-9); // -- Expand to disable its selection in Step Manager comparison
343  }
344  }
345 
346 
347  // ----------------------------------------------
348  // Returns the fGhostSafety as the proposedSafety
349  // The SteppingManager will take care of keeping
350  // the smallest one.
351  // ----------------------------------------------
352  return returnedStep;
353 }
354 
357 AlongStepDoIt(const G4Track& track,
358  const G4Step&)
359 {
361  return &fDummyParticleChange;
362 }
363 
364 
365 //--------------------------------------------
366 // At Rest parameterisation:
367 //--------------------------------------------
368 // AtRestGetPhysiscalInteractionLength:
369 //--------------------------------------------
370 G4double
374 {
375  const G4VPhysicalVolume* currentVolume(0);
377  else currentVolume = track.GetVolume();
380  {
381  // Ask for trigger:
384  {
385  // Dirty trick to take control over stepping. Does anyone will ever use that ?
386  *condition = NotForced;
387  return -1.0;
388  }
389  }
390 
391  // -- no fast simulation occuring there:
392  *condition = NotForced;
393  return DBL_MAX;
394 }
395 
396 //-----------------------------------------------
397 // AtRestDoIt:
398 //-----------------------------------------------
400 {
402 }
403 
404 
406 {
407  /* G4cout << " >>>>> Trigger Status : ";
408  switch(fFastSimulationManager->GetTriggerStatus())
409  {
410  case NoModel:
411  G4cout << "NoModel" << G4endl;
412  break;
413  case OnBoundaryButLeaving:
414  G4cout << "OnBoundaryButLeaving" << G4endl;
415  break;
416  case OneModelTrigger:
417  G4cout << "OneModelTrigger" << G4endl;
418  break;
419  case NoModelTrigger:
420  G4cout << "NoModelTrigger" << G4endl;
421  break;
422  case Undefined:
423  G4cout << "Undefined" << G4endl;
424  break;
425  default:
426  G4cout << " Bizarre..." << G4endl;
427  break;
428  }*/
429 }
430 
431 
void PrepareNewTrack(const G4ThreeVector &position, const G4ThreeVector &direction, G4VPhysicalVolume *massStartVol=0)
static G4PathFinder * GetInstance()
Definition: G4PathFinder.cc:57
G4double condition(const G4ErrorSymMatrix &m)
void RemoveFSMP(G4FastSimulationManagerProcess *)
G4VPhysicalVolume * IsWorldExisting(const G4String &worldName)
virtual void Initialize(const G4Track &)
G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
G4VParticleChange * InvokePostStepDoIt()
G4int verboseLevel
Definition: G4VProcess.hh:368
G4bool AtRestGetFastSimulationManagerTrigger(const G4Track &, const G4Navigator *a=0)
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
ELimited
const G4ThreeVector & GetPosition() const
G4Navigator * GetNavigatorForTracking() const
static G4GlobalFastSimulationManager * GetGlobalFastSimulationManager()
G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &step)
#define G4ThreadLocal
Definition: tls.hh:89
G4GLOB_DLL std::ostream G4cout
G4int GetCurrentStepNumber() const
const G4String & GetName() const
G4VPhysicalVolume * GetLocatedVolume(G4int navId) const
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
G4FastSimulationManager * GetFastSimulationManager() const
void DeActivateNavigator(G4Navigator *aNavigator)
G4double ComputeStep(const G4FieldTrack &pFieldTrack, G4double pCurrentProposedStepLength, G4int navigatorId, G4int stepNo, G4double &pNewSafety, ELimited &limitedStep, G4FieldTrack &EndState, G4VPhysicalVolume *currentVolume)
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
Definition: G4Step.hh:76
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4TransportationManager * GetTransportationManager()
G4int ActivateNavigator(G4Navigator *aNavigator)
const G4ThreeVector & GetMomentumDirection() const
G4LogicalVolume * GetLogicalVolume() const
G4VParticleChange * InvokeAtRestDoIt()
G4Navigator * GetNavigator(const G4String &worldName)
G4FastSimulationManagerProcess(const G4String &processName="G4FastSimulationManagerProcess", G4ProcessType theType=fParameterisation)
G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
G4VPhysicalVolume * GetVolume() const
#define G4endl
Definition: G4ios.hh:61
G4TrackStatus GetTrackStatus() const
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
virtual G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4double pProposedMaxLength=DBL_MAX, const G4bool keepState=true)
G4ForceCondition
G4bool PostStepGetFastSimulationManagerTrigger(const G4Track &, const G4Navigator *a=0)
G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
#define DBL_MAX
Definition: templates.hh:83
G4VPhysicalVolume * GetWorldVolume() const
static void Update(G4FieldTrack *, const G4Track *)
void AddFSMP(G4FastSimulationManagerProcess *)
G4GPILSelection
G4ProcessType