Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 101152 2016-11-08 08:07:39Z 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 ( nullptr ),
58  fIsTrackingTime ( false ),
59  fIsFirstStep ( false ),
60  fGhostNavigator ( nullptr ),
61  fGhostNavigatorIndex ( -1 ),
62  fIsGhostGeometry ( false ),
63  fGhostSafety ( -1.0 ),
64  fFieldTrack ( '0' ),
65  fFastSimulationManager( nullptr ),
66  fFastSimulationTrigger( false )
67 {
68  // -- set Process Sub Type
69  SetProcessSubType(static_cast<int>(FASTSIM_ManagerProcess));
70 
71 
72  fPathFinder = G4PathFinder::GetInstance();
73  fTransportationManager = G4TransportationManager::GetTransportationManager();
74 
75  SetWorldVolume(fTransportationManager->GetNavigatorForTracking()->GetWorldVolume()->GetName());
76  if (verboseLevel>0) G4cout << "G4FastSimulationManagerProcess `" << GetProcessName()
77  << "' is created, and will message geometry with world volume `"
78  << fWorldVolume->GetName() << "'." << G4endl;
80 }
81 
82 
85  const G4String& worldVolumeName,
86  G4ProcessType theType) :
87  G4VProcess(processName,theType),
88  fWorldVolume ( nullptr ),
89  fIsTrackingTime ( false ),
90  fIsFirstStep ( false ),
91  fGhostNavigator ( nullptr ),
92  fGhostNavigatorIndex ( -1 ),
93  fIsGhostGeometry ( false ),
94  fGhostSafety ( -1.0 ),
95  fFieldTrack ( '0' ),
96  fFastSimulationManager( nullptr ),
97  fFastSimulationTrigger( false )
98 {
99  // -- set Process Sub Type
100  SetProcessSubType(static_cast<int>(FASTSIM_ManagerProcess));
101 
102 
103  fPathFinder = G4PathFinder::GetInstance();
104  fTransportationManager = G4TransportationManager::GetTransportationManager();
105 
106  SetWorldVolume(worldVolumeName);
107  if (verboseLevel>0) G4cout << "G4FastSimulationManagerProcess `" << GetProcessName()
108  << "' is created, and will message geometry with world volume `"
109  << fWorldVolume->GetName() << "'." << G4endl;
111 }
112 
113 
116  G4VPhysicalVolume* worldVolume,
117  G4ProcessType theType) :
118  G4VProcess(processName,theType),
119  fWorldVolume ( nullptr ),
120  fIsTrackingTime ( false ),
121  fIsFirstStep ( false ),
122  fGhostNavigator ( nullptr ),
123  fGhostNavigatorIndex ( -1 ),
124  fIsGhostGeometry ( false ),
125  fGhostSafety ( -1.0 ),
126  fFieldTrack ( '0' ),
127  fFastSimulationManager( nullptr ),
128  fFastSimulationTrigger( false )
129 {
130  // -- set Process Sub Type
131  SetProcessSubType(static_cast<int>(FASTSIM_ManagerProcess));
132 
133 
134  fPathFinder = G4PathFinder::GetInstance();
135  fTransportationManager = G4TransportationManager::GetTransportationManager();
136 
137  SetWorldVolume(worldVolume);
138  if (verboseLevel>0) G4cout << "G4FastSimulationManagerProcess `" << GetProcessName()
139  << "' is created, and will message geometry with world volume `"
140  << fWorldVolume->GetName() << "'." << G4endl;
142 }
143 
144 
146 {
148 }
149 
150 
151 // -----------------------
152 // User access methods:
153 // -----------------------
155 {
156  if (fIsTrackingTime)
157  {
159  ed << "G4FastSimulationManagerProcess `" << GetProcessName()
160  << "': changing of world volume at tracking time is not allowed." << G4endl;
161  G4Exception("G4FastSimulationManagerProcess::SetWorldVolume(const G4String)",
162  "FastSim002",
163  JustWarning, ed,
164  "Call ignored.");
165  }
166  else
167  {
168  G4VPhysicalVolume* newWorld = fTransportationManager->IsWorldExisting(newWorldName);
169  if (newWorld == 0)
170  {
171  G4ExceptionDescription tellWhatIsWrong;
172  tellWhatIsWrong << "Volume newWorldName = `" << newWorldName
173  << "' is not a parallel world nor the mass world volume."
174  << G4endl;
175  G4Exception("G4FastSimulationManagerProcess::SetWorldVolume(const G4String)",
176  "FastSim003",
178  tellWhatIsWrong);
179  }
180  if (verboseLevel>0)
181  {
182  if (fWorldVolume) G4cout << "G4FastSimulationManagerProcess `" << GetProcessName()
183  << "': changing world volume from '" << fWorldVolume->GetName()
184  << "' to `" << newWorld << "'." << G4endl;
185  else G4cout << "G4FastSimulationManagerProcess `" << GetProcessName()
186  << "': setting world volume from to `"<< newWorld->GetName() << "'." << G4endl;
187  }
188  fWorldVolume = newWorld;
189  }
190 }
191 
192 
194 {
195  if (newWorld) SetWorldVolume(newWorld->GetName());
196  else
197  {
198  G4ExceptionDescription tellWhatIsWrong;
199  tellWhatIsWrong << "Null pointer passed for world volume." << G4endl;
200  G4Exception("G4FastSimulationManagerProcess::SetWorldVolume(const G4VPhysicalVolume* newWorld)",
201  "FastSim004",
203  tellWhatIsWrong);
204  }
205 }
206 
207 
208 // --------------------
209 // Start/End tracking:
210 // --------------------
212 {
213  fIsTrackingTime = true;
214  fIsFirstStep = true;
215 
216  // -- fetch the navigator (and its index) and activate it:
218  fGhostNavigator = transportationManager->GetNavigator(fWorldVolume);
219  fIsGhostGeometry = (fGhostNavigator != transportationManager->GetNavigatorForTracking());
220  if (fIsGhostGeometry) fGhostNavigatorIndex = transportationManager->ActivateNavigator(fGhostNavigator);
221  else fGhostNavigatorIndex = -1;
222 
223  fPathFinder->PrepareNewTrack(track->GetPosition(), track->GetMomentumDirection());
224 }
225 
226 
227 void
230 {
231  fIsTrackingTime = false;
232  if ( fIsGhostGeometry ) fTransportationManager->DeActivateNavigator(fGhostNavigator);
233 }
234 
235 
236 // ------------------------------------------
237 // PostStepGetPhysicalInteractionLength():
238 // ------------------------------------------
239 G4double
242  G4double,
244 {
245  // -- Get current volume, and check for presence of fast simulation manager.
246  // -- For the case of the navigator for tracking (fGhostNavigatorIndex == 0)
247  // -- we use the track volume. This allows the code to be valid for both
248  // -- cases where the PathFinder is used (G4CoupledTranportation) or not
249  // -- (G4Transportation).
250  const G4VPhysicalVolume* currentVolume(0);
251  if ( fIsGhostGeometry ) currentVolume = fPathFinder->GetLocatedVolume(fGhostNavigatorIndex);
252  else currentVolume = track.GetVolume();
253 
254  if ( currentVolume )
255  {
256  fFastSimulationManager = currentVolume->GetLogicalVolume()->GetFastSimulationManager();
257  if( fFastSimulationManager )
258  {
259  // Ask for trigger:
260  fFastSimulationTrigger = fFastSimulationManager->PostStepGetFastSimulationManagerTrigger(track, fGhostNavigator);
261  if( fFastSimulationTrigger )
262  {
263  // Take control over stepping:
264  *condition = ExclusivelyForced;
265  return 0.0;
266  }
267  }
268  }
269 
270  // -- no fast simulation occuring there:
271  *condition = NotForced;
272  return DBL_MAX;
273 }
274 
275 //------------------------------------
276 // PostStepDoIt()
277 //------------------------------------
281  const G4Step&)
282 {
283  G4VParticleChange* finalState = fFastSimulationManager->InvokePostStepDoIt();
284 
285  // If the particle is still alive, suspend it to force physics re-initialisation:
286  if (finalState->GetTrackStatus() != fStopAndKill) finalState->ProposeTrackStatus(fSuspend);
287 
288  return finalState;
289 }
290 
291 
292 G4double
295  G4double previousStepSize,
296  G4double currentMinimumStep,
297  G4double& proposedSafety,
298  G4GPILSelection* selection)
299 {
300 
301  *selection = NotCandidateForSelection;
302  G4double returnedStep = DBL_MAX;
303 
304  // ---------------------------------------------------
305  // -- Below code valid for ghost geometry, otherwise
306  // -- useless for fast simulation attached to mass
307  // -- geometry. Warn user in case along used for
308  // -- mass geometry ?
309  // --------------------------------------------------
310  if ( fIsGhostGeometry )
311  {
312  static G4ThreadLocal G4FieldTrack *endTrack_G4MT_TLS_ = 0 ;
313  if (!endTrack_G4MT_TLS_) endTrack_G4MT_TLS_ = new G4FieldTrack ('0') ;
314  G4FieldTrack &endTrack = *endTrack_G4MT_TLS_;
315 
316  static G4ThreadLocal ELimited *eLimited_G4MT_TLS_ = 0 ;
317  if (!eLimited_G4MT_TLS_) eLimited_G4MT_TLS_ = new ELimited ;
318  ELimited &eLimited = *eLimited_G4MT_TLS_;
319 
320  if (previousStepSize > 0.) fGhostSafety -= previousStepSize;
321  if (fGhostSafety < 0.) fGhostSafety = 0.0;
322 
323  // ------------------------------------------
324  // Determination of the proposed step length:
325  // ------------------------------------------
326  if (currentMinimumStep <= fGhostSafety && currentMinimumStep > 0.)
327  {
328  // -- No chance to limit the step, as proposed move inside safety
329  returnedStep = currentMinimumStep;
330  proposedSafety = fGhostSafety - currentMinimumStep;
331  }
332  else
333  {
334  // -- Proposed move exceeds safety, need to state
335  G4FieldTrackUpdator::Update(&fFieldTrack, &track);
336  returnedStep = fPathFinder->ComputeStep(fFieldTrack,
337  currentMinimumStep,
338  fGhostNavigatorIndex,
339  track.GetCurrentStepNumber(),
340  fGhostSafety,
341  eLimited,
342  endTrack,
343  track.GetVolume());
344 
345  if(eLimited == kDoNot) fGhostSafety = fGhostNavigator->ComputeSafety(endTrack.GetPosition()); // -- step no limited by ghost
346  proposedSafety = fGhostSafety;
347  if (eLimited == kUnique || eLimited == kSharedOther) *selection = CandidateForSelection;
348  else if (eLimited == kSharedTransport) returnedStep *= (1.0 + 1.0e-9); // -- Expand to disable its selection in Step Manager comparison
349  }
350  }
351 
352 
353  // ----------------------------------------------
354  // Returns the fGhostSafety as the proposedSafety
355  // The SteppingManager will take care of keeping
356  // the smallest one.
357  // ----------------------------------------------
358  return returnedStep;
359 }
360 
363 AlongStepDoIt(const G4Track& track,
364  const G4Step&)
365 {
366  fDummyParticleChange.Initialize(track);
367  return &fDummyParticleChange;
368 }
369 
370 
371 //--------------------------------------------
372 // At Rest parameterisation:
373 //--------------------------------------------
374 // AtRestGetPhysiscalInteractionLength:
375 //--------------------------------------------
376 G4double
380 {
381  const G4VPhysicalVolume* currentVolume(0);
382  if ( fIsGhostGeometry ) currentVolume = fPathFinder->GetLocatedVolume(fGhostNavigatorIndex);
383  else currentVolume = track.GetVolume();
384  fFastSimulationManager = currentVolume->GetLogicalVolume()->GetFastSimulationManager();
385  if( fFastSimulationManager )
386  {
387  // Ask for trigger:
388  fFastSimulationTrigger = fFastSimulationManager->AtRestGetFastSimulationManagerTrigger(track, fGhostNavigator);
389  if( fFastSimulationTrigger )
390  {
391  // Dirty trick to take control over stepping. Does anyone will ever use that ?
392  *condition = NotForced;
393  return -1.0;
394  }
395  }
396 
397  // -- no fast simulation occuring there:
398  *condition = NotForced;
399  return DBL_MAX;
400 }
401 
402 //-----------------------------------------------
403 // AtRestDoIt:
404 //-----------------------------------------------
406 {
407  return fFastSimulationManager->InvokeAtRestDoIt();
408 }
409 
410 
412 {
413  /* G4cout << " >>>>> Trigger Status : ";
414  switch(fFastSimulationManager->GetTriggerStatus())
415  {
416  case NoModel:
417  G4cout << "NoModel" << G4endl;
418  break;
419  case OnBoundaryButLeaving:
420  G4cout << "OnBoundaryButLeaving" << G4endl;
421  break;
422  case OneModelTrigger:
423  G4cout << "OneModelTrigger" << G4endl;
424  break;
425  case NoModelTrigger:
426  G4cout << "NoModelTrigger" << G4endl;
427  break;
428  case Undefined:
429  G4cout << "Undefined" << G4endl;
430  break;
431  default:
432  G4cout << " Bizarre..." << G4endl;
433  break;
434  }*/
435 }
436 
437 
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