Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4ITPathFinder.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: G4ITPathFinder.cc 73734 2013-09-09 15:01:46Z gcosmo $
28 // GEANT4 tag $ Name: $
29 //
30 // class G4ITPathFinder Implementation
31 //
32 // Original author: John Apostolakis, April 2006
33 //
34 // --------------------------------------------------------------------
35 
36 #include <iomanip>
37 
38 #include "G4ITPathFinder.hh"
39 
40 #include "G4SystemOfUnits.hh"
41 #include "G4GeometryTolerance.hh"
42 #include "G4ITNavigator.hh"
43 #include "G4PropagatorInField.hh"
45 #include "G4ITMultiNavigator.hh"
46 #include "G4ITSafetyHelper.hh"
47 
48 // to ease comparaison with G4PathFinder
49 #define State(X) fpTrackState->X
50 #define fNewTrack State(fNewTrack)
51 #define fLimitedStep State(fLimitedStep)
52 #define fLimitTruth State(fLimitTruth)
53 #define fCurrentStepSize State(fCurrentStepSize)
54 #define fNoGeometriesLimiting State(fNoGeometriesLimiting)
55 #define fPreSafetyLocation State(fPreSafetyLocation)
56 #define fPreSafetyMinValue State(fPreSafetyMinValue)
57 #define fPreSafetyValues State(fPreSafetyValues)
58 #define fPreStepLocation State(fPreStepLocation)
59 #define fMinSafety_PreStepPt State(fMinSafety_PreStepPt)
60 #define fCurrentPreStepSafety State(fCurrentPreStepSafety)
61 #define fPreStepCenterRenewed State(fPreStepCenterRenewed)
62 #define fMinStep State(fMinStep)
63 #define fTrueMinStep State(fTrueMinStep)
64 #define fLocatedVolume State(fLocatedVolume)
65 #define fLastLocatedPosition State(fLastLocatedPosition)
66 #define fEndState State(fEndState)
67 #define fFieldExertedForce State(fFieldExertedForce)
68 #define fRelocatedPoint State(fRelocatedPoint)
69 #define fSafetyLocation State(fSafetyLocation)
70 #define fMinSafety_atSafLocation State(fMinSafety_atSafLocation)
71 #define fNewSafetyComputed State(fNewSafetyComputed)
72 #define fLastStepNo State(fLastStepNo)
73 #define fCurrentStepNo State(fCurrentStepNo)
74 
75 // Initialise the static instance of the singleton
76 //
77 G4ThreadLocal G4ITPathFinder* G4ITPathFinder::fpPathFinder=0;
78 
79 // ----------------------------------------------------------------------------
80 // GetInstance()
81 //
82 // Retrieve the static instance of the singleton
83 //
85 {
86  if( ! fpPathFinder )
87  {
88  fpPathFinder = new G4ITPathFinder;
89  }
90  return fpPathFinder;
91 }
92 
93 // ----------------------------------------------------------------------------
94 // Constructor
95 //
97  fVerboseLevel(0)
98 {
99  fpMultiNavigator= new G4ITMultiNavigator();
100 
102  // fpFieldPropagator = fpTransportManager->GetPropagatorInField();
103 
105 
106  fNoActiveNavigators= 0;
107  for( G4int num=0; num< G4ITNavigator::fMaxNav; ++num )
108  {
109  fpNavigator[num] = 0;
110  }
111 }
112 
113 // ----------------------------------------------------------------------------
114 // Destructor
115 //
117 {
118  delete fpMultiNavigator;
119  if (fpPathFinder) { delete fpPathFinder; fpPathFinder=0; }
120 }
121 
122 // ----------------------------------------------------------------------------
123 //
124 void
126 {
127 /*
128  G4ITNavigator *navigatorForPropagation=0, *massNavigator=0;
129 
130  massNavigator= fpTransportManager->GetNavigatorForTracking();
131 */
132  if( enableChoice )
133  {
134  // navigatorForPropagation= fpMultiNavigator;
135 
136  // Enable SafetyHelper to use PF
137  //
138  fpTransportManager->GetSafetyHelper()->EnableParallelNavigation(true);
139  }
140  else
141  {
142  // navigatorForPropagation= massNavigator;
143 
144  // Disable SafetyHelper to use PF
145  //
146  fpTransportManager->GetSafetyHelper()->EnableParallelNavigation(false);
147  }
148  // fpFieldPropagator->SetNavigatorForPropagating(navigatorForPropagation);
149 }
150 
151 // ----------------------------------------------------------------------------
152 //
153 G4double
154 G4ITPathFinder::ComputeStep( const G4FieldTrack &InitialFieldTrack,
155  G4double proposedStepLength,
156  G4int navigatorNo,
157  G4int stepNo, // find next step
158  G4double &pNewSafety, // for this geom
159  ELimited &limitedStep,
160  G4FieldTrack &EndState,
161  G4VPhysicalVolume* /*currentVolume*/)
162 {
163  G4double possibleStep= -1.0;
164 
165 #ifdef G4DEBUG_PATHFINDER
166  if( fVerboseLevel > 2 )
167  {
168  G4cout << " -------------------------" << G4endl;
169  G4cout << " G4ITPathFinder::ComputeStep - entered " << G4endl;
170  G4cout << " - stepNo = " << std::setw(4) << stepNo << " "
171  << " navigatorId = " << std::setw(2) << navigatorNo << " "
172  << " proposed step len = " << proposedStepLength << " " << G4endl;
173  G4cout << " PF::ComputeStep requested step "
174  << " from " << InitialFieldTrack.GetPosition()
175  << " dir " << InitialFieldTrack.GetMomentumDirection() << G4endl;
176  }
177 #endif
178 #ifdef G4VERBOSE
179  if( navigatorNo >= fNoActiveNavigators )
180  {
181  std::ostringstream message;
182  message << "Bad Navigator ID !" << G4endl
183  << " Requested Navigator ID = " << navigatorNo << G4endl
184  << " Number of active navigators = " << fNoActiveNavigators;
185  G4Exception("G4ITPathFinder::ComputeStep()", "GeomNav0002",
186  FatalException, message);
187  }
188 #endif
189 
190  if( fNewTrack || (stepNo != fLastStepNo) )
191  {
192  // This is a new track or a new step, so we must make the step
193  // ( else we can simply retrieve its results for this Navigator Id )
194 
195  G4FieldTrack currentState= InitialFieldTrack;
196 
197  fCurrentStepNo = stepNo;
198 
199  // Check whether a process shifted the position
200  // since the last step -- by physics processes
201  //
202  G4ThreeVector newPosition = InitialFieldTrack.GetPosition();
203  G4ThreeVector moveVector= newPosition - fLastLocatedPosition;
204  G4double moveLenSq= moveVector.mag2();
205  if( moveLenSq > kCarTolerance * kCarTolerance )
206  {
207  G4ThreeVector newDirection = InitialFieldTrack.GetMomentumDirection();
208 #ifdef G4DEBUG_PATHFINDER
209  if( fVerboseLevel > 2 )
210  {
211  G4double moveLen= std::sqrt( moveLenSq );
212  G4cout << " G4ITPathFinder::ComputeStep : Point moved since last step "
213  << " -- at step # = " << stepNo << G4endl
214  << " by " << moveLen << " to " << newPosition << G4endl;
215  }
216 #endif
217  MovePoint(); // Unintentional changed -- ????
218 
219  // Relocate to cope with this move -- else could abort !?
220  //
221  Locate( newPosition, newDirection );
222  }
223 
224  // Check whether the particle have an (EM) field force exerting upon it
225  //
226  /*
227  G4double particleCharge= currentState.GetCharge();
228 
229  G4FieldManager* fieldMgr=0;
230  G4bool fieldExertsForce = false ;
231  if( (particleCharge != 0.0) )
232  {
233  fieldMgr= fpFieldPropagator->FindAndSetFieldManager( currentVolume );
234 
235  // Protect for case where field manager has no field (= field is zero)
236  //
237  fieldExertsForce = (fieldMgr != 0)
238  && (fieldMgr->GetDetectorField() != 0);
239  }
240  fFieldExertedForce = fieldExertsForce; // Store for use in later calls
241  // referring to this 'step'.
242 
243  fNoGeometriesLimiting= -1; // At start of track, no process limited step
244  if( fieldExertsForce )
245  {
246  DoNextCurvedStep( currentState, proposedStepLength, currentVolume );
247  //--------------
248  }else{
249  */
250  DoNextLinearStep( currentState, proposedStepLength );
251  //--------------
252 // }
253  fLastStepNo= stepNo;
254 
255 #ifdef G4DEBUG_PATHFINDER
256  if ( (fNoGeometriesLimiting < 0)
257  || (fNoGeometriesLimiting > fNoActiveNavigators) )
258  {
259  std::ostringstream message;
260  message << "Number of geometries limiting the step not set." << G4endl
261  << " Number of geometries limiting step = "
263  G4Exception("G4ITPathFinder::ComputeStep()",
264  "GeomNav0002", FatalException, message);
265  }
266 #endif
267  }
268 #ifdef G4DEBUG_PATHFINDER
269  else
270  {
271  if( proposedStepLength < fTrueMinStep ) // For 2nd+ geometry
272  {
273  std::ostringstream message;
274  message << "Problem in step size request." << G4endl
275  << " Error can be caused by incorrect process ordering."
276  << " Being requested to make a step which is shorter"
277  << " than the minimum Step " << G4endl
278  << " already computed for any Navigator/geometry during"
279  << " this tracking-step: " << G4endl
280  << " This can happen due to an error in process ordering."
281  << G4endl
282  << " Check that all physics processes are registered"
283  << G4endl
284  << " before all processes with a navigator/geometry."
285  << G4endl
286  << " If using pre-packaged physics list and/or"
287  << G4endl
288  << " functionality, please report this error."
289  << G4endl << G4endl
290  << " Additional information for problem: " << G4endl
291  << " Steps request/proposed = " << proposedStepLength
292  << G4endl
293  << " MinimumStep (true) = " << fTrueMinStep
294  << G4endl
295  << " MinimumStep (navraw) = " << fMinStep
296  << G4endl
297  << " Navigator raw return value" << G4endl
298  << " Requested step now = " << proposedStepLength
299  << G4endl
300  << " Difference min-req = "
301  << fTrueMinStep-proposedStepLength << G4endl
302  << " -- Step info> stepNo= " << stepNo
303  << " last= " << fLastStepNo
304  << " newTr= " << fNewTrack << G4endl;
305  G4Exception("G4ITPathFinder::ComputeStep()",
306  "GeomNav0003", FatalException, message);
307  }
308  else
309  {
310  // This is neither a new track nor a new step -- just another
311  // client accessing information for the current track, step
312  // We will simply retrieve the results of the synchronous
313  // stepping for this Navigator Id below.
314  //
315  if( fVerboseLevel > 1 )
316  {
317  G4cout << " G4P::CS -> Not calling DoNextLinearStep: "
318  << " stepNo= " << stepNo << " last= " << fLastStepNo
319  << " new= " << fNewTrack << " Step already done" << G4endl;
320  }
321  }
322  }
323 #endif
324 
325  fNewTrack= false;
326 
327  // Prepare the information to return
328 
329  pNewSafety = fCurrentPreStepSafety[ navigatorNo ];
330  limitedStep = fLimitedStep[ navigatorNo ];
331  fRelocatedPoint= false;
332 
333  possibleStep= std::min(proposedStepLength, fCurrentStepSize[ navigatorNo ]);
334  EndState = fEndState; // now corrected for smaller step, if needed
335 
336 #ifdef G4DEBUG_PATHFINDER
337  if( fVerboseLevel > 0 )
338  {
339  G4cout << " G4ITPathFinder::ComputeStep returns "
340  << fCurrentStepSize[ navigatorNo ]
341  << " for Navigator " << navigatorNo
342  << " Limited step = " << limitedStep
343  << " Safety(mm) = " << pNewSafety / mm
344  << G4endl;
345  }
346 #endif
347 
348  return possibleStep;
349 }
350 
351 // ----------------------------------------------------------------------
352 
353 void
355  const G4ThreeVector& direction,
356  G4VPhysicalVolume* massStartVol)
357 {
358  // Key purposes:
359  // - Check and cache set of active navigators
360  // - Reset state for new track
361 
362  G4int num=0;
363 
365  // Switch PropagatorInField to use MultiNavigator
366 
367  fpTransportManager->GetSafetyHelper()->InitialiseHelper();
368  // Reinitialise state of safety helper -- avoid problems with overlaps
369 
370  fNewTrack= true;
371  this->MovePoint(); // Signal further that the last status is wiped
372 
373  // Message the G4NavigatorPanel / Dispatcher to find active navigators
374  //
375  std::vector<G4ITNavigator*>::iterator pNavigatorIter;
376 
377  fNoActiveNavigators= fpTransportManager-> GetNoActiveNavigators();
378  if( fNoActiveNavigators > G4ITNavigator::fMaxNav )
379  {
380  std::ostringstream message;
381  message << "Too many active Navigators / worlds." << G4endl
382  << " Transportation Manager has "
383  << fNoActiveNavigators << " active navigators." << G4endl
384  << " This is more than the number allowed = "
385  << G4ITNavigator::fMaxNav << " !";
386  G4Exception("G4ITPathFinder::PrepareNewTrack()", "GeomNav0002",
387  FatalException, message);
388  }
389 
390  fpMultiNavigator->PrepareNavigators();
391  //------------------------------------
392 
393  pNavigatorIter= fpTransportManager->GetActiveNavigatorsIterator();
394  for( num=0; num< fNoActiveNavigators; ++pNavigatorIter,++num )
395  {
396  // Keep information in C-array ... for creating touchables - at least
397 
398  fpNavigator[num] = *pNavigatorIter;
399  fLimitTruth[num] = false;
400  fLimitedStep[num] = kDoNot;
401  fCurrentStepSize[num] = 0.0;
402  fLocatedVolume[num] = 0;
403  }
404  fNoGeometriesLimiting= 0; // At start of track, no process limited step
405 
406  // In case of one geometry, the tracking will have done the locating!!
407 
408  if( fNoActiveNavigators > 1 )
409  {
410  Locate( position, direction, false );
411  }
412  else
413  {
414  // Update state -- depending on the tracking's call to Mass Navigator
415 
417  fLocatedVolume[0]= massStartVol; // This information must be given
418  // by transportation
419  fLimitedStep[0] = kDoNot;
420  fCurrentStepSize[0] = 0.0;
421  }
422 
423  // Reset Safety Information -- as in case of overlaps this can cause
424  // inconsistencies ...
425  //
427 
428  for( num=0; num< fNoActiveNavigators; ++num )
429  {
430  fPreSafetyValues[num]= 0.0;
431  fNewSafetyComputed[num]= 0.0;
432  fCurrentPreStepSafety[num] = 0.0;
433  }
434 
435  // The first location for each Navigator must be non-relative
436  // or else call ResetStackAndState() for each Navigator
437 
438  fRelocatedPoint= false;
439 }
440 
442  const G4ThreeVector& NewVector,
443  const G4String& Quantity ) const
444 {
445  G4ThreeVector moveVec = ( NewVector - OldVector );
446 
447  G4int prc= G4cerr.precision(12);
448  std::ostringstream message;
449  message << "Endpoint moved between value returned by ComputeStep()"
450  << " and call to Locate(). " << G4endl
451  << " Change of " << Quantity << " is "
452  << moveVec.mag() / mm << " mm long" << G4endl
453  << " and its vector is "
454  << (1.0/mm) * moveVec << " mm " << G4endl
455  << " Endpoint of ComputeStep() was " << OldVector << G4endl
456  << " and current position to locate is " << NewVector;
457  G4Exception("G4ITPathFinder::ReportMove()", "GeomNav1002",
458  JustWarning, message);
459  G4cerr.precision(prc);
460 }
461 
462 void
464  const G4ThreeVector& direction,
465  G4bool relative)
466 {
467  // Locate the point in each geometry
468 
469  std::vector<G4ITNavigator*>::iterator pNavIter=
470  fpTransportManager->GetActiveNavigatorsIterator();
471 
472  G4ThreeVector lastEndPosition= fEndState.GetPosition();
473  G4ThreeVector moveVec = (position - lastEndPosition );
474  G4double moveLenSq= moveVec.mag2();
475  if( (!fNewTrack) && (!fRelocatedPoint)
476  && ( moveLenSq> 10*kCarTolerance*kCarTolerance ) )
477  {
478  ReportMove( position, lastEndPosition, "Position" );
479  }
481 
482 #ifdef G4DEBUG_PATHFINDER
483  if( fVerboseLevel > 2 )
484  {
485  G4cout << G4endl;
486  G4cout << " G4ITPathFinder::Locate : entered " << G4endl;
487  G4cout << " -------------------- -------" << G4endl;
488  G4cout << " Locating at position " << position
489  << " with direction " << direction
490  << " relative= " << relative << G4endl;
491  if ( (fVerboseLevel > 1) || ( moveLenSq > 0.0) )
492  {
493  G4cout << " lastEndPosition = " << lastEndPosition
494  << " moveVec = " << moveVec
495  << " newTr = " << fNewTrack
496  << " relocated = " << fRelocatedPoint << G4endl;
497  }
498 
499  G4cout << " Located at " << position ;
500  if( fNoActiveNavigators > 1 ) { G4cout << G4endl; }
501  }
502 #endif
503 
504  for ( G4int num=0; num< fNoActiveNavigators ; ++pNavIter,++num )
505  {
506  // ... who limited the step ....
507 
508  if( fLimitTruth[num] ) { (*pNavIter)->SetGeometricallyLimitedStep(); }
509 
510  G4VPhysicalVolume *pLocated=
511  (*pNavIter)->LocateGlobalPointAndSetup( position, &direction,
512  relative,
513  false);
514  // Set the state related to the location
515  //
516  fLocatedVolume[num] = pLocated;
517 
518  // Clear state related to the step
519  //
520  fLimitedStep[num] = kDoNot;
521  fCurrentStepSize[num] = 0.0;
522 
523 #ifdef G4DEBUG_PATHFINDER
524  if( fVerboseLevel > 2 )
525  {
526  G4cout << " - In world " << num << " geomLimStep= " << fLimitTruth[num]
527  << " gives volume= " << pLocated ;
528  if( pLocated )
529  {
530  G4cout << " name = '" << pLocated->GetName() << "'";
531  G4cout << " - CopyNo= " << pLocated->GetCopyNo();
532  }
533  G4cout << G4endl;
534  }
535 #endif
536  }
537 
538  fRelocatedPoint= false;
539 }
540 
542 {
543  // Locate the point in each geometry
544 
545  std::vector<G4ITNavigator*>::iterator pNavIter=
546  fpTransportManager->GetActiveNavigatorsIterator();
547 
548 #ifdef G4DEBUG_PATHFINDER
549 
550  // Check that this relocation does not violate safety
551  // - at endpoint (computed from start point) AND
552  // - at last safety location (likely just called)
553 
554  G4ThreeVector lastEndPosition= fEndState.GetPosition();
555 
556  // Calculate end-point safety ...
557  //
558  G4double DistanceStartEnd= (lastEndPosition - fPreStepLocation).mag();
559  G4double endPointSafety_raw = fMinSafety_PreStepPt - DistanceStartEnd;
560  G4double endPointSafety_Est1 = std::max( 0.0, endPointSafety_raw );
561 
562  // ... and check move from endpoint against this endpoint safety
563  //
564  G4ThreeVector moveVecEndPos = position - lastEndPosition;
565  G4double moveLenEndPosSq = moveVecEndPos.mag2();
566 
567  // Check that move from endpoint of last step is within safety
568  // -- or check against last location or relocation ??
569  //
570  G4ThreeVector moveVecSafety= position - fSafetyLocation;
571  G4double moveLenSafSq= moveVecSafety.mag2();
572 
573  G4double distCheckEnd_sq= ( moveLenEndPosSq - endPointSafety_Est1
574  *endPointSafety_Est1 );
575  G4double distCheckSaf_sq= ( moveLenSafSq - fMinSafety_atSafLocation
577 
578  G4bool longMoveEnd = distCheckEnd_sq > 0.0;
579  G4bool longMoveSaf = distCheckSaf_sq > 0.0;
580 
581  G4double revisedSafety= 0.0;
582 
583  if( (!fNewTrack) && ( longMoveEnd && longMoveSaf ) )
584  {
585  // Recompute ComputeSafety for end position
586  //
587  revisedSafety= ComputeSafety(lastEndPosition);
588 
589  const G4double kRadTolerance =
591  const G4double cErrorTolerance=1e-12;
592  // Maximum relative error from roundoff of arithmetic
593 
594  G4double distCheckRevisedEnd= moveLenEndPosSq-revisedSafety*revisedSafety;
595 
596  G4bool longMoveRevisedEnd= ( distCheckRevisedEnd > 0. ) ;
597 
598  G4double moveMinusSafety= 0.0;
599  G4double moveLenEndPosition= std::sqrt( moveLenEndPosSq );
600  moveMinusSafety = moveLenEndPosition - revisedSafety;
601 
602  if ( longMoveRevisedEnd && (moveMinusSafety > 0.0 )
603  && ( revisedSafety > 0.0 ) )
604  {
605  // Take into account possibility of roundoff error causing
606  // this apparent move further than safety
607 
608  if( fVerboseLevel > 0 )
609  {
610  G4cout << " G4PF:Relocate> Ratio to revised safety is "
611  << std::fabs(moveMinusSafety)/revisedSafety << G4endl;
612  }
613 
614  G4double absMoveMinusSafety= std::fabs(moveMinusSafety);
615  G4bool smallRatio= absMoveMinusSafety < kRadTolerance * revisedSafety ;
616  G4double maxCoordPos = std::max(
617  std::max( std::fabs(position.x()),
618  std::fabs(position.y())),
619  std::fabs(position.z()) );
620  G4bool smallValue= absMoveMinusSafety < cErrorTolerance * maxCoordPos;
621  if( ! (smallRatio || smallValue) )
622  {
623  G4cout << " G4PF:Relocate> Ratio to revised safety is "
624  << std::fabs(moveMinusSafety)/revisedSafety << G4endl;
625  G4cout << " Difference of move and safety is not very small."
626  << G4endl;
627  }
628  else
629  {
630  moveMinusSafety = 0.0;
631  longMoveRevisedEnd = false; // Numerical issue -- not too long!
632 
633  G4cout << " Difference of move & safety is very small in magnitude, "
634  << absMoveMinusSafety << G4endl;
635  if( smallRatio )
636  {
637  G4cout << " ratio to safety " << revisedSafety
638  << " is " << absMoveMinusSafety / revisedSafety
639  << "smaller than " << kRadTolerance << " of safety ";
640  }
641  else
642  {
643  G4cout << " as fraction " << absMoveMinusSafety / maxCoordPos
644  << " of position vector max-coord " << maxCoordPos
645  << " smaller than " << cErrorTolerance ;
646  }
647  G4cout << " -- reset moveMinusSafety to "
648  << moveMinusSafety << G4endl;
649  }
650  }
651 
652  if ( longMoveEnd && longMoveSaf
653  && longMoveRevisedEnd && (moveMinusSafety>0.0) )
654  {
655  G4int oldPrec= G4cout.precision(9);
656  std::ostringstream message;
657  message << "ReLocation is further than end-safety value." << G4endl
658  << " Moved from last endpoint by " << moveLenEndPosition
659  << " compared to end safety (from preStep point) = "
660  << endPointSafety_Est1 << G4endl
661  << " --> last PreSafety Location was " << fPreSafetyLocation
662  << G4endl
663  << " safety value = " << fPreSafetyMinValue << G4endl
664  << " --> last PreStep Location was " << fPreStepLocation
665  << G4endl
666  << " safety value = " << fMinSafety_PreStepPt << G4endl
667  << " --> last EndStep Location was " << lastEndPosition
668  << G4endl
669  << " safety value = " << endPointSafety_Est1
670  << " raw-value = " << endPointSafety_raw << G4endl
671  << " --> Calling again at this endpoint, we get "
672  << revisedSafety << " as safety value." << G4endl
673  << " --> last position for safety " << fSafetyLocation
674  << G4endl
675  << " its safety value = " << fMinSafety_atSafLocation
676  << G4endl
677  << " move from safety location = "
678  << std::sqrt(moveLenSafSq) << G4endl
679  << " again= " << moveVecSafety.mag() << G4endl
680  << " safety - Move-from-end= "
681  << revisedSafety - moveLenEndPosition
682  << " (negative is Bad.)" << G4endl
683  << " Debug: distCheckRevisedEnd = "
684  << distCheckRevisedEnd;
685  ReportMove( lastEndPosition, position, "Position" );
686  G4Exception("G4ITPathFinder::ReLocate", "GeomNav0003",
687  FatalException, message);
688  G4cout.precision(oldPrec);
689  }
690  }
691 
692  if( fVerboseLevel > 2 )
693  {
694  G4cout << G4endl;
695  G4cout << " G4ITPathFinder::ReLocate : entered " << G4endl;
696  G4cout << " ---------------------- -------" << G4endl;
697  G4cout << " *Re*Locating at position " << position << G4endl;
698  // << " with direction " << direction
699  // << " relative= " << relative << G4endl;
700  if ( (fVerboseLevel > -1) || ( moveLenEndPosSq > 0.0) )
701  {
702  G4cout << " lastEndPosition = " << lastEndPosition
703  << " moveVec from step-end = " << moveVecEndPos
704  << " is new Track = " << fNewTrack
705  << " relocated = " << fRelocatedPoint << G4endl;
706  }
707  }
708 #endif // G4DEBUG_PATHFINDER
709 
710  for ( G4int num=0; num< fNoActiveNavigators ; ++pNavIter,++num )
711  {
712  // ... none limited the step
713 
714  (*pNavIter)->LocateGlobalPointWithinVolume( position );
715 
716  // Clear state related to the step
717  //
718  fLimitedStep[num] = kDoNot;
719  fCurrentStepSize[num] = 0.0;
720  fLimitTruth[num] = false;
721  }
722 
724  fRelocatedPoint= false;
725 
726 #ifdef G4DEBUG_PATHFINDER
727  if( fVerboseLevel > 2 )
728  {
729  G4cout << " G4ITPathFinder::ReLocate : exiting "
730  << " at position " << fLastLocatedPosition << G4endl << G4endl;
731  }
732 #endif
733 }
734 
735 // -----------------------------------------------------------------------------
736 
738 {
739  // Recompute safety for the relevant point
740 
741  G4double minSafety= kInfinity;
742 
743  std::vector<G4ITNavigator*>::iterator pNavigatorIter;
744  pNavigatorIter= fpTransportManager->GetActiveNavigatorsIterator();
745 
746  for( G4int num=0; num<fNoActiveNavigators; ++pNavigatorIter,++num )
747  {
748  G4double safety = (*pNavigatorIter)->ComputeSafety( position,true );
749  if( safety < minSafety ) { minSafety = safety; }
750  fNewSafetyComputed[num]= safety;
751  }
752 
754  fMinSafety_atSafLocation = minSafety;
755 
756 #ifdef G4DEBUG_PATHFINDER
757  if( fVerboseLevel > 1 )
758  {
759  G4cout << " G4ITPathFinder::ComputeSafety - returns "
760  << minSafety << " at location " << position << G4endl;
761  }
762 #endif
763  return minSafety;
764 }
765 
766 
767 // -----------------------------------------------------------------------------
768 
771 {
772 #ifdef G4DEBUG_PATHFINDER
773  if( fVerboseLevel > 2 )
774  {
775  G4cout << "G4ITPathFinder::CreateTouchableHandle : navId = "
776  << navId << " -- " << GetNavigator(navId) << G4endl;
777  }
778 #endif
779 
780  G4TouchableHistory* touchHist;
781  touchHist= GetNavigator(navId) -> CreateTouchableHistory();
782 
783  G4VPhysicalVolume* locatedVolume= fLocatedVolume[navId];
784  if( locatedVolume == 0 )
785  {
786  // Workaround to ensure that the touchable is fixed !! // TODO: fix
787 
788  touchHist->UpdateYourself( locatedVolume, touchHist->GetHistory() );
789  }
790 
791 #ifdef G4DEBUG_PATHFINDER
792  if( fVerboseLevel > 2 )
793  {
794  G4String VolumeName("None");
795  if( locatedVolume ) { VolumeName= locatedVolume->GetName(); }
796  G4cout << " Touchable History created at address " << touchHist
797  << " volume = " << locatedVolume << " name= " << VolumeName
798  << G4endl;
799  }
800 #endif
801 
802  return G4TouchableHandle(touchHist);
803 }
804 
805 G4double
807  G4double proposedStepLength )
808 {
809  std::vector<G4ITNavigator*>::iterator pNavigatorIter;
810  G4double safety= 0.0, step=0.0;
811  G4double minSafety= kInfinity, minStep;
812 
813  const G4int IdTransport= 0; // Id of Mass Navigator !!
814  G4int num=0;
815 
816 #ifdef G4DEBUG_PATHFINDER
817  if( fVerboseLevel > 2 )
818  {
819  G4cout << " G4ITPathFinder::DoNextLinearStep : entered " << G4endl;
820  G4cout << " Input field track= " << initialState << G4endl;
821  G4cout << " Requested step= " << proposedStepLength << G4endl;
822  }
823 #endif
824 
825  G4ThreeVector initialPosition= initialState.GetPosition();
826  G4ThreeVector initialDirection= initialState.GetMomentumDirection();
827 
828  G4ThreeVector OriginShift = initialPosition - fPreSafetyLocation;
829  G4double MagSqShift = OriginShift.mag2() ;
830  G4double MagShift; // Only given value if it larger than minimum safety
831 
832  // Potential optimisation using Maximum Value of safety!
833  // if( MagSqShift >= sqr(fPreSafetyMaxValue ) ){
834  // MagShift= kInfinity; // Not a useful value -- all will not use/ignore
835  // else
836  // MagShift= std::sqrt(MagSqShift) ;
837 
838  MagShift= std::sqrt(MagSqShift) ;
839 
840 #ifdef G4PATHFINDER_OPTIMISATION
841 
842  G4double fullSafety; // For all geometries, for prestep point
843 
844  if( MagSqShift >= sqr(fPreSafetyMinValue ) )
845  {
846  fullSafety = 0.0 ;
847  }
848  else
849  {
850  fullSafety = fPreSafetyMinValue - MagShift;
851  }
852  if( proposedStepLength < fullSafety )
853  {
854  // Move is smaller than all safeties
855  // -> so we do not have to move the safety center
856 
857  fPreStepCenterRenewed= false;
858 
859  for( num=0; num< fNoActiveNavigators; ++num )
860  {
862  safety = std::max( 0.0, fPreSafetyValues[num] - MagShift);
863  minSafety= std::min( safety, minSafety );
864  fCurrentPreStepSafety[num]= safety;
865  }
866  minStep= kInfinity;
867 
868 #ifdef G4DEBUG_PATHFINDER
869  if( fVerboseLevel > 2 )
870  {
871  G4cout << "G4ITPathFinder::DoNextLinearStep : Quick Stepping. " << G4endl
872  << " proposedStepLength " << proposedStepLength
873  << " < (full) safety = " << fullSafety
874  << " at " << initialPosition
875  << G4endl;
876  }
877 #endif
878  }
879  else
880 #endif // End of G4PATHFINDER_OPTIMISATION 1
881  {
882  // Move is larger than at least one of the safeties
883  // -> so we must move the safety center!
884 
885  fPreStepCenterRenewed= true;
886  pNavigatorIter= fpTransportManager-> GetActiveNavigatorsIterator();
887 
888  minStep= kInfinity; // Not proposedStepLength;
889 
890  for( num=0; num< fNoActiveNavigators; ++pNavigatorIter,++num )
891  {
892  safety = std::max( 0.0, fPreSafetyValues[num] - MagShift);
893 
894 #ifdef G4PATHFINDER_OPTIMISATION
895  if( proposedStepLength <= safety ) // Should be just < safety ?
896  {
897  // The Step is guaranteed to be taken
898 
899  step= kInfinity; // ComputeStep Would return this
900 
901 #ifdef G4DEBUG_PATHFINDER
902  G4cout.precision(8);
903  G4cout << "G4ITNavigator::ComputeStep> small proposed step = "
904  << proposedStepLength
905  << " <= safety = " << safety << " for nav " << num
906  << " Step fully taken. " << G4endl;
907 #endif
908  }
909  else
910 #endif // End of G4PATHFINDER_OPTIMISATION 2
911  {
912 #ifdef G4DEBUG_PATHFINDER
913  G4double previousSafety= safety;
914 #endif
915  step= (*pNavigatorIter)->ComputeStep( initialPosition,
916  initialDirection,
917  proposedStepLength,
918  safety );
919  minStep = std::min( step, minStep);
920 
921  // TODO: consider whether/how to reduce the proposed step
922  // to the latest minStep value - to reduce calculations
923 
924 #ifdef G4DEBUG_PATHFINDER
925  if( fVerboseLevel > 0)
926  {
927  G4cout.precision(8);
928  G4cout << "G4ITNavigator::ComputeStep> long proposed step = "
929  << proposedStepLength
930  << " > safety = " << previousSafety
931  << " for nav " << num
932  << " . New safety = " << safety << " step= " << step
933  << G4endl;
934  }
935 #endif
936  }
937  fCurrentStepSize[num] = step;
938 
939  // Save safety value, must be done for all geometries "together"
940  // (even if not recomputed using call to ComputeStep)
941  // since they share the fPreSafetyLocation
942 
943  fPreSafetyValues[num]= safety;
944  fCurrentPreStepSafety[num]= safety;
945 
946  minSafety= std::min( safety, minSafety );
947 
948 #ifdef G4DEBUG_PATHFINDER
949  if( fVerboseLevel > 2 )
950  {
951  G4cout << "G4ITPathFinder::DoNextLinearStep : Navigator ["
952  << num << "] -- step size " << step << G4endl;
953  }
954 #endif
955  }
956 
957  // Only change these when safety is recalculated
958  // it is good/relevant only for safety calculations
959 
960  fPreSafetyLocation= initialPosition;
961  fPreSafetyMinValue= minSafety;
962  } // end of else for if( proposedStepLength <= fullSafety)
963 
964  // For use in Relocation, need PreStep point location, min-safety
965  //
966  fPreStepLocation= initialPosition;
967  fMinSafety_PreStepPt= minSafety;
968 
969  fMinStep= minStep;
970 
971  if( fMinStep == kInfinity )
972  {
973  minStep = proposedStepLength; // Use this below for endpoint !!
974  }
975  fTrueMinStep = minStep;
976 
977  // Set the EndState
978 
979  G4ThreeVector endPosition;
980 
981  fEndState= initialState;
982  endPosition= initialPosition + minStep * initialDirection ;
983 
984 #ifdef G4DEBUG_PATHFINDER
985  if( fVerboseLevel > 1 )
986  {
987  G4cout << "G4ITPathFinder::DoNextLinearStep : "
988  << " initialPosition = " << initialPosition
989  << " and endPosition = " << endPosition<< G4endl;
990  }
991 #endif
992 
993  fEndState.SetPosition( endPosition );
994  fEndState.SetProperTimeOfFlight( -1.000 ); // Not defined YET
995 
996  if( fNoActiveNavigators == 1 )
997  {
998  G4bool transportLimited = (fMinStep!= kInfinity);
999  fLimitTruth[IdTransport] = transportLimited;
1000  fLimitedStep[IdTransport] = transportLimited ? kUnique : kDoNot;
1001 
1002  // Set fNoGeometriesLimiting - as WhichLimited does
1003  fNoGeometriesLimiting = transportLimited ? 1 : 0;
1004  }
1005  else
1006  {
1007  WhichLimited();
1008  }
1009 
1010 #ifdef G4DEBUG_PATHFINDER
1011  if( fVerboseLevel > 2 )
1012  {
1013  G4cout << " G4ITPathFinder::DoNextLinearStep : exits returning "
1014  << minStep << G4endl;
1015  G4cout << " Endpoint values = " << fEndState << G4endl;
1016  G4cout << G4endl;
1017  }
1018 #endif
1019 
1020  return minStep;
1021 }
1022 
1024 {
1025  // Flag which processes limited the step
1026 
1027  G4int num=-1, last=-1;
1028  G4int noLimited=0;
1029  ELimited shared= kSharedOther;
1030 
1031  const G4int IdTransport= 0; // Id of Mass Navigator !!
1032 
1033  // Assume that [IdTransport] is Mass / Transport
1034  //
1035  G4bool transportLimited = (fCurrentStepSize[IdTransport] == fMinStep)
1036  && ( fMinStep!= kInfinity) ;
1037 
1038  if( transportLimited ) {
1039  shared= kSharedTransport;
1040  }
1041 
1042  for ( num= 0; num < fNoActiveNavigators; num++ )
1043  {
1044  G4bool limitedStep;
1045 
1046  G4double step= fCurrentStepSize[num];
1047 
1048  limitedStep = ( std::fabs(step - fMinStep) < kCarTolerance )
1049  && ( step != kInfinity);
1050 
1051  fLimitTruth[ num ] = limitedStep;
1052  if( limitedStep )
1053  {
1054  noLimited++;
1055  fLimitedStep[num] = shared;
1056  last= num;
1057  }
1058  else
1059  {
1060  fLimitedStep[num] = kDoNot;
1061  }
1062  }
1063  fNoGeometriesLimiting= noLimited; // Save # processes limiting step
1064 
1065  if( (last > -1) && (noLimited == 1 ) )
1066  {
1067  fLimitedStep[ last ] = kUnique;
1068  }
1069 
1070 #ifdef G4DEBUG_PATHFINDER
1071  if( fVerboseLevel > 1 )
1072  {
1073  PrintLimited(); // --> for tracing
1074  if( fVerboseLevel > 4 ) {
1075  G4cout << " G4ITPathFinder::WhichLimited - exiting. " << G4endl;
1076  }
1077  }
1078 #endif
1079 }
1080 
1082 {
1083  // Report results -- for checking
1084 
1085  G4cout << "G4ITPathFinder::PrintLimited reports: " ;
1086  G4cout << " Minimum step (true)= " << fTrueMinStep
1087  << " reported min = " << fMinStep
1088  << G4endl;
1089  if( (fCurrentStepNo <= 2) || (fVerboseLevel>=2) )
1090  {
1091  G4cout << std::setw(5) << " Step#" << " "
1092  << std::setw(5) << " NavId" << " "
1093  << std::setw(12) << " step-size " << " "
1094  << std::setw(12) << " raw-size " << " "
1095  << std::setw(12) << " pre-safety " << " "
1096  << std::setw(15) << " Limited / flag" << " "
1097  << std::setw(15) << " World " << " "
1098  << G4endl;
1099  }
1100  G4int num;
1101  for ( num= 0; num < fNoActiveNavigators; num++ )
1102  {
1103  G4double rawStep = fCurrentStepSize[num];
1104  G4double stepLen = fCurrentStepSize[num];
1105  if( stepLen > fTrueMinStep )
1106  {
1107  stepLen = fTrueMinStep; // did not limit (went as far as asked)
1108  }
1109  G4int oldPrec= G4cout.precision(9);
1110 
1111  G4cout << std::setw(5) << fCurrentStepNo << " "
1112  << std::setw(5) << num << " "
1113  << std::setw(12) << stepLen << " "
1114  << std::setw(12) << rawStep << " "
1115  << std::setw(12) << fCurrentPreStepSafety[num] << " "
1116  << std::setw(5) << (fLimitTruth[num] ? "YES" : " NO") << " ";
1117  G4String limitedStr= LimitedString(fLimitedStep[num]);
1118  G4cout << " " << std::setw(15) << limitedStr << " ";
1119  G4cout.precision(oldPrec);
1120 
1121  G4ITNavigator *pNav= GetNavigator( num );
1122  G4String WorldName( "Not-Set" );
1123  if (pNav)
1124  {
1125  G4VPhysicalVolume *pWorld= pNav->GetWorldVolume();
1126  if( pWorld )
1127  {
1128  WorldName = pWorld->GetName();
1129  }
1130  }
1131  G4cout << " " << WorldName ;
1132  G4cout << G4endl;
1133  }
1134 
1135  if( fVerboseLevel > 4 )
1136  {
1137  G4cout << " G4ITPathFinder::PrintLimited - exiting. " << G4endl;
1138  }
1139 }
1140 
1141 G4double
1143  G4double proposedStepLength,
1144  G4VPhysicalVolume* /*pCurrentPhysicalVolume*/ )
1145 {
1146  const G4double toleratedRelativeError= 1.0e-10;
1147  G4double minStep= kInfinity, newSafety=0.0;
1148  G4int numNav;
1149  G4FieldTrack fieldTrack= initialState;
1150  G4ThreeVector startPoint= initialState.GetPosition();
1151 
1152 #ifdef G4DEBUG_PATHFINDER
1153  G4int prc= G4cout.precision(9);
1154  if( fVerboseLevel > 2 )
1155  {
1156  G4cout << " G4ITPathFinder::DoNextCurvedStep ****** " << G4endl;
1157  G4cout << " Initial value of field track is " << fieldTrack
1158  << " and proposed step= " << proposedStepLength << G4endl;
1159  }
1160 #endif
1161 
1162  fPreStepCenterRenewed= true; // Always update PreSafety with PreStep point
1163 
1164  if( fNoActiveNavigators > 1 )
1165  {
1166  // Calculate the safety values before making the step
1167 
1168  G4double minSafety= kInfinity, safety;
1169  for( numNav=0; numNav < fNoActiveNavigators; ++numNav )
1170  {
1171  safety= fpNavigator[numNav]->ComputeSafety( startPoint, false );
1172  fPreSafetyValues[numNav]= safety;
1173  fCurrentPreStepSafety[numNav]= safety;
1174  minSafety = std::min( safety, minSafety );
1175  }
1176 
1177  // Save safety value, related position
1178 
1179  fPreSafetyLocation= startPoint;
1180  fPreSafetyMinValue= minSafety;
1181  fPreStepLocation= startPoint;
1182  fMinSafety_PreStepPt= minSafety;
1183  }
1184 
1185  /*
1186  // Allow Propagator In Field to do the hard work, calling G4MultiNavigator
1187  //
1188  minStep= fpFieldPropagator->ComputeStep( fieldTrack,
1189  proposedStepLength,
1190  newSafety,
1191  pCurrentPhysicalVolume );
1192 */
1193  // fieldTrack now contains the endpoint information
1194  //
1195  fEndState= fieldTrack;
1196  fMinStep= minStep;
1197  fTrueMinStep = std::min( minStep, proposedStepLength );
1198 
1199  if( fNoActiveNavigators== 1 )
1200  {
1201  // Update the 'PreSafety' sphere - as any ComputeStep was called
1202  // (must be done anyway in field)
1203 
1204  fPreSafetyValues[0]= newSafety;
1205  fPreSafetyLocation= startPoint;
1206  fPreSafetyMinValue= newSafety;
1207 
1208  // Update the current 'PreStep' point's values - mandatory
1209  //
1210  fCurrentPreStepSafety[0]= newSafety;
1211  fPreStepLocation= startPoint;
1212  fMinSafety_PreStepPt= newSafety;
1213  }
1214 
1215 #ifdef G4DEBUG_PATHFINDER
1216  if( fVerboseLevel > 2 )
1217  {
1218  G4cout << "G4ITPathFinder::DoNextCurvedStep : " << G4endl
1219  << " initialState = " << initialState << G4endl
1220  << " and endState = " << fEndState << G4endl;
1221  G4cout << "G4ITPathFinder::DoNextCurvedStep : "
1222  << " minStep = " << minStep
1223  << " proposedStepLength " << proposedStepLength
1224  << " safety = " << newSafety << G4endl;
1225  }
1226 #endif
1227  G4double currentStepSize; // = 0.0;
1228  if( minStep < proposedStepLength ) // if == , then a boundary found at end ??
1229  {
1230  // Recover the remaining information from MultiNavigator
1231  // especially regarding which Navigator limited the step
1232 
1233  G4int noLimited= 0; // No geometries limiting step
1234  for( numNav=0; numNav < fNoActiveNavigators; ++numNav )
1235  {
1236  G4double finalStep, lastPreSafety=0.0, minStepLast;
1237  ELimited didLimit;
1238  G4bool limited;
1239 
1240  finalStep= fpMultiNavigator->ObtainFinalStep( numNav, lastPreSafety,
1241  minStepLast, didLimit );
1242 
1243  // Calculate the step for this geometry, using the
1244  // final step (the only one which can differ.)
1245 
1246  currentStepSize = fTrueMinStep;
1247  G4double diffStep= 0.0;
1248  if( (minStepLast != kInfinity) )
1249  {
1250  diffStep = (finalStep-minStepLast);
1251  if ( std::abs(diffStep) <= toleratedRelativeError * finalStep )
1252  {
1253  diffStep = 0.0;
1254  }
1255  currentStepSize += diffStep;
1256  }
1257  fCurrentStepSize[numNav] = currentStepSize;
1258 
1259  // TODO: could refine the way to obtain safeties for > 1 geometries
1260  // - for pre step safety
1261  // notify MultiNavigator about new set of sub-steps
1262  // allow it to return this value in ObtainFinalStep
1263  // instead of lastPreSafety (or as well?)
1264  // - for final step start (available)
1265  // get final Step start from MultiNavigator
1266  // and corresponding safety values
1267  // and/or ALSO calculate ComputeSafety at endpoint
1268  // endSafety= fpNavigator[numNav]->ComputeSafety( endPoint );
1269 
1270  fLimitedStep[numNav] = didLimit;
1271  fLimitTruth[numNav] = limited = (didLimit != kDoNot );
1272  if( limited ) { noLimited++; }
1273 
1274 #ifdef G4DEBUG_PATHFINDER
1275  G4bool StepError= (currentStepSize < 0)
1276  || ( (minStepLast != kInfinity) && (diffStep < 0) ) ;
1277  if( StepError || (fVerboseLevel > 2) )
1278  {
1279  G4String limitedString= LimitedString( fLimitedStep[numNav] );
1280 
1281  G4cout << " G4ITPathFinder::ComputeStep. Geometry " << numNav
1282  << " step= " << fCurrentStepSize[numNav]
1283  << " from final-step= " << finalStep
1284  << " fTrueMinStep= " << fTrueMinStep
1285  << " minStepLast= " << minStepLast
1286  << " limited = " << (fLimitTruth[numNav] ? "YES" : " NO")
1287  << " ";
1288  G4cout << " status = " << limitedString << " #= " << didLimit
1289  << G4endl;
1290 
1291  if( StepError )
1292  {
1293  std::ostringstream message;
1294  message << "Incorrect calculation of step size for one navigator"
1295  << G4endl
1296  << " currentStepSize = " << currentStepSize
1297  << ", diffStep= " << diffStep << G4endl
1298  << "ERROR in computing step size for this navigator.";
1299  G4Exception("G4ITPathFinder::DoNextCurvedStep",
1300  "GeomNav0003", FatalException, message);
1301  }
1302  }
1303 #endif
1304  } // for num Navigators
1305 
1306  fNoGeometriesLimiting= noLimited; // Save # processes limiting step
1307  }
1308  else if ( (minStep == proposedStepLength)
1309  || (minStep == kInfinity)
1310  || ( std::abs(minStep-proposedStepLength)
1311  < toleratedRelativeError * proposedStepLength ) )
1312  {
1313  // In case the step was not limited, use default responses
1314  // --> all Navigators
1315  // Also avoid problems in case of G4ITPathFinder using safety to optimise
1316  // - it is possible that the Navigators were not called
1317  // if the safety was already satisfactory.
1318  // (In that case calling ObtainFinalStep gives invalid results.)
1319 
1320  currentStepSize= minStep;
1321  for( numNav=0; numNav < fNoActiveNavigators; ++numNav )
1322  {
1323  fCurrentStepSize[numNav] = minStep;
1324  // Safety for endpoint ?? // Can eventuall improve it -- see TODO above
1325  fLimitedStep[numNav] = kDoNot;
1326  fLimitTruth[numNav] = false;
1327  }
1328  fNoGeometriesLimiting= 0; // Save # processes limiting step
1329  }
1330  else // (minStep > proposedStepLength) and not (minStep == kInfinity)
1331  {
1332  std::ostringstream message;
1333  message << "Incorrect calculation of step size for one navigator." << G4endl
1334  << " currentStepSize = " << minStep << " is larger than "
1335  << " proposed StepSize = " << proposedStepLength << ".";
1336  G4Exception("G4ITPathFinder::DoNextCurvedStep()",
1337  "GeomNav0003", FatalException, message);
1338  }
1339 
1340 #ifdef G4DEBUG_PATHFINDER
1341  if( fVerboseLevel > 2 )
1342  {
1343  G4cout << " Exiting G4ITPathFinder::DoNextCurvedStep " << G4endl;
1344  PrintLimited();
1345  }
1346  G4cout.precision(prc);
1347 #endif
1348 
1349  return minStep;
1350 }
1351 
1353 {
1354  static G4String StrDoNot("DoNot"),
1355  StrUnique("Unique"),
1356  StrUndefined("Undefined"),
1357  StrSharedTransport("SharedTransport"),
1358  StrSharedOther("SharedOther");
1359 
1360  G4String* limitedStr;
1361  switch ( lim )
1362  {
1363  case kDoNot: limitedStr= &StrDoNot; break;
1364  case kUnique: limitedStr = &StrUnique; break;
1365  case kSharedTransport: limitedStr= &StrSharedTransport; break;
1366  case kSharedOther: limitedStr = &StrSharedOther; break;
1367  default: limitedStr = &StrUndefined; break;
1368  }
1369  return *limitedStr;
1370 }
1371 
1373 {
1376  for( G4int nav=0; nav < fNoActiveNavigators; ++nav )
1377  {
1379  }
1380 }
1381 
G4ITNavigator * GetNavigator(G4int n) const
G4ITSafetyHelper * GetSafetyHelper() const
void ReLocate(const G4ThreeVector &position)
static G4ITPathFinder * GetInstance()
#define fCurrentPreStepSafety
#define fCurrentStepSize
static constexpr double mm
Definition: G4SIunits.hh:115
static const G4double kInfinity
Definition: geomdefs.hh:42
#define fNoGeometriesLimiting
#define fNewTrack
double x() const
#define fPreSafetyValues
void EnableParallelNavigation(G4bool parallel)
void EnableParallelNavigation(G4bool enableChoice=true)
void Locate(const G4ThreeVector &position, const G4ThreeVector &direction, G4bool relativeSearch=true)
G4TouchableHandle CreateTouchableHandle(G4int navId) const
ELimited
G4double GetSurfaceTolerance() const
void PrepareNewTrack(const G4ThreeVector &position, const G4ThreeVector &direction, G4VPhysicalVolume *massStartVol=0)
#define fMinStep
#define G4ThreadLocal
Definition: tls.hh:89
int G4int
Definition: G4Types.hh:78
#define fPreStepCenterRenewed
double z() const
G4ReferenceCountedHandle< G4VTouchable > G4TouchableHandle
#define fLocatedVolume
#define fPreSafetyMinValue
static G4ITTransportationManager * GetTransportationManager()
G4double DoNextLinearStep(const G4FieldTrack &FieldTrack, G4double proposedStepLength)
#define position
Definition: xmlparse.cc:622
G4ThreeVector GetPosition() const
G4GLOB_DLL std::ostream G4cout
const G4String & GetName() const
G4double GetRadialTolerance() const
bool G4bool
Definition: G4Types.hh:79
void ReportMove(const G4ThreeVector &OldV, const G4ThreeVector &NewV, const G4String &Quantity) const
#define fLastLocatedPosition
G4double ComputeSafety(const G4ThreeVector &globalPoint)
#define fMinSafety_atSafLocation
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4double ObtainFinalStep(G4int navigatorId, G4double &pNewSafety, G4double &minStepLast, ELimited &limitedStep)
G4double DoNextCurvedStep(const G4FieldTrack &FieldTrack, G4double proposedStepLength, G4VPhysicalVolume *pCurrentPhysVolume)
#define fCurrentStepNo
std::vector< G4ITNavigator * >::iterator GetActiveNavigatorsIterator()
#define fMinSafety_PreStepPt
G4String & LimitedString(ELimited lim)
#define fLimitTruth
T max(const T t1, const T t2)
brief Return the largest of the two arguments
#define fPreSafetyLocation
#define fLimitedStep
double y() const
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
double mag2() const
#define G4endl
Definition: G4ios.hh:61
#define fPreStepLocation
void UpdateYourself(G4VPhysicalVolume *pPhysVol, const G4NavigationHistory *history=0)
#define fSafetyLocation
#define fLastStepNo
const G4NavigationHistory * GetHistory() const
T sqr(const T &x)
Definition: templates.hh:145
#define fEndState
G4double ComputeStep(const G4FieldTrack &pFieldTrack, G4double pCurrentProposedStepLength, G4int navigatorId, G4int stepNo, G4double &pNewSafety, ELimited &limitedStep, G4FieldTrack &EndState, G4VPhysicalVolume *currentVolume)
double G4double
Definition: G4Types.hh:76
#define fRelocatedPoint
double mag() const
void PushPostSafetyToPreSafety()
static G4GeometryTolerance * GetInstance()
#define fTrueMinStep
G4GLOB_DLL std::ostream G4cerr
G4ThreeVector GetMomentumDirection() const
#define fNewSafetyComputed