Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4MultiLevelLocator.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 G4MultiLevelLocator implementation
27//
28// 27.10.08 - Tatiana Nikitina.
29// 04.10.11 - John Apostolakis, revised convergence to use Surface Normal
30// ---------------------------------------------------------------------------
31
32#include <iomanip>
33
34#include "G4ios.hh"
38
40 : G4VIntersectionLocator(theNavigator)
41{
42 // In case of too slow progress in finding Intersection Point
43 // intermediates Points on the Track must be stored.
44 // Initialise the array of Pointers [max_depth+1] to do this
45
46 G4ThreeVector zeroV(0.0,0.0,0.0);
47 for ( auto idepth=0; idepth<max_depth+1; ++idepth )
48 {
49 ptrInterMedFT[ idepth ] = new G4FieldTrack( zeroV, zeroV, 0., 0., 0., 0.);
50 }
51
52 if (fCheckMode)
53 {
54 // Trial values Loose Medium Tight
55 // To happen: Infrequent Occasional Often
56 SetMaxSteps(150); // 300 85 25
57 SetWarnSteps(80); // 250 60 15
58 }
59}
60
62{
63 for ( auto idepth=0; idepth<max_depth+1; ++idepth )
64 {
65 delete ptrInterMedFT[idepth];
66 }
67#ifdef G4DEBUG_FIELD
69#endif
70}
71
72
73// --------------------------------------------------------------------------
74// G4bool G4PropagatorInField::LocateIntersectionPoint(
75// const G4FieldTrack& CurveStartPointVelocity, // A
76// const G4FieldTrack& CurveEndPointVelocity, // B
77// const G4ThreeVector& TrialPoint, // E
78// G4FieldTrack& IntersectedOrRecalculated // Output
79// G4bool& recalculated ) // Out
80// --------------------------------------------------------------------------
81//
82// Function that returns the intersection of the true path with the surface
83// of the current volume (either the external one or the inner one with one
84// of the daughters:
85//
86// A = Initial point
87// B = another point
88//
89// Both A and B are assumed to be on the true path:
90//
91// E is the first point of intersection of the chord AB with
92// a volume other than A (on the surface of A or of a daughter)
93//
94// Convention of Use :
95// i) If it returns "true", then IntersectionPointVelocity is set
96// to the approximate intersection point.
97// ii) If it returns "false", no intersection was found.
98// Potential reasons:
99// a) no segment found an intersection
100// b) too many steps were required - after that it abandoned the effort
101// and is returning how far it could go. (New - 29 Oct 2015)
102// (If so, it must set 'recalculated' to true.)
103// TODO/idea: add a new flag: 'unfinished' to identify these cases.
104//
105// IntersectedOrRecalculated means different things:
106// a) if it is the same curve lenght along, it is a revision of the
107// original enpdoint due to the need for re-integration.
108// b) if it is at a shorter curve length, it is 'end of what it could do'
109// i.e. as far as it could go, because it took too many steps!
110// Note: IntersectedOrRecalculated is valid only if 'recalculated' is
111// 'true'.
112// --------------------------------------------------------------------------
113// NOTE: implementation taken from G4PropagatorInField
114//
116 const G4FieldTrack& CurveStartPointVelocity, // A
117 const G4FieldTrack& CurveEndPointVelocity, // B
118 const G4ThreeVector& TrialPoint, // E
119 G4FieldTrack& IntersectedOrRecalculatedFT, // Output
120 G4bool& recalculatedEndPoint, // Out
121 G4double& previousSafety, // In/Out
122 G4ThreeVector& previousSftOrigin) // In/Out
123{
124 // Find Intersection Point ( A, B, E ) of true path AB - start at E.
125 const char* MethodName= "G4MultiLevelLocator::EstimateIntersectionPoint()";
126
127 G4bool found_approximate_intersection = false;
128 G4bool there_is_no_intersection = false;
129
130 G4FieldTrack CurrentA_PointVelocity = CurveStartPointVelocity;
131 G4FieldTrack CurrentB_PointVelocity = CurveEndPointVelocity;
132 G4ThreeVector CurrentE_Point = TrialPoint;
133 G4bool validNormalAtE = false;
134 G4ThreeVector NormalAtEntry;
135
136 G4FieldTrack ApproxIntersecPointV(CurveEndPointVelocity); // FT-Def-Construct
137 G4bool validIntersectP= true; // Is it current ?
138 G4double NewSafety = 0.0;
139 G4bool last_AF_intersection = false;
140
141 auto integrDriver = GetChordFinderFor()->GetIntegrationDriver();
142 G4bool driverReIntegrates = integrDriver->DoesReIntegrate();
143
144 G4bool first_section = true;
145 recalculatedEndPoint = false;
146 G4bool restoredFullEndpoint = false;
147
148 unsigned int substep_no = 0;
149
150 // Statistics for substeps
151 static G4ThreadLocal unsigned int max_no_seen= 0;
152
153#ifdef G4DEBUG_FIELD
154 unsigned int trigger_substepno_print = 0;
155 const G4double tolerance = 1.0e-8 * CLHEP::mm;
156 unsigned int biggest_depth = 0;
157 // using kInitialisingCL = G4LocatorChangeRecord::kInitialisingCL;
158#endif
159
160 // Log the location, iteration of changes in A,B
161 //----------------------------------------------
162 static G4ThreadLocal G4LocatorChangeLogger endChangeA("StartPointA"),
163 endChangeB("EndPointB"), recApproxPoint("ApproxPoint"),
164 pointH_logger("Trial points 'E': position, normal");
165 unsigned int eventCount = 0;
166
167 if (fCheckMode)
168 {
169 // Clear previous call's data
170 endChangeA.clear();
171 endChangeB.clear();
172 recApproxPoint.clear();
173 pointH_logger.clear();
174
175 // Record the initialisation
176 ++eventCount;
177 endChangeA.AddRecord( G4LocatorChangeRecord::kInitialisingCL, substep_no,
178 eventCount, CurrentA_PointVelocity );
179 endChangeB.AddRecord( G4LocatorChangeRecord::kInitialisingCL, substep_no,
180 eventCount, CurrentB_PointVelocity );
181 }
182
183 //--------------------------------------------------------------------------
184 // Algorithm for the case if progress in founding intersection is too slow.
185 // Process is defined too slow if after N=param_substeps advances on the
186 // path, it will be only 'fraction_done' of the total length.
187 // In this case the remaining length is divided in two half and
188 // the loop is restarted for each half.
189 // If progress is still too slow, the division in two halfs continue
190 // until 'max_depth'.
191 //--------------------------------------------------------------------------
192
193 const G4int param_substeps = 5; // Test value for the maximum number
194 // of substeps
195 const G4double fraction_done = 0.3;
196
197 G4bool Second_half = false; // First half or second half of divided step
198
199 // We need to know this for the 'final_section':
200 // real 'final_section' or first half 'final_section'
201 // In algorithm it is considered that the 'Second_half' is true
202 // and it becomes false only if we are in the first-half of level
203 // depthness or if we are in the first section
204
205 unsigned int depth = 0; // Depth counts subdivisions of initial step made
206 ++fNumCalls;
207
208 NormalAtEntry = GetSurfaceNormal(CurrentE_Point, validNormalAtE);
209
210 if (fCheckMode)
211 {
212 pointH_logger.AddRecord( G4LocatorChangeRecord::kInitialisingCL,
213 substep_no, eventCount,
214 G4FieldTrack( CurrentE_Point,0.,NormalAtEntry,0.,
215 0., 1.,G4ThreeVector(0.),0.,0.) );
216 #if (G4DEBUG_FIELD>1)
217 G4ThreeVector StartPosition = CurveStartPointVelocity.GetPosition();
218 if( (TrialPoint - StartPosition).mag2() < sqr(tolerance))
219 {
220 ReportImmediateHit( MethodName, StartPosition, TrialPoint,
221 tolerance, fNumCalls);
222 }
223 #endif
224 }
225
226 // Intermediates Points on the Track = Subdivided Points must be stored.
227 // Give the initial values to 'InterMedFt'
228 // Important is 'ptrInterMedFT[0]', it saves the 'EndCurvePoint'
229 //
230 *ptrInterMedFT[0] = CurveEndPointVelocity;
231 for ( auto idepth=1; idepth<max_depth+1; ++idepth )
232 {
233 *ptrInterMedFT[idepth] = CurveStartPointVelocity;
234 }
235
236 // Final_section boolean store
237 //
238 G4bool fin_section_depth[max_depth];
239 for ( auto idepth=0; idepth<max_depth; ++idepth )
240 {
241 fin_section_depth[idepth] = true;
242 }
243 // 'SubStartPoint' is needed to calculate the length of the divided step
244 //
245 G4FieldTrack SubStart_PointVelocity = CurveStartPointVelocity;
246
247 do // Loop checking, 07.10.2016, J.Apostolakis
248 {
249 unsigned int substep_no_p = 0;
250 G4bool sub_final_section = false; // the same as final_section,
251 // but for 'sub_section'
252 SubStart_PointVelocity = CurrentA_PointVelocity;
253
254 do // Loop checking, 07.10.2016, J.Apostolakis
255 { // REPEAT param
256 G4ThreeVector Point_A = CurrentA_PointVelocity.GetPosition();
257 G4ThreeVector Point_B = CurrentB_PointVelocity.GetPosition();
258
259#ifdef G4DEBUG_FIELD
260 const G4double lenA = CurrentA_PointVelocity.GetCurveLength() ;
261 const G4double lenB = CurrentB_PointVelocity.GetCurveLength() ;
262 G4double curv_lenAB = lenB - lenA;
263 G4double distAB = (Point_B - Point_A).mag();
264 if( curv_lenAB < distAB * ( 1. - 10.*fiEpsilonStep ) )
265 {
266 G4cerr << "ERROR> (Start) Point A coincides with or has gone past (end) point B"
267 << "MLL: iters = " << substep_no << G4endl;
268 G4long op=G4cerr.precision(6);
269 G4cerr << " Difference = " << distAB - curv_lenAB
270 << " exceeds limit of relative dist (10*epsilon)= " << 10*fiEpsilonStep
271 << " i.e. limit = " << 10 * fiEpsilonStep * distAB << G4endl;
272 G4cerr.precision(9);
273 G4cerr << " Len A, B = " << lenA << " " << lenB << G4endl
274 << " Position A: " << Point_A << G4endl
275 << " Position B: " << Point_B << G4endl;
276 G4cerr.precision(op);
277 // G4LocatorChangeRecord::ReportVector(G4cerr, "endPointB", endChangeB );
278 // G4cerr<<"EndPoints A(start) and B(end): combined changes " << G4endl;
279 if (fCheckMode) {
280 G4LocatorChangeLogger::ReportEndChanges(G4cerr, endChangeA, endChangeB);
281 }
282 }
283#endif
284 if( !validIntersectP ){
286 errmsg << "Assertion FAILURE - invalid (stale) Interection point. Substep: "
287 << substep_no << " call: " << fNumCalls << G4endl;
288 if (fCheckMode)
289 G4LocatorChangeRecord::ReportEndChanges(errmsg, endChangeA, endChangeB );
290 G4Exception("G4MultiLevelLocator::EstimateIntersectionPoint", "GeomNav0004",
291 JustWarning, errmsg);
292 }
293
294 // F = a point on true AB path close to point E
295 // (the closest if possible)
296 //
297 ApproxIntersecPointV = GetChordFinderFor()
298 ->ApproxCurvePointV( CurrentA_PointVelocity,
299 CurrentB_PointVelocity,
300 CurrentE_Point,
302 // The above method is the key & most intuitive part ...
303
304#ifdef G4DEBUG_FIELD
306 substep_no, eventCount, ApproxIntersecPointV ) );
307 G4double lenIntsc= ApproxIntersecPointV.GetCurveLength();
308 G4double checkVsEnd= lenB - lenIntsc;
309
310 if( lenIntsc > lenB )
311 {
312 std::ostringstream errmsg;
313 errmsg.precision(17);
314 G4double ratio = checkVsEnd / lenB;
315 G4double ratioTol = std::fabs(ratio) / tolerance;
316 errmsg << "Intermediate F point is past end B point" << G4endl
317 << " l( intersection ) = " << lenIntsc << G4endl
318 << " l( endpoint ) = " << lenB << G4endl;
319 errmsg.precision(8);
320 errmsg << " l_end - l_inters = " << checkVsEnd << G4endl
321 << " / l_end = " << ratio << G4endl
322 << " ratio / tolerance = " << ratioTol << G4endl;
323 if( ratioTol < 1.0 )
324 G4Exception(MethodName, "GeomNav0003", JustWarning, errmsg );
325 else
326 G4Exception(MethodName, "GeomNav0003", FatalException, errmsg );
327 }
328#endif
329
330 G4ThreeVector CurrentF_Point= ApproxIntersecPointV.GetPosition();
331
332 // First check whether EF is small - then F is a good approx. point
333 // Calculate the length and direction of the chord AF
334 //
335 G4ThreeVector ChordEF_Vector = CurrentF_Point - CurrentE_Point;
336
337 G4ThreeVector NewMomentumDir = ApproxIntersecPointV.GetMomentumDir();
338 G4double MomDir_dot_Norm = NewMomentumDir.dot( NormalAtEntry );
339
340#ifdef G4DEBUG_FIELD
341 if( fVerboseLevel > 3 )
342 {
343 G4ThreeVector ChordAB = Point_B - Point_A;
344 G4double ABchord_length = ChordAB.mag();
345 G4double MomDir_dot_ABchord;
346 MomDir_dot_ABchord = (1.0 / ABchord_length)
347 * NewMomentumDir.dot( ChordAB );
348 G4VIntersectionLocator::ReportTrialStep( substep_no, ChordAB,
349 ChordEF_Vector, NewMomentumDir, NormalAtEntry, validNormalAtE );
350 G4cout << " dot( MomentumDir, ABchord_unit ) = "
351 << MomDir_dot_ABchord << G4endl;
352 }
353#endif
354 G4bool adequate_angle =
355 ( MomDir_dot_Norm >= 0.0 ) // Can use ( > -epsilon) instead
356 || (! validNormalAtE ); // Invalid, cannot use this criterion
357 G4double EF_dist2 = ChordEF_Vector.mag2();
358 if ( ( EF_dist2 <= sqr(fiDeltaIntersection) && ( adequate_angle ) )
359 || ( EF_dist2 <= kCarTolerance*kCarTolerance ) )
360 {
361 found_approximate_intersection = true;
362
363 // Create the "point" return value
364 //
365 IntersectedOrRecalculatedFT = ApproxIntersecPointV;
366 IntersectedOrRecalculatedFT.SetPosition( CurrentE_Point );
367
369 {
370 // Try to Get Correction of IntersectionPoint using SurfaceNormal()
371 //
372 G4ThreeVector IP;
373 G4ThreeVector MomentumDir=ApproxIntersecPointV.GetMomentumDirection();
374 G4bool goodCorrection = AdjustmentOfFoundIntersection(Point_A,
375 CurrentE_Point, CurrentF_Point, MomentumDir,
376 last_AF_intersection, IP, NewSafety,
377 previousSafety, previousSftOrigin );
378 if ( goodCorrection )
379 {
380 IntersectedOrRecalculatedFT = ApproxIntersecPointV;
381 IntersectedOrRecalculatedFT.SetPosition(IP);
382 }
383 }
384 // Note: in order to return a point on the boundary,
385 // we must return E. But it is F on the curve.
386 // So we must "cheat": we are using the position at point E
387 // and the velocity at point F !!!
388 //
389 // This must limit the length we can allow for displacement!
390 }
391 else // E is NOT close enough to the curve (ie point F)
392 {
393 // Check whether any volumes are encountered by the chord AF
394 // ---------------------------------------------------------
395 // First relocate to restore any Voxel etc information
396 // in the Navigator before calling ComputeStep()
397 //
399
400 G4ThreeVector PointG; // Candidate intersection point
401 G4double stepLengthAF;
402 G4bool Intersects_FB = false;
403 G4bool Intersects_AF = IntersectChord( Point_A, CurrentF_Point,
404 NewSafety, previousSafety,
405 previousSftOrigin,
406 stepLengthAF,
407 PointG );
408 last_AF_intersection = Intersects_AF;
409 if( Intersects_AF )
410 {
411 // G is our new Candidate for the intersection point.
412 // It replaces "E" and we will see if it's good enough.
413 CurrentB_PointVelocity = ApproxIntersecPointV; // B <- F
414 CurrentE_Point = PointG; // E <- G
415
416 validIntersectP = true; // 'E' has been updated.
417
418 G4bool validNormalLast;
419 NormalAtEntry = GetSurfaceNormal( PointG, validNormalLast );
420 validNormalAtE = validNormalLast;
421
422 // As we move point B, must take care in case the current
423 // AF has no intersection to try current FB!!
424 fin_section_depth[depth] = false;
425
426 if (fCheckMode)
427 {
428 ++eventCount;
429 endChangeB.push_back(
431 substep_no, eventCount, CurrentB_PointVelocity) );
432 }
433#ifdef G4VERBOSE
434 if( fVerboseLevel > 3 )
435 {
436 G4cout << "G4PiF::LI> Investigating intermediate point"
437 << " at s=" << ApproxIntersecPointV.GetCurveLength()
438 << " on way to full s="
439 << CurveEndPointVelocity.GetCurveLength() << G4endl;
440 }
441#endif
442 }
443 else // not Intersects_AF
444 {
445 // In this case:
446 // There is NO intersection of AF with a volume boundary.
447 // We must continue the search in the segment FB!
448 //
450
451 G4double stepLengthFB;
452 G4ThreeVector PointH;
453
454 // Check whether any volumes are encountered by the chord FB
455 // ---------------------------------------------------------
456
457 Intersects_FB = IntersectChord( CurrentF_Point, Point_B,
458 NewSafety, previousSafety,
459 previousSftOrigin,
460 stepLengthFB,
461 PointH );
462 if( Intersects_FB )
463 {
464 // There is an intersection of FB with a volume boundary
465 // H <- First Intersection of Chord FB
466
467 // H is our new Candidate for the intersection point.
468 // It replaces "E" and we will repeat the test to see if
469 // it is a good enough approximate point for us.
470
471 // Note that F must be in volume volA (the same as A)
472 // (otherwise AF would meet a volume boundary!)
473 // A <- F
474 // E <- H
475 //
476 CurrentA_PointVelocity = ApproxIntersecPointV;
477 CurrentE_Point = PointH;
478
479 validIntersectP = true; // 'E' has been updated.
480
481 G4bool validNormalH;
482 NormalAtEntry = GetSurfaceNormal( PointH, validNormalH );
483 validNormalAtE = validNormalH;
484
485 if (fCheckMode)
486 {
487 ++eventCount;
488 endChangeA.push_back(
490 substep_no, eventCount, CurrentA_PointVelocity) );
491 G4FieldTrack intersectH_pn('0'); // Point and normal
492 // nothing else will be valid
493 intersectH_pn.SetPosition( PointH );
494 intersectH_pn.SetMomentum( NormalAtEntry );
495 pointH_logger.AddRecord(G4LocatorChangeRecord::kIntersectsFB,
496 substep_no, eventCount, intersectH_pn );
497 }
498 }
499 else // not Intersects_FB
500 {
501 validIntersectP = false; // Intersections are now stale
502 if( fin_section_depth[depth] )
503 {
504 // If B is the original endpoint, this means that whatever
505 // volume(s) intersected the original chord, none touch the
506 // smaller chords we have used.
507 // The value of 'IntersectedOrRecalculatedFT' returned is
508 // likely not valid
509
510 // Check on real final_section or SubEndSection
511 //
512 if( ((Second_half)&&(depth==0)) || (first_section) )
513 {
514 there_is_no_intersection = true; // real final_section
515 }
516 else
517 {
518 // end of subsection, not real final section
519 // exit from the and go to the depth-1 level
520 substep_no_p = param_substeps+2; // exit from the loop
521
522 // but 'Second_half' is still true because we need to find
523 // the 'CurrentE_point' for the next loop
524 Second_half = true;
525 sub_final_section = true;
526 }
527 }
528 else
529 {
530 CurrentA_PointVelocity = CurrentB_PointVelocity; // Got to B
531 CurrentB_PointVelocity = (depth==0) ? CurveEndPointVelocity
532 : *ptrInterMedFT[depth] ;
533 SubStart_PointVelocity = CurrentA_PointVelocity;
534 restoredFullEndpoint = true;
535
536 validIntersectP = false; // 'E' has NOT been updated.
537
538 if (fCheckMode)
539 {
540 ++eventCount;
541 endChangeA.push_back(
544 substep_no, eventCount, CurrentA_PointVelocity) );
545 endChangeB.push_back(
548 substep_no, eventCount, CurrentB_PointVelocity) );
549 }
550 }
551 } // Endif (Intersects_FB)
552 } // Endif (Intersects_AF)
553
554 G4int errorEndPt = 0; // Default: no error (if not calling CheckAnd...
555
556 G4bool recalculatedB= false;
557 if( driverReIntegrates )
558 {
559 G4FieldTrack RevisedB_FT = CurrentB_PointVelocity;
560 recalculatedB= CheckAndReEstimateEndpoint(CurrentA_PointVelocity,
561 CurrentB_PointVelocity,
562 RevisedB_FT,
563 errorEndPt );
564 if( recalculatedB )
565 {
566 CurrentB_PointVelocity = RevisedB_FT; // Use it !
567 // Do not invalidate intersection F -- it is still roughly OK.
568 //
569 // The best course would be to invalidate (reset validIntersectP)
570 // BUT if we invalidate it, we must re-estimate it somewhere! E.g.
571 // validIntersectP = false; // 'E' has NOT been updated.
572
573 if ( (fin_section_depth[depth]) // real final section
574 &&( first_section || ((Second_half)&&(depth==0)) ) )
575 {
576 recalculatedEndPoint = true;
577 IntersectedOrRecalculatedFT = RevisedB_FT;
578 // So that we can return it, if it is the endpoint!
579 }
580 // else
581 // Move forward the other points
582 // - or better flag it, so that they are re-computed when next used
583 // [ Implementation: a counter for # of recomputations
584 // => avoids extra work]
585 }
586 if (fCheckMode)
587 {
588 ++eventCount;
589 endChangeB.push_back(
591 substep_no, eventCount, RevisedB_FT ) );
592 }
593 }
594 else
595 {
596 if( CurrentB_PointVelocity.GetCurveLength() < CurrentA_PointVelocity.GetCurveLength() )
597 errorEndPt = 2;
598 }
599
600 if( errorEndPt > 1 ) // errorEndPt = 1 is milder, just: len(B)=len(A)
601 {
602 std::ostringstream errmsg;
603
605 CurveStartPointVelocity, CurveEndPointVelocity,
606 NewSafety, fiEpsilonStep,
607 CurrentA_PointVelocity, CurrentB_PointVelocity,
608 SubStart_PointVelocity, CurrentE_Point,
609 ApproxIntersecPointV, substep_no, substep_no_p, depth);
610
611 if (fCheckMode) {
612 G4LocatorChangeRecord::ReportEndChanges(errmsg, endChangeA, endChangeB );
613 }
614
615 errmsg << G4endl << " * Location: " << MethodName
616 << "- After EndIf(Intersects_AF)" << G4endl;
617 errmsg << " * Bool flags: Recalculated = " << recalculatedB
618 << " Intersects_AF = " << Intersects_AF
619 << " Intersects_FB = " << Intersects_FB << G4endl;
620 errmsg << " * Number of calls to MLL:EIP= " << fNumCalls << G4endl;
621 G4Exception(MethodName, "GeomNav0003", FatalException, errmsg);
622 }
623 if( restoredFullEndpoint )
624 {
625 fin_section_depth[depth] = restoredFullEndpoint;
626 restoredFullEndpoint = false;
627 }
628 } // EndIf ( E is close enough to the curve, ie point F. )
629 // tests ChordAF_Vector.mag() <= maximum_lateral_displacement
630
631#ifdef G4DEBUG_FIELD
632 if( trigger_substepno_print == 0)
633 {
634 trigger_substepno_print= fWarnSteps - 20;
635 }
636
637 if( substep_no >= trigger_substepno_print )
638 {
639 G4cout << "Difficulty in converging in " << MethodName
640 << " Substep no = " << substep_no << G4endl;
641 if( substep_no == trigger_substepno_print )
642 {
643 G4cout << " Start: ";
644 printStatus( CurveStartPointVelocity, CurveEndPointVelocity,
645 -1.0, NewSafety, 0 );
646 if( fCheckMode ) {
647 G4LocatorChangeRecord::ReportEndChanges(G4cout, endChangeA, endChangeB );
648 } else {
649 G4cout << " ** For more information enable 'check mode' in G4MultiLevelLocator "
650 << "-- (it saves and can output change records) " << G4endl;
651 }
652 }
653 G4cout << " Point A: ";
654 printStatus( CurrentA_PointVelocity, CurrentA_PointVelocity,
655 -1.0, NewSafety, substep_no-1 );
656 G4cout << " Point B: ";
657 printStatus( CurrentA_PointVelocity, CurrentB_PointVelocity,
658 -1.0, NewSafety, substep_no );
659 }
660#endif
661 ++substep_no;
662 ++substep_no_p;
663
664 } while ( ( ! found_approximate_intersection )
665 && ( ! there_is_no_intersection )
666 && validIntersectP // New condition: must refresh intersection !!
667 && ( substep_no_p <= param_substeps) ); // UNTIL found or
668 // failed param substep
669
670 if( (!found_approximate_intersection) && (!there_is_no_intersection) )
671 {
672 G4double did_len = std::abs( CurrentA_PointVelocity.GetCurveLength()
673 - SubStart_PointVelocity.GetCurveLength());
674 G4double all_len = std::abs( CurrentB_PointVelocity.GetCurveLength()
675 - SubStart_PointVelocity.GetCurveLength());
676
677 G4double distAB = -1;
678 //
679 // Is progress is too slow, and is it possible to go deeper?
680 // If so, then *halve the step*
681 // ==============
682 if( (did_len < fraction_done*all_len)
683 && (depth<max_depth) && (!sub_final_section) )
684 {
685#ifdef G4DEBUG_FIELD
686 static G4ThreadLocal unsigned int numSplits=0; // For debugging only
687 biggest_depth = std::max(depth, biggest_depth);
688 ++numSplits;
689#endif
690 Second_half = false;
691 ++depth;
692 first_section = false;
693
694 G4double Sub_len = (all_len-did_len)/(2.);
695 G4FieldTrack midPoint = CurrentA_PointVelocity;
696 G4bool fullAdvance=
697 integrDriver->AccurateAdvance(midPoint, Sub_len, fiEpsilonStep);
698
699 ++fNumAdvanceTrials;
700 if( fullAdvance ) { ++fNumAdvanceFull; }
701
702 G4double lenAchieved=
703 midPoint.GetCurveLength()-CurrentA_PointVelocity.GetCurveLength();
704
705 const G4double adequateFraction = (1.0-CLHEP::perThousand);
706 G4bool goodAdvance = (lenAchieved >= adequateFraction * Sub_len);
707 if ( goodAdvance ) { ++fNumAdvanceGood; }
708
709#ifdef G4DEBUG_FIELD
710 else // !goodAdvance
711 {
712 G4cout << "MLL> AdvanceChordLimited not full at depth=" << depth
713 << " total length achieved = " << lenAchieved << " of "
714 << Sub_len << " fraction= ";
715 if (Sub_len != 0.0 ) { G4cout << lenAchieved / Sub_len; }
716 else { G4cout << "DivByZero"; }
717 G4cout << " Good-enough fraction = " << adequateFraction;
718 G4cout << " Number of call to mll = " << fNumCalls
719 << " iteration " << substep_no
720 << " inner = " << substep_no_p << G4endl;
721 G4cout << " Epsilon of integration = " << fiEpsilonStep << G4endl;
722 G4cout << " State at start is = " << CurrentA_PointVelocity
723 << G4endl
724 << " at end (midpoint)= " << midPoint << G4endl;
725 G4cout << " Particle mass = " << midPoint.GetRestMass() << G4endl;
726
727 G4EquationOfMotion *equation = integrDriver->GetEquationOfMotion();
728 ReportFieldValue( CurrentA_PointVelocity, "start", equation );
729 ReportFieldValue( midPoint, "midPoint", equation );
730 G4cout << " Original Start = "
731 << CurveStartPointVelocity << G4endl;
732 G4cout << " Original End = "
733 << CurveEndPointVelocity << G4endl;
734 G4cout << " Original TrialPoint = "
735 << TrialPoint << G4endl;
736 G4cout << " (this is STRICT mode) "
737 << " num Splits= " << numSplits;
738 G4cout << G4endl;
739 }
740#endif
741
742 *ptrInterMedFT[depth] = midPoint;
743 CurrentB_PointVelocity = midPoint;
744
745 if (fCheckMode)
746 {
747 ++eventCount;
748 endChangeB.push_back(
750 substep_no, eventCount, midPoint) );
751 }
752
753 // Adjust 'SubStartPoint' to calculate the 'did_length' in next loop
754 //
755 SubStart_PointVelocity = CurrentA_PointVelocity;
756
757 // Find new trial intersection point needed at start of the loop
758 //
759 G4ThreeVector Point_A = CurrentA_PointVelocity.GetPosition();
760 G4ThreeVector Point_B = CurrentB_PointVelocity.GetPosition();
761
762 G4ThreeVector PointGe;
764 G4bool Intersects_AB = IntersectChord(Point_A, Point_B,
765 NewSafety, previousSafety,
766 previousSftOrigin, distAB,
767 PointGe);
768 if( Intersects_AB )
769 {
770 last_AF_intersection = Intersects_AB;
771 CurrentE_Point = PointGe;
772 fin_section_depth[depth] = true;
773
774 validIntersectP = true; // 'E' has been updated.
775
776 G4bool validNormalAB;
777 NormalAtEntry = GetSurfaceNormal( PointGe, validNormalAB );
778 validNormalAtE = validNormalAB;
779 }
780 else
781 {
782 // No intersection found for first part of curve
783 // (CurrentA,InterMedPoint[depth]). Go to the second part
784 //
785 Second_half = true;
786
787 validIntersectP= false; // No new 'E' chord intersection found
788 }
789 } // if did_len
790
791 G4bool unfinished = Second_half;
792 while ( unfinished && (depth>0) ) // Loop checking, 07.10.2016, JA
793 {
794 // Second part of curve (InterMed[depth],Intermed[depth-1]))
795 // On the depth-1 level normally we are on the 'second_half'
796
797 // Find new trial intersection point needed at start of the loop
798 //
799 SubStart_PointVelocity = *ptrInterMedFT[depth];
800 CurrentA_PointVelocity = *ptrInterMedFT[depth];
801 CurrentB_PointVelocity = *ptrInterMedFT[depth-1];
802
803 if (fCheckMode)
804 {
805 ++eventCount;
807 substep_no, eventCount, CurrentA_PointVelocity);
808 endChangeA.push_back( chngPop_a );
810 substep_no, eventCount, CurrentB_PointVelocity);
811 endChangeB.push_back( chngPop_b );
812 }
813
814 // Ensure that the new endpoints are not further apart in space
815 // than on the curve due to different errors in the integration
816 //
817 G4int errorEndPt = -1;
818 G4bool recalculatedB= false;
819 if( driverReIntegrates )
820 {
821 // Ensure that the new endpoints are not further apart in space
822 // than on the curve due to different errors in the integration
823 //
824 G4FieldTrack RevisedEndPointFT = CurrentB_PointVelocity;
825 recalculatedB =
826 CheckAndReEstimateEndpoint( CurrentA_PointVelocity,
827 CurrentB_PointVelocity,
828 RevisedEndPointFT,
829 errorEndPt );
830 if( recalculatedB )
831 {
832 CurrentB_PointVelocity = RevisedEndPointFT; // Use it !
833
834 if ( depth == 1 )
835 {
836 recalculatedEndPoint = true;
837 IntersectedOrRecalculatedFT = RevisedEndPointFT;
838 // So that we can return it, if it is the endpoint!
839 }
840 }
841 else
842 {
843 if( CurrentB_PointVelocity.GetCurveLength() < CurrentA_PointVelocity.GetCurveLength() )
844 errorEndPt = 2;
845 }
846
847 if (fCheckMode)
848 {
849 ++eventCount;
850 endChangeB.push_back(
852 substep_no, eventCount, RevisedEndPointFT));
853 }
854 }
855 if( errorEndPt > 1 ) // errorEndPt = 1 is milder, just: len(B)=len(A)
856 {
857 std::ostringstream errmsg;
858 ReportReversedPoints(errmsg,
859 CurveStartPointVelocity, CurveEndPointVelocity,
860 NewSafety, fiEpsilonStep,
861 CurrentA_PointVelocity, CurrentB_PointVelocity,
862 SubStart_PointVelocity, CurrentE_Point,
863 ApproxIntersecPointV, substep_no, substep_no_p, depth);
864 errmsg << " * Location: " << MethodName << "- Second-Half" << G4endl;
865 errmsg << " * Recalculated = " << recalculatedB << G4endl; // false
866 G4Exception(MethodName, "GeomNav0003", FatalException, errmsg);
867 }
868
869 G4ThreeVector Point_A = CurrentA_PointVelocity.GetPosition();
870 G4ThreeVector Point_B = CurrentB_PointVelocity.GetPosition();
871 G4ThreeVector PointGi;
873 G4bool Intersects_AB = IntersectChord(Point_A, Point_B, NewSafety,
874 previousSafety,
875 previousSftOrigin, distAB,
876 PointGi);
877 if( Intersects_AB )
878 {
879 last_AF_intersection = Intersects_AB;
880 CurrentE_Point = PointGi;
881
882 validIntersectP = true; // 'E' has been updated.
883 NormalAtEntry = GetSurfaceNormal( PointGi, validNormalAtE );
884 }
885 else
886 {
887 validIntersectP = false; // No new 'E' chord intersection found
888 if( depth == 1)
889 {
890 there_is_no_intersection = true;
891 }
892 }
893 depth--;
894 fin_section_depth[depth] = true;
895 unfinished = !validIntersectP;
896 }
897#ifdef G4DEBUG_FIELD
898 if( ! ( validIntersectP || there_is_no_intersection ) )
899 {
900 // What happened ??
901 G4cout << "MLL - WARNING Potential FAILURE: Conditions not met!"
902 << G4endl
903 << " Depth = " << depth << G4endl
904 << " Num Substeps= " << substep_no << G4endl;
905 G4cout << " Found intersection= " << found_approximate_intersection
906 << G4endl;
907 G4cout << " Progress report: -- " << G4endl;
909 CurveStartPointVelocity, CurveEndPointVelocity,
910 substep_no, CurrentA_PointVelocity,
911 CurrentB_PointVelocity,
912 NewSafety, depth );
913 G4cout << G4endl;
914 }
915#endif
916 } // if(!found_aproximate_intersection)
917
918 assert( validIntersectP || there_is_no_intersection
919 || found_approximate_intersection);
920
921 } while ( ( ! found_approximate_intersection )
922 && ( ! there_is_no_intersection )
923 && ( substep_no <= fMaxSteps) ); // UNTIL found or failed
924
925 if( substep_no > max_no_seen )
926 {
927 max_no_seen = substep_no;
928#ifdef G4DEBUG_FIELD
929 if( max_no_seen > fWarnSteps )
930 {
931 trigger_substepno_print = max_no_seen-20; // Want to see last 20 steps
932 }
933#endif
934 }
935
936 if( !there_is_no_intersection && !found_approximate_intersection )
937 {
938 if( substep_no >= fMaxSteps)
939 {
940 // Since we cannot go further (yet), we return as far as we have gone
941
942 recalculatedEndPoint = true;
943 IntersectedOrRecalculatedFT = CurrentA_PointVelocity;
944 found_approximate_intersection = false;
945
946 std::ostringstream message;
947 message << G4endl;
948 message << "Convergence is requiring too many substeps: "
949 << substep_no << " ( limit = "<< fMaxSteps << ")"
950 << G4endl
951 << " Abandoning effort to intersect. " << G4endl << G4endl;
952 message << " Number of calls to MLL: " << fNumCalls;
953 message << " iteration = " << substep_no <<G4endl << G4endl;
954
955 message.precision( 12 );
956 G4double done_len = CurrentA_PointVelocity.GetCurveLength();
957 G4double full_len = CurveEndPointVelocity.GetCurveLength();
958 message << " Undertaken only length: " << done_len
959 << " out of " << full_len << " required." << G4endl
960 << " Remaining length = " << full_len - done_len;
961
962 message << " Start and end-point of requested Step:" << G4endl;
963 printStatus( CurveStartPointVelocity, CurveEndPointVelocity,
964 -1.0, NewSafety, 0, message, -1 );
965 message << " Start and end-point of current Sub-Step:" << G4endl;
966 printStatus( CurrentA_PointVelocity, CurrentA_PointVelocity,
967 -1.0, NewSafety, substep_no-1, message, -1 );
968 printStatus( CurrentA_PointVelocity, CurrentB_PointVelocity,
969 -1.0, NewSafety, substep_no, message, -1 );
970
971 G4Exception(MethodName, "GeomNav0003", JustWarning, message);
972 }
973 else if( substep_no >= fWarnSteps)
974 {
975 std::ostringstream message;
976 message << "Many substeps while trying to locate intersection."
977 << G4endl
978 << " Undertaken length: "
979 << CurrentB_PointVelocity.GetCurveLength()
980 << " - Needed: " << substep_no << " substeps." << G4endl
981 << " Warning number = " << fWarnSteps
982 << " and maximum substeps = " << fMaxSteps;
983 G4Exception(MethodName, "GeomNav1002", JustWarning, message);
984 }
985 }
986
987 return (!there_is_no_intersection) && found_approximate_intersection;
988 // Success or failure
989}
990
992{
993 G4cout << " Number of calls = " << fNumCalls << G4endl;
994 G4cout << " Number of split level ('advances'): "
995 << fNumAdvanceTrials << G4endl;
996 G4cout << " Number of full advances: "
997 << fNumAdvanceGood << G4endl;
998 G4cout << " Number of good advances: "
999 << fNumAdvanceFull << G4endl;
1000}
1001
1002void G4MultiLevelLocator::ReportFieldValue( const G4FieldTrack& locationPV,
1003 const char* nameLoc,
1004 const G4EquationOfMotion* equation )
1005{
1006 enum { maxNumFieldComp = 24 };
1007
1008 G4ThreeVector position = locationPV.GetPosition();
1009 G4double startPoint[4] = { position.x(), position.y(), position.z(),
1010 locationPV.GetLabTimeOfFlight() };
1011 G4double FieldVec[maxNumFieldComp]; // 24 ;
1012 for (auto i=0; i<maxNumFieldComp; ++i )
1013 {
1014 FieldVec[i] = 0.0;
1015 }
1016 equation->GetFieldValue( startPoint, FieldVec);
1017 G4cout << " B-field value (" << nameLoc << ")= "
1018 << FieldVec[0] << " " << FieldVec[1] << " " << FieldVec[2];
1019 G4double Emag2= G4ThreeVector( FieldVec[3],
1020 FieldVec[4],
1021 FieldVec[5] ).mag2();
1022 if( Emag2 > 0.0 )
1023 {
1024 G4cout << " Electric = " << FieldVec[3] << " "
1025 << FieldVec[4] << " "
1026 << FieldVec[5]<< G4endl;
1027 }
1028 return;
1029}
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
const size_t max_depth
Definition: G4Octree.hh:54
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
long G4long
Definition: G4Types.hh:87
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 mag2() const
double dot(const Hep3Vector &) const
double mag() const
G4FieldTrack ApproxCurvePointV(const G4FieldTrack &curveAPointVelocity, const G4FieldTrack &curveBPointVelocity, const G4ThreeVector &currentEPoint, G4double epsStep)
G4VIntegrationDriver * GetIntegrationDriver()
void GetFieldValue(const G4double Point[4], G4double Field[]) const
const G4ThreeVector & GetMomentumDir() const
void SetMomentum(G4ThreeVector nMomDir)
G4double GetCurveLength() const
G4ThreeVector GetMomentumDirection() const
G4ThreeVector GetPosition() const
void SetPosition(G4ThreeVector nPos)
G4double GetRestMass() const
G4double GetLabTimeOfFlight() const
static std::ostream & ReportEndChanges(std::ostream &os, const G4LocatorChangeLogger &startA, const G4LocatorChangeLogger &endB)
static std::ostream & ReportEndChanges(std::ostream &os, const std::vector< G4LocatorChangeRecord > &startA, const std::vector< G4LocatorChangeRecord > &endB)
G4bool EstimateIntersectionPoint(const G4FieldTrack &curveStartPointTangent, const G4FieldTrack &curveEndPointTangent, const G4ThreeVector &trialPoint, G4FieldTrack &intersectPointTangent, G4bool &recalculatedEndPoint, G4double &fPreviousSafety, G4ThreeVector &fPreviousSftOrigin)
G4MultiLevelLocator(G4Navigator *theNavigator)
void SetMaxSteps(unsigned int valMax)
void SetWarnSteps(unsigned int valWarn)
virtual void LocateGlobalPointWithinVolume(const G4ThreeVector &position)
Definition: G4Navigator.cc:600
virtual G4bool DoesReIntegrate() const =0
G4Navigator * GetNavigatorFor()
G4ChordFinder * GetChordFinderFor()
G4ThreeVector GetSurfaceNormal(const G4ThreeVector &CurrentInt_Point, G4bool &validNormal)
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)
G4double GetEpsilonStepFor()
void ReportImmediateHit(const char *MethodName, const G4ThreeVector &StartPosition, const G4ThreeVector &TrialPoint, G4double tolerance, unsigned long int numCalls)
G4bool GetAdjustementOfFoundIntersection()
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)
T sqr(const T &x)
Definition: templates.hh:128
#define G4ThreadLocal
Definition: tls.hh:77