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