Geant4-11
G4VIntersectionLocator.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// Class G4VIntersectionLocator implementation
27//
28// 27.10.08 - John Apostolakis, Tatiana Nikitina.
29// ---------------------------------------------------------------------------
30
31#include <iomanip>
32#include <sstream>
33
34#include "globals.hh"
35#include "G4ios.hh"
36#include "G4AutoDelete.hh"
37#include "G4SystemOfUnits.hh"
40
42//
43// Constructor
44//
46 : fiNavigator(theNavigator)
47{
49
50 if( fiNavigator->GetExternalNavigation() == nullptr )
51 {
53 }
54 else // Must clone the navigator, together with External Navigation
55 {
57 }
58}
59
61//
62// Destructor.
63//
65{
66 delete fHelpingNavigator;
67 delete fpTouchable;
68}
69
71//
72// Dump status of propagator to cout (old method)
73//
74void
76 const G4FieldTrack& CurrentFT,
77 G4double requestStep,
78 G4double safety,
79 G4int stepNo)
80{
81 std::ostringstream os;
82 printStatus( StartFT,CurrentFT,requestStep,safety,stepNo,os,fVerboseLevel);
83 G4cout << os.str();
84}
85
87//
88// Dumps status of propagator.
89//
90void
92 const G4FieldTrack& CurrentFT,
93 G4double requestStep,
94 G4double safety,
95 G4int stepNo,
96 std::ostream& os,
97 G4int verboseLevel)
98{
99 // const G4int verboseLevel= fVerboseLevel;
100 const G4ThreeVector StartPosition = StartFT.GetPosition();
101 const G4ThreeVector StartUnitVelocity = StartFT.GetMomentumDir();
102 const G4ThreeVector CurrentPosition = CurrentFT.GetPosition();
103 const G4ThreeVector CurrentUnitVelocity = CurrentFT.GetMomentumDir();
104
105 G4double step_len = CurrentFT.GetCurveLength() - StartFT.GetCurveLength();
106 G4int oldprc; // cout/cerr precision settings
107
108 if( ((stepNo == 0) && (verboseLevel <3)) || (verboseLevel >= 3) )
109 {
110 oldprc = os.precision(4);
111 os << std::setw( 6) << " "
112 << std::setw( 25) << " Current Position and Direction" << " "
113 << G4endl;
114 os << std::setw( 5) << "Step#"
115 << std::setw(10) << " s " << " "
116 << std::setw(10) << "X(mm)" << " "
117 << std::setw(10) << "Y(mm)" << " "
118 << std::setw(10) << "Z(mm)" << " "
119 << std::setw( 7) << " N_x " << " "
120 << std::setw( 7) << " N_y " << " "
121 << std::setw( 7) << " N_z " << " " ;
122 os << std::setw( 7) << " Delta|N|" << " "
123 << std::setw( 9) << "StepLen" << " "
124 << std::setw(12) << "StartSafety" << " "
125 << std::setw( 9) << "PhsStep" << " ";
126 os << G4endl;
127 os.precision(oldprc);
128 }
129 if((stepNo == 0) && (verboseLevel <=3))
130 {
131 // Recurse to print the start values
132 //
133 printStatus( StartFT, StartFT, -1.0, safety, -1, os, verboseLevel);
134 }
135 if( verboseLevel <= 3 )
136 {
137 if( stepNo >= 0)
138 {
139 os << std::setw( 4) << stepNo << " ";
140 }
141 else
142 {
143 os << std::setw( 5) << "Start" ;
144 }
145 oldprc = os.precision(8);
146 os << std::setw(10) << CurrentFT.GetCurveLength() << " ";
147 os << std::setw(10) << CurrentPosition.x() << " "
148 << std::setw(10) << CurrentPosition.y() << " "
149 << std::setw(10) << CurrentPosition.z() << " ";
150 os.precision(4);
151 os << std::setw( 7) << CurrentUnitVelocity.x() << " "
152 << std::setw( 7) << CurrentUnitVelocity.y() << " "
153 << std::setw( 7) << CurrentUnitVelocity.z() << " ";
154 os.precision(3);
155 os << std::setw( 7)
156 << CurrentFT.GetMomentum().mag()- StartFT.GetMomentum().mag()
157 << " ";
158 os << std::setw( 9) << step_len << " ";
159 os << std::setw(12) << safety << " ";
160 if( requestStep != -1.0 )
161 {
162 os << std::setw( 9) << requestStep << " ";
163 }
164 else
165 {
166 os << std::setw( 9) << "Init/NotKnown" << " ";
167 }
168 os << G4endl;
169 os.precision(oldprc);
170 }
171 else // if( verboseLevel > 3 )
172 {
173 // Multi-line output
174
175 os << "Step taken was " << step_len
176 << " out of PhysicalStep= " << requestStep << G4endl;
177 os << "Final safety is: " << safety << G4endl;
178 os << "Chord length = " << (CurrentPosition-StartPosition).mag()
179 << G4endl;
180 os << G4endl;
181 }
182}
183
185//
186// ReEstimateEndPoint.
187//
189ReEstimateEndpoint( const G4FieldTrack& CurrentStateA,
190 const G4FieldTrack& EstimatedEndStateB,
191 G4double , // linearDistSq, // NOT used
193#ifdef G4DEBUG_FIELD
194 curveDist
195#endif
196 )
197{
198 G4FieldTrack newEndPoint( CurrentStateA );
199 auto integrDriver = GetChordFinderFor()->GetIntegrationDriver();
200
201 G4FieldTrack retEndPoint( CurrentStateA );
202 G4bool goodAdvance;
203 G4int itrial = 0;
204 const G4int no_trials = 20;
205
206
207 G4double endCurveLen= EstimatedEndStateB.GetCurveLength();
208
209 do // Loop checking, 07.10.2016, JA
210 {
211 G4double currentCurveLen = newEndPoint.GetCurveLength();
212 G4double advanceLength = endCurveLen - currentCurveLen ;
213 if (std::abs(advanceLength)<kCarTolerance)
214 {
215 goodAdvance=true;
216 }
217 else
218 {
219 goodAdvance = integrDriver->AccurateAdvance(newEndPoint, advanceLength,
221 }
222 }
223 while( !goodAdvance && (++itrial < no_trials) );
224
225 if( goodAdvance )
226 {
227 retEndPoint = newEndPoint;
228 }
229 else
230 {
231 retEndPoint = EstimatedEndStateB; // Could not improve without major work !!
232 }
233
234 // All the work is done
235 // below are some diagnostics only -- before the return!
236 //
237 const G4String MethodName("G4VIntersectionLocator::ReEstimateEndpoint()");
238
239#ifdef G4VERBOSE
240 G4int latest_good_trials = 0;
241 if( itrial > 1)
242 {
243 if( fVerboseLevel > 0 )
244 {
245 G4cout << MethodName << " called - goodAdv= " << goodAdvance
246 << " trials = " << itrial
247 << " previous good= " << latest_good_trials
248 << G4endl;
249 }
250 latest_good_trials = 0;
251 }
252 else
253 {
254 ++latest_good_trials;
255 }
256#endif
257
258#ifdef G4DEBUG_FIELD
259 G4double lengthDone = newEndPoint.GetCurveLength()
260 - CurrentStateA.GetCurveLength();
261 if( !goodAdvance )
262 {
263 if( fVerboseLevel >= 3 )
264 {
265 G4cout << MethodName << "> AccurateAdvance failed " ;
266 G4cout << " in " << itrial << " integration trials/steps. " << G4endl;
267 G4cout << " It went only " << lengthDone << " instead of " << curveDist
268 << " -- a difference of " << curveDist - lengthDone << G4endl;
269 G4cout << " ReEstimateEndpoint> Reset endPoint to original value!"
270 << G4endl;
271 }
272 }
273 G4double linearDist = ( EstimatedEndStateB.GetPosition()
274 - CurrentStateA.GetPosition() ).mag();
275 static G4int noInaccuracyWarnings = 0;
276 G4int maxNoWarnings = 10;
277 if ( (noInaccuracyWarnings < maxNoWarnings )
278 || (fVerboseLevel > 1) )
279 {
280 G4ThreeVector move = newEndPoint.GetPosition()
281 - EstimatedEndStateB.GetPosition();
282 std::ostringstream message;
283 message.precision(12);
284 message << " Integration inaccuracy requires"
285 << " an adjustment in the step's endpoint." << G4endl
286 << " Two mid-points are further apart than their"
287 << " curve length difference" << G4endl
288 << " Dist = " << linearDist
289 << " curve length = " << curveDist << G4endl;
290 message << " Correction applied is " << move.mag() << G4endl
291 << " Old Estimated B position= "
292 << EstimatedEndStateB.GetPosition() << G4endl
293 << " Recalculated Position= "
294 << newEndPoint.GetPosition() << G4endl
295 << " Change ( new - old ) = " << move;
296 G4Exception("G4VIntersectionLocator::ReEstimateEndpoint()",
297 "GeomNav1002", JustWarning, message);
298 }
299/*
300#else
301 // Statistics on the RMS value of the corrections
302
303 static G4ThreadLocal G4int noCorrections = 0;
304 ++noCorrections;
305 if( goodAdvance )
306 {
307 static G4ThreadLocal G4double sumCorrectionsSq;
308 sumCorrectionsSq += (EstimatedEndStateB.GetPosition() -
309 newEndPoint.GetPosition()).mag2();
310 }
311*/
312#endif
313
314 return retEndPoint;
315}
316
317
319//
320// ReEstimateEndPoint.
321//
322// The following values are returned in curveError
323// 0: Normal - no problem
324// 1: Unexpected co-incidence - milder mixup
325// 2: Real mixup - EndB is NOT past StartA
326// ( ie. StartA's curve-lengh is bigger than EndB's)
327
328
330CheckAndReEstimateEndpoint( const G4FieldTrack& CurrentStartA,
331 const G4FieldTrack& EstimatedEndB,
332 G4FieldTrack& RevisedEndPoint,
333 G4int& curveError)
334{
335 G4double linDistSq, curveDist;
336
337 G4bool recalculated = false;
338 curveError= 0;
339
340 linDistSq = ( EstimatedEndB.GetPosition()
341 - CurrentStartA.GetPosition() ).mag2();
342 curveDist = EstimatedEndB.GetCurveLength()
343 - CurrentStartA.GetCurveLength();
344 if( (curveDist>=0.0)
345 && (curveDist*curveDist *(1.0+2.0*fiEpsilonStep ) < linDistSq ) )
346 {
347 G4FieldTrack newEndPointFT = EstimatedEndB; // Unused
348
349 if (curveDist>0.0)
350 {
351 // Re-integrate to obtain a new B
352 RevisedEndPoint = ReEstimateEndpoint( CurrentStartA,
353 EstimatedEndB,
354 linDistSq,
355 curveDist );
356 recalculated = true;
357 }
358 else
359 {
360 // Zero length -> no advance!
361 newEndPointFT = CurrentStartA;
362 recalculated = true;
363 curveError = 1; // Unexpected co-incidence - milder mixup
364
365 G4Exception("G4MultiLevelLocator::EstimateIntersectionPoint()",
366 "GeomNav1002", JustWarning,
367 "A & B are at equal distance in 2nd half. A & B will coincide." );
368 }
369 }
370
371 // Sanity check
372 //
373 if( curveDist < 0.0 )
374 {
375 curveError = 2; // Real mixup
376 }
377 return recalculated;
378}
379
381//
382// Method for finding SurfaceNormal of Intersecting Solid
383//
385GetLocalSurfaceNormal(const G4ThreeVector& CurrentE_Point, G4bool& validNormal)
386{
387 G4ThreeVector Normal(G4ThreeVector(0.0,0.0,0.0));
388 G4VPhysicalVolume* located;
389
390 validNormal = false;
392 located = fHelpingNavigator->LocateGlobalPointAndSetup( CurrentE_Point );
393
394 delete fpTouchable;
396
397 // To check if we can use GetGlobalExitNormal()
398 //
399 G4ThreeVector localPosition = fpTouchable->GetHistory()
400 ->GetTopTransform().TransformPoint(CurrentE_Point);
401
402 // Issue: in the case of coincident surfaces, this version does not recognise
403 // which side you are located onto (can return vector with wrong sign.)
404 // TO-DO: use direction (of chord) to identify volume we will be "entering"
405
406 if( located != 0)
407 {
408 G4LogicalVolume* pLogical= located->GetLogicalVolume();
409 G4VSolid* pSolid;
410
411 if( (pLogical != nullptr) && ( (pSolid=pLogical->GetSolid()) != nullptr ) )
412 {
413 if ( ( pSolid->Inside(localPosition)==kSurface )
414 || ( pSolid->DistanceToOut(localPosition) < 1000.0 * kCarTolerance ) )
415 {
416 Normal = pSolid->SurfaceNormal(localPosition);
417 validNormal = true;
418
419#ifdef G4DEBUG_FIELD
420 if( std::fabs(Normal.mag2() - 1.0 ) > CLHEP::perThousand)
421 {
422 G4cerr << "PROBLEM in G4VIntersectionLocator::GetLocalSurfaceNormal."
423 << G4endl;
424 G4cerr << " Normal is not unit - mag=" << Normal.mag() << G4endl;
425 G4cerr << " at trial local point " << CurrentE_Point << G4endl;
426 G4cerr << " Solid is " << *pSolid << G4endl;
427 }
428#endif
429 }
430 }
431 }
432 return Normal;
433}
434
436//
437// Adjustment of Found Intersection
438//
440AdjustmentOfFoundIntersection( const G4ThreeVector& CurrentA_Point,
441 const G4ThreeVector& CurrentE_Point,
442 const G4ThreeVector& CurrentF_Point,
443 const G4ThreeVector& MomentumDir,
444 const G4bool IntersectAF,
445 G4ThreeVector& IntersectionPoint, // I/O
446 G4double& NewSafety, // I/O
449{
450 G4double dist,lambda;
451 G4ThreeVector Normal, NewPoint, Point_G;
452 G4bool goodAdjust = false, Intersects_FP = false, validNormal = false;
453
454 // Get SurfaceNormal of Intersecting Solid
455 //
456 Normal = GetGlobalSurfaceNormal(CurrentE_Point,validNormal);
457 if(!validNormal) { return false; }
458
459 // Intersection between Line and Plane
460 //
461 G4double n_d_m = Normal.dot(MomentumDir);
462 if ( std::abs(n_d_m)>kCarTolerance )
463 {
464#ifdef G4VERBOSE
465 if ( fVerboseLevel>1 )
466 {
467 G4Exception("G4VIntersectionLocator::AdjustmentOfFoundIntersection()",
468 "GeomNav0003", JustWarning,
469 "No intersection. Parallels lines!");
470 }
471#endif
472 lambda =- Normal.dot(CurrentF_Point-CurrentE_Point)/n_d_m;
473
474 // New candidate for Intersection
475 //
476 NewPoint = CurrentF_Point+lambda*MomentumDir;
477
478 // Distance from CurrentF to Calculated Intersection
479 //
480 dist = std::abs(lambda);
481
482 if ( dist<kCarTolerance*0.001 ) { return false; }
483
484 // Calculation of new intersection point on the path.
485 //
486 if ( IntersectAF ) // First part intersects
487 {
488 G4double stepLengthFP;
489 G4ThreeVector Point_P = CurrentA_Point;
491 Intersects_FP = IntersectChord( Point_P, NewPoint, NewSafety,
493 stepLengthFP, Point_G );
494
495 }
496 else // Second part intersects
497 {
498 G4double stepLengthFP;
500 Intersects_FP = IntersectChord( CurrentF_Point, NewPoint, NewSafety,
502 stepLengthFP, Point_G );
503 }
504 if ( Intersects_FP )
505 {
506 goodAdjust = true;
507 IntersectionPoint = Point_G;
508 }
509 }
510
511 return goodAdjust;
512}
513
515//
516// GetSurfaceNormal.
517//
519GetSurfaceNormal(const G4ThreeVector& CurrentInt_Point,
520 G4bool& validNormal)
521{
522 G4ThreeVector NormalAtEntry; // ( -10. , -10., -10. );
523
524 G4ThreeVector NormalAtEntryLast, NormalAtEntryGlobal, diffNormals;
525 G4bool validNormalLast;
526
527 // Relies on a call to Navigator::ComputeStep in IntersectChord before
528 // this call
529 //
530 NormalAtEntryLast = GetLastSurfaceNormal( CurrentInt_Point, validNormalLast );
531 // May return valid=false in cases, including
532 // - if the candidate volume was not found (eg exiting world), or
533 // - a replica was involved -- determined the step size.
534 // (This list is not complete.)
535
536#ifdef G4DEBUG_FIELD
537 if ( validNormalLast
538 && ( std::fabs(NormalAtEntryLast.mag2() - 1.0) > perThousand ) )
539 {
540 std::ostringstream message;
541 message << "PROBLEM: Normal is not unit - magnitude = "
542 << NormalAtEntryLast.mag() << G4endl;
543 message << " at trial intersection point " << CurrentInt_Point << G4endl;
544 message << " Obtained from Get *Last* Surface Normal.";
545 G4Exception("G4VIntersectionLocator::GetSurfaceNormal()",
546 "GeomNav1002", JustWarning, message);
547 }
548#endif
549
550 if( validNormalLast )
551 {
552 NormalAtEntry = NormalAtEntryLast;
553 }
554 validNormal = validNormalLast;
555
556 return NormalAtEntry;
557}
558
560//
561// GetGlobalSurfaceNormal.
562//
564GetGlobalSurfaceNormal(const G4ThreeVector& CurrentE_Point,
565 G4bool& validNormal)
566{
567 G4ThreeVector localNormal = GetLocalSurfaceNormal(CurrentE_Point,validNormal);
568 G4AffineTransform localToGlobal = // Must use the same Navigator !!
570 G4ThreeVector globalNormal = localToGlobal.TransformAxis( localNormal );
571
572#ifdef G4DEBUG_FIELD
573 if( validNormal && ( std::fabs(globalNormal.mag2() - 1.0) > perThousand ) )
574 {
575 std::ostringstream message;
576 message << "**************************************************************"
577 << G4endl;
578 message << " Bad Normal in G4VIntersectionLocator::GetGlobalSurfaceNormal "
579 << G4endl;
580 message << " * Constituents: " << G4endl;
581 message << " Local Normal= " << localNormal << G4endl;
582 message << " Transform: " << G4endl
583 << " Net Translation= " << localToGlobal.NetTranslation()
584 << G4endl
585 << " Net Rotation = " << localToGlobal.NetRotation()
586 << G4endl;
587 message << " * Result: " << G4endl;
588 message << " Global Normal= " << localNormal << G4endl;
589 message << "**************************************************************";
590 G4Exception("G4VIntersectionLocator::GetGlobalSurfaceNormal()",
591 "GeomNav1002", JustWarning, message);
592 }
593#endif
594
595 return globalNormal;
596}
597
599//
600// GetLastSurfaceNormal.
601//
603GetLastSurfaceNormal( const G4ThreeVector& intersectPoint,
604 G4bool& normalIsValid) const
605{
606 G4ThreeVector normalVec;
607 G4bool validNorm;
608 normalVec = fiNavigator->GetGlobalExitNormal( intersectPoint, &validNorm );
609 normalIsValid = validNorm;
610
611 return normalVec;
612}
613
615//
616// ReportTrialStep.
617//
619 const G4ThreeVector& ChordAB_v,
620 const G4ThreeVector& ChordEF_v,
621 const G4ThreeVector& NewMomentumDir,
622 const G4ThreeVector& NormalAtEntry,
623 G4bool validNormal )
624{
625 G4double ABchord_length = ChordAB_v.mag();
626 G4double MomDir_dot_Norm = NewMomentumDir.dot( NormalAtEntry );
627 G4double MomDir_dot_ABchord;
628 MomDir_dot_ABchord = (1.0 / ABchord_length) * NewMomentumDir.dot( ChordAB_v );
629
630 std::ostringstream outStream;
631 outStream << std::setw(6) << " Step# "
632 << std::setw(17) << " |ChordEF|(mag)" << " "
633 << std::setw(18) << " uMomentum.Normal" << " "
634 << std::setw(18) << " uMomentum.ABdir " << " "
635 << std::setw(16) << " AB-dist " << " "
636 << " Chord Vector (EF) "
637 << G4endl;
638 outStream.precision(7);
639 outStream << " " << std::setw(5) << step_no
640 << " " << std::setw(18) << ChordEF_v.mag()
641 << " " << std::setw(18) << MomDir_dot_Norm
642 << " " << std::setw(18) << MomDir_dot_ABchord
643 << " " << std::setw(12) << ABchord_length
644 << " " << ChordEF_v
645 << G4endl;
646 outStream << " MomentumDir= " << " " << NewMomentumDir
647 << " Normal at Entry E= " << NormalAtEntry
648 << " AB chord = " << ChordAB_v
649 << G4endl;
650 G4cout << outStream.str();
651
652 if( ( std::fabs(NormalAtEntry.mag2() - 1.0) > perThousand ) )
653 {
654 std::ostringstream message;
655 message << "Normal is not unit - mag= " << NormalAtEntry.mag() << G4endl
656 << " ValidNormalAtE = " << validNormal;
657 G4Exception("G4VIntersectionLocator::ReportTrialStep()",
658 "GeomNav1002", JustWarning, message);
659 }
660 return;
661}
662
664//
665// LocateGlobalPointWithinVolumeAndCheck.
666//
667// Locate point using navigator: updates state of Navigator
668// By default, it assumes that the point is inside the current volume,
669// and returns true.
670// In check mode, checks whether the point is *inside* the volume.
671// If it is inside, it returns true
672// If not, issues a warning and returns false.
673//
676{
677 G4bool good = true;
679 const G4String
680 MethodName("G4VIntersectionLocator::LocateGlobalPointWithinVolumeAndCheck()");
681
682 if( fCheckMode )
683 {
684 G4bool navCheck= nav->IsCheckModeActive(); // Recover original value
685 nav->CheckMode(true);
686
687 // Identify the current volume
688
690 G4VPhysicalVolume* motherPhys = startTH->GetVolume();
691 G4VSolid* motherSolid = startTH->GetSolid();
693 G4int motherCopyNo = motherPhys->GetCopyNo();
694
695 // Let's check that the point is inside the current solid
696 G4ThreeVector localPosition = transform.TransformPoint(position);
697 EInside inMother = motherSolid->Inside( localPosition );
698 if( inMother != kInside )
699 {
700 std::ostringstream message;
701 message << "Position located "
702 << ( inMother == kSurface ? " on Surface " : " outside " )
703 << "expected volume" << G4endl
704 << " Safety (from Outside) = "
705 << motherSolid->DistanceToIn(localPosition);
706 G4Exception("G4VIntersectionLocator::LocateGlobalPointWithinVolumeAndCheck()",
707 "GeomNav1002", JustWarning, message);
708 }
709
710 // 1. Simple next step - quick relocation and check result.
711 // nav->LocateGlobalPointWithinVolume( position );
712
713 // 2. Full relocation - to cross-check answer !
715 if( (nextPhysical != motherPhys)
716 || (nextPhysical->GetCopyNo() != motherCopyNo )
717 )
718 {
719 G4Exception("G4VIntersectionLocator::LocateGlobalPointWithinVolumeAndCheck()",
720 "GeomNav1002", JustWarning,
721 "Position located outside expected volume.");
722 }
723 nav->CheckMode(navCheck); // Recover original value
724 }
725 else
726 {
728 }
729 return good;
730}
731
733//
734// LocateGlobalPointWithinVolumeCheckAndReport.
735//
738 const G4String& CodeLocationInfo,
739 G4int /* CheckMode */)
740{
741 // Save value of Check mode first
742 G4bool oldCheck = GetCheckMode();
743
745 if( !ok )
746 {
747 std::ostringstream message;
748 message << "Failed point location." << G4endl
749 << " Code Location info: " << CodeLocationInfo;
750 G4Exception("G4VIntersectionLocator::LocateGlobalPointWithinVolumeCheckAndReport()",
751 "GeomNav1002", JustWarning, message);
752 }
753
754 SetCheckMode( oldCheck );
755}
756
758//
759// ReportReversedPoints.
760//
762ReportReversedPoints( std::ostringstream& msg,
763 const G4FieldTrack& StartPointVel,
764 const G4FieldTrack& EndPointVel,
765 G4double NewSafety, G4double epsStep,
766 const G4FieldTrack& A_PtVel,
767 const G4FieldTrack& B_PtVel,
768 const G4FieldTrack& SubStart_PtVel,
769 const G4ThreeVector& E_Point,
770 const G4FieldTrack& ApproxIntersecPointV,
771 G4int substep_no, G4int substep_no_p, G4int depth )
772{
773 // Expect that 'msg' can hold the name of the calling method
774
775 // FieldTrack 'points' A and B have been tangled
776 // Whereas A should be before B, it is found that curveLen(B) < curveLen(A)
777 G4int verboseLevel= 5;
778 G4double curveDist = B_PtVel.GetCurveLength() - A_PtVel.GetCurveLength();
779 G4VIntersectionLocator::printStatus( A_PtVel, B_PtVel,
780 -1.0, NewSafety, substep_no, msg, verboseLevel );
781 msg << "Error in advancing propagation." << G4endl
782 << " The final curve point is NOT further along"
783 << " than the original!" << G4endl
784 << " Going *backwards* from len(A) = " << A_PtVel.GetCurveLength()
785 << " to len(B) = " << B_PtVel.GetCurveLength() << G4endl
786 << " Curve distance is " << curveDist / CLHEP::millimeter << " mm "
787 << G4endl
788 << " Point A' (start) is " << A_PtVel << G4endl
789 << " Point B' (end) is " << B_PtVel << G4endl;
790 msg << " fEpsStep= " << epsStep << G4endl << G4endl;
791
792 G4int oldprc = msg.precision(20);
793 msg << " In full precision, the position, momentum, E_kin, length, rest mass "
794 << " ... are: " << G4endl;
795 msg << " Point A[0] (Curve start) is " << StartPointVel << G4endl
796 << " Point S (Sub start) is " << SubStart_PtVel
797 << " Point A' (Current start) is " << A_PtVel << G4endl
798 << " Point E (Trial Point) is " << E_Point << G4endl
799 << " Point F (Intersection) is " << ApproxIntersecPointV << G4endl
800 << " Point B' (Current end) is " << B_PtVel << G4endl
801 << " Point B[0] (Curve end) is " << EndPointVel << G4endl
802 << G4endl
803 << " LocateIntersection parameters are : " << G4endl
804 << " Substep no (total) = " << substep_no << G4endl
805 << " Substep no = " << substep_no_p << " at depth= " << depth;
806 msg.precision(oldprc);
807}
808
810//
811// ReportProgress.
812//
814 const G4FieldTrack& StartPointVel,
815 const G4FieldTrack& EndPointVel,
816 G4int substep_no,
817 const G4FieldTrack& A_PtVel,
818 const G4FieldTrack& B_PtVel,
819 G4double safetyLast,
820 G4int depth )
821
822{
823 oss << "ReportProgress: Current status of intersection search: " << G4endl;
824 if( depth > 0 ) oss << " Depth= " << depth;
825 oss << " Substep no = " << substep_no << G4endl;
826 G4int verboseLevel = 5;
827 G4double safetyPrev = -1.0; // Add as argument ?
828
829 printStatus( StartPointVel, EndPointVel, -1.0, -1.0, -1,
830 oss, verboseLevel);
831 oss << " * Start and end-point of requested Step:" << G4endl;
832 oss << " ** State of point A: ";
833 printStatus( A_PtVel, A_PtVel, -1.0, safetyPrev, substep_no-1,
834 oss, verboseLevel);
835 oss << " ** State of point B: ";
836 printStatus( A_PtVel, B_PtVel, -1.0, safetyLast, substep_no,
837 oss, verboseLevel);
838}
839
841//
842// ReportImmediateHit.
843//
844void
846 const G4ThreeVector& StartPosition,
847 const G4ThreeVector& TrialPoint,
848 G4double tolerance,
849 unsigned long int numCalls )
850{
851 static G4ThreadLocal unsigned int occurredOnTop= 0;
852 static G4ThreadLocal G4ThreeVector* ptrLast = nullptr;
853 if( ptrLast == nullptr )
854 {
855 ptrLast= new G4ThreeVector( DBL_MAX, DBL_MAX, DBL_MAX );
856 G4AutoDelete::Register(ptrLast);
857 }
858 G4ThreeVector &lastStart= *ptrLast;
859
860 if( (TrialPoint - StartPosition).mag2() < tolerance*tolerance)
861 {
862 static G4ThreadLocal unsigned int numUnmoved = 0;
863 static G4ThreadLocal unsigned int numStill = 0; // Still at same point
864
865 G4cout << "Intersection F == start A in " << MethodName;
866 G4cout << "Start Point: " << StartPosition << G4endl;
867 G4cout << " Start-Trial: " << TrialPoint - StartPosition;
868 G4cout << " Start-last: " << StartPosition - lastStart;
869
870 if( (StartPosition - lastStart).mag() < tolerance )
871 {
872 // We are at position of last 'Start' position - ie unmoved
873 ++numUnmoved;
874 ++numStill;
875 G4cout << " { Unmoved: " << " still#= " << numStill
876 << " total # = " << numUnmoved << " } - ";
877 }
878 else
879 {
880 numStill = 0;
881 }
882 G4cout << " Occurred: " << ++occurredOnTop;
883 G4cout << " out of total calls= " << numCalls;
884 G4cout << G4endl;
885 lastStart = StartPosition;
886 }
887} // End of ReportImmediateHit()
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
#define fPreviousSftOrigin
#define fPreviousSafety
static constexpr double perThousand
Definition: G4SIunits.hh:326
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
double z() const
double x() const
double mag2() const
double y() const
double dot(const Hep3Vector &) const
double mag() const
G4ThreeVector NetTranslation() const
G4RotationMatrix NetRotation() const
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4ThreeVector TransformAxis(const G4ThreeVector &axis) const
G4VIntegrationDriver * GetIntegrationDriver()
const G4ThreeVector & GetMomentumDir() const
G4double GetCurveLength() const
G4ThreeVector GetPosition() const
G4ThreeVector GetMomentum() const
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()
G4VSolid * GetSolid() const
const G4AffineTransform & GetTopTransform() const
G4TouchableHistory * CreateTouchableHistory() const
G4VExternalNavigation * GetExternalNavigation() const
void CheckMode(G4bool mode)
virtual void LocateGlobalPointWithinVolume(const G4ThreeVector &position)
Definition: G4Navigator.cc:614
virtual G4ThreeVector GetGlobalExitNormal(const G4ThreeVector &point, G4bool *valid)
virtual G4TouchableHistoryHandle CreateTouchableHistoryHandle() const
const G4AffineTransform GetLocalToGlobalTransform() const
virtual G4VPhysicalVolume * LocateGlobalPointAndSetup(const G4ThreeVector &point, const G4ThreeVector *direction=nullptr, const G4bool pRelativeSearch=true, const G4bool ignoreDirection=true)
Definition: G4Navigator.cc:132
G4Navigator * Clone() const
void SetWorldVolume(G4VPhysicalVolume *pWorld)
G4bool IsCheckModeActive() const
const G4AffineTransform & GetGlobalToLocalTransform() const
const G4NavigationHistory * GetHistory() const
G4ThreeVector GetGlobalSurfaceNormal(const G4ThreeVector &CurrentE_Point, G4bool &validNormal)
G4Navigator * GetNavigatorFor()
G4ChordFinder * GetChordFinderFor()
G4VIntersectionLocator(G4Navigator *theNavigator)
G4TouchableHistory * fpTouchable
G4ThreeVector GetLastSurfaceNormal(const G4ThreeVector &intersectPoint, G4bool &validNormal) const
G4ThreeVector GetSurfaceNormal(const G4ThreeVector &CurrentInt_Point, G4bool &validNormal)
void LocateGlobalPointWithinVolumeCheckAndReport(const G4ThreeVector &pos, const G4String &CodeLocationInfo, G4int CheckMode)
void ReportTrialStep(G4int step_no, const G4ThreeVector &ChordAB_v, const G4ThreeVector &ChordEF_v, const G4ThreeVector &NewMomentumDir, const G4ThreeVector &NormalAtEntry, G4bool validNormal)
G4bool IntersectChord(const G4ThreeVector &StartPointA, const G4ThreeVector &EndPointB, G4double &NewSafety, G4double &PreviousSafety, G4ThreeVector &PreviousSftOrigin, G4double &LinearStepLength, G4ThreeVector &IntersectionPoint, G4bool *calledNavigator=nullptr)
G4bool CheckAndReEstimateEndpoint(const G4FieldTrack &CurrentStartA, const G4FieldTrack &EstimatedEndB, G4FieldTrack &RevisedEndPoint, G4int &errorCode)
void ReportProgress(std::ostream &oss, const G4FieldTrack &StartPointVel, const G4FieldTrack &EndPointVel, G4int substep_no, const G4FieldTrack &A_PtVel, const G4FieldTrack &B_PtVel, G4double safetyLast, G4int depth=-1)
G4ThreeVector GetLocalSurfaceNormal(const G4ThreeVector &CurrentE_Point, G4bool &validNormal)
G4double GetEpsilonStepFor()
void ReportImmediateHit(const char *MethodName, const G4ThreeVector &StartPosition, const G4ThreeVector &TrialPoint, G4double tolerance, unsigned long int numCalls)
G4FieldTrack ReEstimateEndpoint(const G4FieldTrack &CurrentStateA, const G4FieldTrack &EstimtdEndStateB, G4double linearDistSq, G4double curveDist)
G4bool LocateGlobalPointWithinVolumeAndCheck(const G4ThreeVector &pos)
void printStatus(const G4FieldTrack &startFT, const G4FieldTrack &currentFT, G4double requestStep, G4double safety, G4int stepNum)
G4bool AdjustmentOfFoundIntersection(const G4ThreeVector &A, const G4ThreeVector &CurrentE_Point, const G4ThreeVector &CurrentF_Point, const G4ThreeVector &MomentumDir, const G4bool IntersectAF, G4ThreeVector &IntersectionPoint, G4double &NewSafety, G4double &fPrevSafety, G4ThreeVector &fPrevSftOrigin)
void ReportReversedPoints(std::ostringstream &ossMsg, const G4FieldTrack &StartPointVel, const G4FieldTrack &EndPointVel, G4double NewSafety, G4double epsStep, const G4FieldTrack &CurrentA_PointVelocity, const G4FieldTrack &CurrentB_PointVelocity, const G4FieldTrack &SubStart_PointVelocity, const G4ThreeVector &CurrentE_Point, const G4FieldTrack &ApproxIntersecPointV, G4int sbstp_no, G4int sbstp_no_p, G4int depth)
G4LogicalVolume * GetLogicalVolume() const
virtual G4int GetCopyNo() const =0
virtual EInside Inside(const G4ThreeVector &p) const =0
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const =0
virtual G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const =0
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const =0
EInside
Definition: geomdefs.hh:67
@ kInside
Definition: geomdefs.hh:70
@ kSurface
Definition: geomdefs.hh:69
static constexpr double perThousand
static constexpr double millimeter
Definition: SystemOfUnits.h:63
void Register(T *inst)
Definition: G4AutoDelete.hh:65
G4bool transform(G4String &input, const G4String &type)
#define DBL_MAX
Definition: templates.hh:62
#define G4ThreadLocal
Definition: tls.hh:77