Geant4-11
Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes | Static Private Attributes
G4BrentLocator Class Reference

#include <G4BrentLocator.hh>

Inheritance diagram for G4BrentLocator:
G4VIntersectionLocator

Public Member Functions

void AddAdjustementOfFoundIntersection (G4bool UseCorrection)
 
void AdjustIntersections (G4bool UseCorrection)
 
G4bool AreIntersectionsAdjusted ()
 
G4bool EstimateIntersectionPoint (const G4FieldTrack &curveStartPointTangent, const G4FieldTrack &curveEndPointTangent, const G4ThreeVector &trialPoint, G4FieldTrack &intersectPointTangent, G4bool &recalculatedEndPoint, G4double &fPreviousSafety, G4ThreeVector &fPreviousSftOrigin)
 
 G4BrentLocator (G4Navigator *theNavigator)
 
G4bool GetAdjustementOfFoundIntersection ()
 
G4bool GetCheckMode ()
 
G4ChordFinderGetChordFinderFor ()
 
G4double GetDeltaIntersectionFor ()
 
G4double GetEpsilonStepFor ()
 
G4NavigatorGetNavigatorFor ()
 
G4int GetVerboseFor ()
 
G4bool IntersectChord (const G4ThreeVector &StartPointA, const G4ThreeVector &EndPointB, G4double &NewSafety, G4double &PreviousSafety, G4ThreeVector &PreviousSftOrigin, G4double &LinearStepLength, G4ThreeVector &IntersectionPoint, G4bool *calledNavigator=nullptr)
 
void printStatus (const G4FieldTrack &startFT, const G4FieldTrack &currentFT, G4double requestStep, G4double safety, G4int stepNum)
 
void SetCheckMode (G4bool value)
 
void SetChordFinderFor (G4ChordFinder *fCFinder)
 
void SetDeltaIntersectionFor (G4double deltaIntersection)
 
void SetEpsilonStepFor (G4double EpsilonStep)
 
void SetNavigatorFor (G4Navigator *fNavigator)
 
void SetSafetyParametersFor (G4bool UseSafety)
 
void SetVerboseFor (G4int fVerbose)
 
 ~G4BrentLocator ()
 

Static Public Member Functions

static void printStatus (const G4FieldTrack &startFT, const G4FieldTrack &currentFT, G4double requestStep, G4double safety, G4int stepNum, std::ostream &oss, G4int verboseLevel)
 

Protected Member Functions

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)
 
G4bool CheckAndReEstimateEndpoint (const G4FieldTrack &CurrentStartA, const G4FieldTrack &EstimatedEndB, G4FieldTrack &RevisedEndPoint, G4int &errorCode)
 
G4ThreeVector GetGlobalSurfaceNormal (const G4ThreeVector &CurrentE_Point, G4bool &validNormal)
 
G4ThreeVector GetSurfaceNormal (const G4ThreeVector &CurrentInt_Point, G4bool &validNormal)
 
G4bool LocateGlobalPointWithinVolumeAndCheck (const G4ThreeVector &pos)
 
void LocateGlobalPointWithinVolumeCheckAndReport (const G4ThreeVector &pos, const G4String &CodeLocationInfo, G4int CheckMode)
 
G4FieldTrack ReEstimateEndpoint (const G4FieldTrack &CurrentStateA, const G4FieldTrack &EstimtdEndStateB, G4double linearDistSq, G4double curveDist)
 
void ReportImmediateHit (const char *MethodName, const G4ThreeVector &StartPosition, const G4ThreeVector &TrialPoint, G4double tolerance, unsigned long int numCalls)
 
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)
 
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)
 
void ReportTrialStep (G4int step_no, const G4ThreeVector &ChordAB_v, const G4ThreeVector &ChordEF_v, const G4ThreeVector &NewMomentumDir, const G4ThreeVector &NormalAtEntry, G4bool validNormal)
 

Protected Attributes

G4bool fCheckMode = false
 
G4NavigatorfHelpingNavigator
 
G4ChordFinderfiChordFinder = nullptr
 
G4double fiDeltaIntersection = -1.0
 
G4double fiEpsilonStep = -1.0
 
G4NavigatorfiNavigator
 
G4bool fiUseSafety = false
 
G4TouchableHistoryfpTouchable = nullptr
 
G4bool fUseNormalCorrection = false
 
G4int fVerboseLevel = 0
 
G4double kCarTolerance
 

Private Member Functions

G4ThreeVector GetLastSurfaceNormal (const G4ThreeVector &intersectPoint, G4bool &validNormal) const
 
G4ThreeVector GetLocalSurfaceNormal (const G4ThreeVector &CurrentE_Point, G4bool &validNormal)
 

Private Attributes

G4FieldTrackptrInterMedFT [max_depth+1]
 

Static Private Attributes

static const G4int max_depth = 4
 

Detailed Description

Definition at line 47 of file G4BrentLocator.hh.

Constructor & Destructor Documentation

◆ G4BrentLocator()

G4BrentLocator::G4BrentLocator ( G4Navigator theNavigator)

Definition at line 37 of file G4BrentLocator.cc.

38 : G4VIntersectionLocator(theNavigator)
39{
40 // In case of too slow progress in finding Intersection Point
41 // intermediates Points on the Track must be stored.
42 // Initialise the array of Pointers [max_depth+1] to do this
43
44 G4ThreeVector zeroV(0.0,0.0,0.0);
45 for (auto idepth=0; idepth<max_depth+1; ++idepth )
46 {
47 ptrInterMedFT[ idepth ] = new G4FieldTrack( zeroV, zeroV, 0., 0., 0., 0.);
48 }
49}
static const G4int max_depth
G4FieldTrack * ptrInterMedFT[max_depth+1]
G4VIntersectionLocator(G4Navigator *theNavigator)

References max_depth, and ptrInterMedFT.

◆ ~G4BrentLocator()

G4BrentLocator::~G4BrentLocator ( )

Definition at line 51 of file G4BrentLocator.cc.

52{
53 for ( auto idepth=0; idepth<max_depth+1; ++idepth )
54 {
55 delete ptrInterMedFT[idepth];
56 }
57}

References max_depth, and ptrInterMedFT.

Member Function Documentation

◆ AddAdjustementOfFoundIntersection()

void G4VIntersectionLocator::AddAdjustementOfFoundIntersection ( G4bool  UseCorrection)
inlineinherited

◆ AdjustIntersections()

void G4VIntersectionLocator::AdjustIntersections ( G4bool  UseCorrection)
inlineinherited

◆ AdjustmentOfFoundIntersection()

G4bool G4VIntersectionLocator::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 
)
protectedinherited

Definition at line 439 of file G4VIntersectionLocator.cc.

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}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
#define fPreviousSftOrigin
#define fPreviousSafety
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
double dot(const Hep3Vector &) const
virtual void LocateGlobalPointWithinVolume(const G4ThreeVector &position)
Definition: G4Navigator.cc:614
G4ThreeVector GetGlobalSurfaceNormal(const G4ThreeVector &CurrentE_Point, G4bool &validNormal)
G4Navigator * GetNavigatorFor()
G4bool IntersectChord(const G4ThreeVector &StartPointA, const G4ThreeVector &EndPointB, G4double &NewSafety, G4double &PreviousSafety, G4ThreeVector &PreviousSftOrigin, G4double &LinearStepLength, G4ThreeVector &IntersectionPoint, G4bool *calledNavigator=nullptr)

References CLHEP::Hep3Vector::dot(), fPreviousSafety, fPreviousSftOrigin, G4VIntersectionLocator::fVerboseLevel, G4Exception(), G4VIntersectionLocator::GetGlobalSurfaceNormal(), G4VIntersectionLocator::GetNavigatorFor(), G4VIntersectionLocator::IntersectChord(), JustWarning, G4VIntersectionLocator::kCarTolerance, G4InuclParticleNames::lambda, and G4Navigator::LocateGlobalPointWithinVolume().

Referenced by EstimateIntersectionPoint(), G4MultiLevelLocator::EstimateIntersectionPoint(), and G4SimpleLocator::EstimateIntersectionPoint().

◆ AreIntersectionsAdjusted()

G4bool G4VIntersectionLocator::AreIntersectionsAdjusted ( )
inlineinherited

◆ CheckAndReEstimateEndpoint()

G4bool G4VIntersectionLocator::CheckAndReEstimateEndpoint ( const G4FieldTrack CurrentStartA,
const G4FieldTrack EstimatedEndB,
G4FieldTrack RevisedEndPoint,
G4int errorCode 
)
protectedinherited

Definition at line 329 of file G4VIntersectionLocator.cc.

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}
G4double GetCurveLength() const
G4ThreeVector GetPosition() const
G4FieldTrack ReEstimateEndpoint(const G4FieldTrack &CurrentStateA, const G4FieldTrack &EstimtdEndStateB, G4double linearDistSq, G4double curveDist)

References G4VIntersectionLocator::fiEpsilonStep, G4Exception(), G4FieldTrack::GetCurveLength(), G4FieldTrack::GetPosition(), JustWarning, and G4VIntersectionLocator::ReEstimateEndpoint().

Referenced by G4MultiLevelLocator::EstimateIntersectionPoint().

◆ EstimateIntersectionPoint()

G4bool G4BrentLocator::EstimateIntersectionPoint ( const G4FieldTrack curveStartPointTangent,
const G4FieldTrack curveEndPointTangent,
const G4ThreeVector trialPoint,
G4FieldTrack intersectPointTangent,
G4bool recalculatedEndPoint,
G4double fPreviousSafety,
G4ThreeVector fPreviousSftOrigin 
)
virtual

Implements G4VIntersectionLocator.

Definition at line 92 of file G4BrentLocator.cc.

101{
102 // Find Intersection Point ( A, B, E ) of true path AB - start at E.
103
104 G4bool found_approximate_intersection = false;
105 G4bool there_is_no_intersection = false;
106
107 G4FieldTrack CurrentA_PointVelocity = CurveStartPointVelocity;
108 G4FieldTrack CurrentB_PointVelocity = CurveEndPointVelocity;
109 G4ThreeVector CurrentE_Point = TrialPoint;
110 G4bool validNormalAtE = false;
111 G4ThreeVector NormalAtEntry;
112
113 G4FieldTrack ApproxIntersecPointV(CurveEndPointVelocity); // FT-Def-Construct
114 G4double NewSafety = 0.0;
115 G4bool last_AF_intersection = false;
116
117 // G4bool final_section= true; // Shows whether current section is last
118 // (i.e. B=full end)
119 G4bool first_section = true;
120 recalculatedEndPoint = false;
121
122 G4bool restoredFullEndpoint = false;
123
124 G4int oldprc; // cout, cerr precision
125 G4int substep_no = 0;
126
127 // Limits for substep number
128 //
129 const G4int max_substeps= 10000; // Test 120 (old value 100 )
130 const G4int warn_substeps= 1000; // 100
131
132 // Statistics for substeps
133 //
134 static G4ThreadLocal G4int max_no_seen= -1;
135
136 // Counter for restarting Bintermed
137 //
138 G4int restartB = 0;
139
140 //--------------------------------------------------------------------------
141 // Algorithm for the case if progress in founding intersection is too slow.
142 // Process is defined too slow if after N=param_substeps advances on the
143 // path, it will be only 'fraction_done' of the total length.
144 // In this case the remaining length is divided in two half and
145 // the loop is restarted for each half.
146 // If progress is still too slow, the division in two halfs continue
147 // until 'max_depth'.
148 //--------------------------------------------------------------------------
149
150 const G4int param_substeps = 50; // Test value for the maximum number
151 // of substeps
152 const G4double fraction_done = 0.3;
153
154 G4bool Second_half = false; // First half or second half of divided step
155
156 NormalAtEntry = GetSurfaceNormal(CurrentE_Point, validNormalAtE);
157
158 // We need to know this for the 'final_section':
159 // real 'final_section' or first half 'final_section'
160 // In algorithm it is considered that the 'Second_half' is true
161 // and it becomes false only if we are in the first-half of level
162 // depthness or if we are in the first section
163
164 G4int depth = 0; // Depth counts how many subdivisions of initial step made
165
166#ifdef G4DEBUG_FIELD
167 const G4double tolerance = 1.0e-8;
168 G4ThreeVector StartPosition = CurveStartPointVelocity.GetPosition();
169 if( (TrialPoint - StartPosition).mag() < tolerance * CLHEP::mm )
170 {
171 G4Exception("G4BrentLocator::EstimateIntersectionPoint()",
172 "GeomNav1002", JustWarning,
173 "Intersection point F is exactly at start point A." );
174 }
175#endif
176
177 // Intermediates Points on the Track = Subdivided Points must be stored.
178 // Give the initial values to 'InterMedFt'
179 // Important is 'ptrInterMedFT[0]', it saves the 'EndCurvePoint'
180 //
181 *ptrInterMedFT[0] = CurveEndPointVelocity;
182 for (auto idepth=1; idepth<max_depth+1; ++idepth )
183 {
184 *ptrInterMedFT[idepth] = CurveStartPointVelocity;
185 }
186
187 //Final_section boolean store
188 G4bool fin_section_depth[max_depth];
189 for (auto idepth=0; idepth<max_depth; ++idepth )
190 {
191 fin_section_depth[idepth] = true;
192 }
193
194 // 'SubStartPoint' is needed to calculate the length of the divided step
195 //
196 G4FieldTrack SubStart_PointVelocity = CurveStartPointVelocity;
197
198 do // Loop checking, 07.10.2016, J.Apostolakis
199 {
200 G4int substep_no_p = 0;
201 G4bool sub_final_section = false; // the same as final_section,
202 // but for 'sub_section'
203 SubStart_PointVelocity = CurrentA_PointVelocity;
204
205 do // Loop checking, 07.10.2016, J.Apostolakis
206 { // REPEAT param
207 G4ThreeVector Point_A = CurrentA_PointVelocity.GetPosition();
208 G4ThreeVector Point_B = CurrentB_PointVelocity.GetPosition();
209
210 // F = a point on true AB path close to point E
211 // (the closest if possible)
212 //
213 if(substep_no_p==0)
214 {
215 ApproxIntersecPointV = GetChordFinderFor()
216 ->ApproxCurvePointV( CurrentA_PointVelocity,
217 CurrentB_PointVelocity,
218 CurrentE_Point,
220 // The above method is the key & most intuitive part ...
221 }
222#ifdef G4DEBUG_FIELD
223 if( ApproxIntersecPointV.GetCurveLength() >
224 CurrentB_PointVelocity.GetCurveLength() * (1.0 + tolerance) )
225 {
226 G4Exception("G4BrentLocator::EstimateIntersectionPoint()",
227 "GeomNav0003", FatalException,
228 "Intermediate F point is past end B point" );
229 }
230#endif
231
232 G4ThreeVector CurrentF_Point = ApproxIntersecPointV.GetPosition();
233
234 // First check whether EF is small - then F is a good approx. point
235 // Calculate the length and direction of the chord AF
236 //
237 G4ThreeVector ChordEF_Vector = CurrentF_Point - CurrentE_Point;
238 G4ThreeVector NewMomentumDir = ApproxIntersecPointV.GetMomentumDir();
239 G4double MomDir_dot_Norm = NewMomentumDir.dot( NormalAtEntry ) ;
240
241#ifdef G4DEBUG_FIELD
242 G4ThreeVector ChordAB = Point_B - Point_A;
243
244 G4VIntersectionLocator::ReportTrialStep( substep_no, ChordAB,
245 ChordEF_Vector, NewMomentumDir, NormalAtEntry, validNormalAtE );
246#endif
247
248 G4bool adequate_angle;
249 adequate_angle = ( MomDir_dot_Norm >= 0.0 ) // Can use -epsilon instead.
250 || (! validNormalAtE ); // Makes criterion invalid
251 G4double EF_dist2 = ChordEF_Vector.mag2();
252 if ( ( EF_dist2 <= sqr(fiDeltaIntersection) && ( adequate_angle ) )
253 || ( EF_dist2 <= kCarTolerance*kCarTolerance ) )
254 {
255 found_approximate_intersection = true;
256
257 // Create the "point" return value
258 //
259 IntersectedOrRecalculatedFT = ApproxIntersecPointV;
260 IntersectedOrRecalculatedFT.SetPosition( CurrentE_Point );
261
263 {
264 // Try to Get Correction of IntersectionPoint using SurfaceNormal()
265 //
266 G4ThreeVector IP;
267 G4ThreeVector MomentumDir=ApproxIntersecPointV.GetMomentumDirection();
268 G4bool goodCorrection = AdjustmentOfFoundIntersection( Point_A,
269 CurrentE_Point, CurrentF_Point, MomentumDir,
270 last_AF_intersection, IP, NewSafety,
272 if ( goodCorrection )
273 {
274 IntersectedOrRecalculatedFT = ApproxIntersecPointV;
275 IntersectedOrRecalculatedFT.SetPosition(IP);
276 }
277 }
278
279 // Note: in order to return a point on the boundary,
280 // we must return E. But it is F on the curve.
281 // So we must "cheat": we are using the position at point E
282 // and the velocity at point F !!!
283 //
284 // This must limit the length we can allow for displacement!
285 }
286 else // E is NOT close enough to the curve (ie point F)
287 {
288 // Check whether any volumes are encountered by the chord AF
289 // ---------------------------------------------------------
290 // First relocate to restore any Voxel etc information
291 // in the Navigator before calling ComputeStep()
292 //
294
295 G4ThreeVector PointG; // Candidate intersection point
296 G4double stepLengthAF;
297 G4bool usedNavigatorAF = false;
298 G4bool Intersects_AF = IntersectChord( Point_A, CurrentF_Point,
299 NewSafety,fPreviousSafety,
301 stepLengthAF,
302 PointG,
303 &usedNavigatorAF);
304 last_AF_intersection = Intersects_AF;
305 if( Intersects_AF )
306 {
307 // G is our new Candidate for the intersection point.
308 // It replaces "E" and we will repeat the test to see if
309 // it is a good enough approximate point for us.
310 // B <- F
311 // E <- G
312 //
313 G4FieldTrack EndPoint = ApproxIntersecPointV;
314 ApproxIntersecPointV = GetChordFinderFor()->ApproxCurvePointS(
315 CurrentA_PointVelocity, CurrentB_PointVelocity,
316 EndPoint,CurrentE_Point, CurrentF_Point,PointG,
317 true, GetEpsilonStepFor() );
318 CurrentB_PointVelocity = EndPoint;
319 CurrentE_Point = PointG;
320
321 // Need to recalculate the Exit Normal at the new PointG
322 // Know that a call was made to Navigator::ComputeStep in
323 // IntersectChord above.
324 //
325 G4bool validNormalLast;
326 NormalAtEntry = GetSurfaceNormal( PointG, validNormalLast );
327 validNormalAtE = validNormalLast;
328
329 // By moving point B, must take care if current
330 // AF has no intersection to try current FB!!
331 //
332 fin_section_depth[depth] = false;
333#ifdef G4VERBOSE
334 if( fVerboseLevel > 3 )
335 {
336 G4cout << "G4PiF::LI> Investigating intermediate point"
337 << " at s=" << ApproxIntersecPointV.GetCurveLength()
338 << " on way to full s="
339 << CurveEndPointVelocity.GetCurveLength() << G4endl;
340 }
341#endif
342 }
343 else // not Intersects_AF
344 {
345 // In this case:
346 // There is NO intersection of AF with a volume boundary.
347 // We must continue the search in the segment FB!
348 //
350
351 G4double stepLengthFB;
352 G4ThreeVector PointH;
353 G4bool usedNavigatorFB = false;
354
355 // Check whether any volumes are encountered by the chord FB
356 // ---------------------------------------------------------
357
358 G4bool Intersects_FB = IntersectChord( CurrentF_Point, Point_B,
359 NewSafety,fPreviousSafety,
361 stepLengthFB,
362 PointH,
363 &usedNavigatorFB);
364 if( Intersects_FB )
365 {
366 // There is an intersection of FB with a volume boundary
367 // H <- First Intersection of Chord FB
368
369 // H is our new Candidate for the intersection point.
370 // It replaces "E" and we will repeat the test to see if
371 // it is a good enough approximate point for us.
372
373 // Note that F must be in volume volA (the same as A)
374 // (otherwise AF would meet a volume boundary!)
375 // A <- F
376 // E <- H
377 //
378 G4FieldTrack InterMed = ApproxIntersecPointV;
379 ApproxIntersecPointV = GetChordFinderFor()->ApproxCurvePointS(
380 CurrentA_PointVelocity,CurrentB_PointVelocity,
381 InterMed,CurrentE_Point,CurrentF_Point,PointH,
382 false,GetEpsilonStepFor());
383 CurrentA_PointVelocity = InterMed;
384 CurrentE_Point = PointH;
385
386 // Need to recalculate the Exit Normal at the new PointG
387 //
388 G4bool validNormalLast;
389 NormalAtEntry = GetSurfaceNormal( PointH, validNormalLast );
390 validNormalAtE = validNormalLast;
391 }
392 else // not Intersects_FB
393 {
394 // There is NO intersection of FB with a volume boundary
395
396 if( fin_section_depth[depth] )
397 {
398 // If B is the original endpoint, this means that whatever
399 // volume(s) intersected the original chord, none touch the
400 // smaller chords we have used.
401 // The value of 'IntersectedOrRecalculatedFT' returned is
402 // likely not valid
403
404 // Check on real final_section or SubEndSection
405 //
406 if( ((Second_half)&&(depth==0)) || (first_section) )
407 {
408 there_is_no_intersection = true; // real final_section
409 }
410 else
411 {
412 // end of subsection, not real final section
413 // exit from the and go to the depth-1 level
414
415 substep_no_p = param_substeps+2; // exit from the loop
416
417 // but 'Second_half' is still true because we need to find
418 // the 'CurrentE_point' for the next loop
419 //
420 Second_half = true;
421 sub_final_section = true;
422 }
423 }
424 else
425 {
426 if( depth==0 )
427 {
428 // We must restore the original endpoint
429 //
430 CurrentA_PointVelocity = CurrentB_PointVelocity; // Got to B
431 CurrentB_PointVelocity = CurveEndPointVelocity;
432 SubStart_PointVelocity = CurrentA_PointVelocity;
433 ApproxIntersecPointV = GetChordFinderFor()
434 ->ApproxCurvePointV( CurrentA_PointVelocity,
435 CurrentB_PointVelocity,
436 CurrentE_Point,
438
439 restoredFullEndpoint = true;
440 ++restartB; // counter
441 }
442 else
443 {
444 // We must restore the depth endpoint
445 //
446 CurrentA_PointVelocity = CurrentB_PointVelocity; // Got to B
447 CurrentB_PointVelocity = *ptrInterMedFT[depth];
448 SubStart_PointVelocity = CurrentA_PointVelocity;
449 ApproxIntersecPointV = GetChordFinderFor()
450 ->ApproxCurvePointV( CurrentA_PointVelocity,
451 CurrentB_PointVelocity,
452 CurrentE_Point,
454 restoredFullEndpoint = true;
455 ++restartB; // counter
456 }
457 }
458 } // Endif (Intersects_FB)
459 } // Endif (Intersects_AF)
460
461 // Ensure that the new endpoints are not further apart in space
462 // than on the curve due to different errors in the integration
463 //
464 G4double linDistSq, curveDist;
465 linDistSq = ( CurrentB_PointVelocity.GetPosition()
466 - CurrentA_PointVelocity.GetPosition() ).mag2();
467 curveDist = CurrentB_PointVelocity.GetCurveLength()
468 - CurrentA_PointVelocity.GetCurveLength();
469
470 // Change this condition for very strict parameters of propagation
471 //
472 if( curveDist*curveDist*(1+2* GetEpsilonStepFor()) < linDistSq )
473 {
474 // Re-integrate to obtain a new B
475 //
476 G4FieldTrack newEndPointFT=
477 ReEstimateEndpoint( CurrentA_PointVelocity,
478 CurrentB_PointVelocity,
479 linDistSq, // to avoid recalculation
480 curveDist );
481 G4FieldTrack oldPointVelB = CurrentB_PointVelocity;
482 CurrentB_PointVelocity = newEndPointFT;
483
484 if ( (fin_section_depth[depth]) // real final section
485 &&( first_section || ((Second_half)&&(depth==0)) ) )
486 {
487 recalculatedEndPoint = true;
488 IntersectedOrRecalculatedFT = newEndPointFT;
489 // So that we can return it, if it is the endpoint!
490 }
491 }
492 if( curveDist < 0.0 )
493 {
494 fVerboseLevel = 5; // Print out a maximum of information
495 printStatus( CurrentA_PointVelocity, CurrentB_PointVelocity,
496 -1.0, NewSafety, substep_no );
497 std::ostringstream message;
498 message << "Error in advancing propagation." << G4endl
499 << " Error in advancing propagation." << G4endl
500 << " Point A (start) is " << CurrentA_PointVelocity
501 << G4endl
502 << " Point B (end) is " << CurrentB_PointVelocity
503 << G4endl
504 << " Curve distance is " << curveDist << G4endl
505 << G4endl
506 << "The final curve point is not further along"
507 << " than the original!" << G4endl;
508
509 if( recalculatedEndPoint )
510 {
511 message << "Recalculation of EndPoint was called with fEpsStep= "
513 }
514 oldprc = G4cerr.precision(20);
515 message << " Point A (Curve start) is " << CurveStartPointVelocity
516 << G4endl
517 << " Point B (Curve end) is " << CurveEndPointVelocity
518 << G4endl
519 << " Point A (Current start) is " << CurrentA_PointVelocity
520 << G4endl
521 << " Point B (Current end) is " << CurrentB_PointVelocity
522 << G4endl
523 << " Point S (Sub start) is " << SubStart_PointVelocity
524 << G4endl
525 << " Point E (Trial Point) is " << CurrentE_Point
526 << G4endl
527 << " Old Point F(Intersection) is " << CurrentF_Point
528 << G4endl
529 << " New Point F(Intersection) is " << ApproxIntersecPointV
530 << G4endl
531 << " LocateIntersection parameters are : Substep no= "
532 << substep_no << G4endl
533 << " Substep depth no= "<< substep_no_p << " Depth= "
534 << depth << G4endl
535 << " Restarted no= "<< restartB << " Epsilon= "
536 << GetEpsilonStepFor() <<" DeltaInters= "
538 G4cerr.precision( oldprc );
539
540 G4Exception("G4BrentLocator::EstimateIntersectionPoint()",
541 "GeomNav0003", FatalException, message);
542 }
543
544 if( restoredFullEndpoint )
545 {
546 fin_section_depth[depth] = restoredFullEndpoint;
547 restoredFullEndpoint = false;
548 }
549 } // EndIf ( E is close enough to the curve, ie point F. )
550 // tests ChordAF_Vector.mag() <= maximum_lateral_displacement
551
552#ifdef G4DEBUG_LOCATE_INTERSECTION
553 G4int trigger_substepno_print= warn_substeps - 20 ;
554
555 if( substep_no >= trigger_substepno_print )
556 {
557 G4cout << "Difficulty in converging in "
558 << "G4BrentLocator::EstimateIntersectionPoint()"
559 << G4endl
560 << " Substep no = " << substep_no << G4endl;
561 if( substep_no == trigger_substepno_print )
562 {
563 printStatus( CurveStartPointVelocity, CurveEndPointVelocity,
564 -1.0, NewSafety, 0);
565 }
566 G4cout << " State of point A: ";
567 printStatus( CurrentA_PointVelocity, CurrentA_PointVelocity,
568 -1.0, NewSafety, substep_no-1, 0);
569 G4cout << " State of point B: ";
570 printStatus( CurrentA_PointVelocity, CurrentB_PointVelocity,
571 -1.0, NewSafety, substep_no);
572 }
573#endif
574 ++substep_no;
575 ++substep_no_p;
576
577 } while ( ( ! found_approximate_intersection )
578 && ( ! there_is_no_intersection )
579 && ( substep_no_p <= param_substeps) ); // UNTIL found or
580 // failed param substep
581 first_section = false;
582
583 if( (!found_approximate_intersection) && (!there_is_no_intersection) )
584 {
585 G4double did_len = std::abs( CurrentA_PointVelocity.GetCurveLength()
586 - SubStart_PointVelocity.GetCurveLength());
587 G4double all_len = std::abs( CurrentB_PointVelocity.GetCurveLength()
588 - SubStart_PointVelocity.GetCurveLength());
589
590 G4double stepLengthAB;
591 G4ThreeVector PointGe;
592
593 // Check if progress is too slow and if it possible to go deeper,
594 // then halve the step if so
595 //
596 if ( ( did_len < fraction_done*all_len )
597 && (depth < max_depth) && (!sub_final_section) )
598 {
599 Second_half=false;
600 ++depth;
601
602 G4double Sub_len = (all_len-did_len)/(2.);
603 G4FieldTrack start = CurrentA_PointVelocity;
604 auto integrDriver =
606 integrDriver->AccurateAdvance(start, Sub_len, GetEpsilonStepFor());
607 *ptrInterMedFT[depth] = start;
608 CurrentB_PointVelocity = *ptrInterMedFT[depth];
609
610 // Adjust 'SubStartPoint' to calculate the 'did_length' in next loop
611 //
612 SubStart_PointVelocity = CurrentA_PointVelocity;
613
614 // Find new trial intersection point needed at start of the loop
615 //
616 G4ThreeVector Point_A = CurrentA_PointVelocity.GetPosition();
617 G4ThreeVector SubE_point = CurrentB_PointVelocity.GetPosition();
618
620 G4bool Intersects_AB = IntersectChord(Point_A, SubE_point,
621 NewSafety, fPreviousSafety,
622 fPreviousSftOrigin,stepLengthAB,
623 PointGe);
624 if( Intersects_AB )
625 {
626 last_AF_intersection = Intersects_AB;
627 CurrentE_Point = PointGe;
628 fin_section_depth[depth] = true;
629
630 // Need to recalculate the Exit Normal at the new PointG
631 //
632 G4bool validNormalAB;
633 NormalAtEntry = GetSurfaceNormal( PointGe, validNormalAB );
634 validNormalAtE = validNormalAB;
635 }
636 else
637 {
638 // No intersection found for first part of curve
639 // (CurrentA,InterMedPoint[depth]). Go to the second part
640 //
641 Second_half = true;
642 }
643 } // if did_len
644
645 if( (Second_half)&&(depth!=0) )
646 {
647 // Second part of curve (InterMed[depth],Intermed[depth-1]) )
648 // On the depth-1 level normally we are on the 'second_half'
649
650 Second_half = true;
651
652 // Find new trial intersection point needed at start of the loop
653 //
654 SubStart_PointVelocity = *ptrInterMedFT[depth];
655 CurrentA_PointVelocity = *ptrInterMedFT[depth];
656 CurrentB_PointVelocity = *ptrInterMedFT[depth-1];
657 // Ensure that the new endpoints are not further apart in space
658 // than on the curve due to different errors in the integration
659 //
660 G4double linDistSq, curveDist;
661 linDistSq = ( CurrentB_PointVelocity.GetPosition()
662 - CurrentA_PointVelocity.GetPosition() ).mag2();
663 curveDist = CurrentB_PointVelocity.GetCurveLength()
664 - CurrentA_PointVelocity.GetCurveLength();
665 if( curveDist*curveDist*(1+2*GetEpsilonStepFor() ) < linDistSq )
666 {
667 // Re-integrate to obtain a new B
668 //
669 G4FieldTrack newEndPointFT =
670 ReEstimateEndpoint( CurrentA_PointVelocity,
671 CurrentB_PointVelocity,
672 linDistSq, // to avoid recalculation
673 curveDist );
674 G4FieldTrack oldPointVelB = CurrentB_PointVelocity;
675 CurrentB_PointVelocity = newEndPointFT;
676 if ( depth==1 )
677 {
678 recalculatedEndPoint = true;
679 IntersectedOrRecalculatedFT = newEndPointFT;
680 // So that we can return it, if it is the endpoint!
681 }
682 }
683
684
685 G4ThreeVector Point_A = CurrentA_PointVelocity.GetPosition();
686 G4ThreeVector SubE_point = CurrentB_PointVelocity.GetPosition();
688 G4bool Intersects_AB = IntersectChord(Point_A, SubE_point, NewSafety,
690 fPreviousSftOrigin,stepLengthAB, PointGe);
691 if( Intersects_AB )
692 {
693 last_AF_intersection = Intersects_AB;
694 CurrentE_Point = PointGe;
695
696 G4bool validNormalAB;
697 NormalAtEntry = GetSurfaceNormal( PointGe, validNormalAB );
698 validNormalAtE = validNormalAB;
699 }
700
701 depth--;
702 fin_section_depth[depth]=true;
703 }
704 } // if(!found_aproximate_intersection)
705
706 } while ( ( ! found_approximate_intersection )
707 && ( ! there_is_no_intersection )
708 && ( substep_no <= max_substeps) ); // UNTIL found or failed
709
710 if( substep_no > max_no_seen )
711 {
712 max_no_seen = substep_no;
713#ifdef G4DEBUG_LOCATE_INTERSECTION
714 if( max_no_seen > warn_substeps )
715 {
716 trigger_substepno_print = max_no_seen-20; // Want to see last 20 steps
717 }
718#endif
719 }
720
721 if( ( substep_no >= max_substeps)
722 && !there_is_no_intersection
723 && !found_approximate_intersection )
724 {
725 G4cout << "ERROR - G4BrentLocator::EstimateIntersectionPoint()" << G4endl
726 << " Start and end-point of requested Step:" << G4endl;
727 printStatus( CurveStartPointVelocity, CurveEndPointVelocity,
728 -1.0, NewSafety, 0);
729 G4cout << " Start and end-point of current Sub-Step:" << G4endl;
730 printStatus( CurrentA_PointVelocity, CurrentA_PointVelocity,
731 -1.0, NewSafety, substep_no-1);
732 printStatus( CurrentA_PointVelocity, CurrentB_PointVelocity,
733 -1.0, NewSafety, substep_no);
734 std::ostringstream message;
735 message << "Too many substeps!" << G4endl
736 << " Convergence is requiring too many substeps: "
737 << substep_no << G4endl
738 << " Abandoning effort to intersect. " << G4endl
739 << " Found intersection = "
740 << found_approximate_intersection << G4endl
741 << " Intersection exists = "
742 << !there_is_no_intersection << G4endl;
743 oldprc = G4cout.precision( 10 );
744 G4double done_len = CurrentA_PointVelocity.GetCurveLength();
745 G4double full_len = CurveEndPointVelocity.GetCurveLength();
746 message << " Undertaken only length: " << done_len
747 << " out of " << full_len << " required." << G4endl
748 << " Remaining length = " << full_len - done_len;
749 G4cout.precision( oldprc );
750
751 G4Exception("G4BrentLocator::EstimateIntersectionPoint()",
752 "GeomNav0003", FatalException, message);
753 }
754 else if( substep_no >= warn_substeps )
755 {
756 oldprc= G4cout.precision( 10 );
757 std::ostringstream message;
758 message << "Many substeps while trying to locate intersection."
759 << G4endl
760 << " Undertaken length: "
761 << CurrentB_PointVelocity.GetCurveLength()
762 << " - Needed: " << substep_no << " substeps." << G4endl
763 << " Warning level = " << warn_substeps
764 << " and maximum substeps = " << max_substeps;
765 G4Exception("G4BrentLocator::EstimateIntersectionPoint()",
766 "GeomNav1002", JustWarning, message);
767 G4cout.precision( oldprc );
768 }
769 return !there_is_no_intersection; // Success or failure
770}
@ FatalException
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
G4FieldTrack ApproxCurvePointV(const G4FieldTrack &curveAPointVelocity, const G4FieldTrack &curveBPointVelocity, const G4ThreeVector &currentEPoint, G4double epsStep)
G4FieldTrack ApproxCurvePointS(const G4FieldTrack &curveAPointVelocity, const G4FieldTrack &curveBPointVelocity, const G4FieldTrack &ApproxCurveV, const G4ThreeVector &currentEPoint, const G4ThreeVector &currentFPoint, const G4ThreeVector &PointG, G4bool first, G4double epsStep)
G4VIntegrationDriver * GetIntegrationDriver()
virtual G4bool AccurateAdvance(G4FieldTrack &track, G4double hstep, G4double eps, G4double hinitial=0)=0
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)
G4double GetEpsilonStepFor()
G4bool GetAdjustementOfFoundIntersection()
G4double GetDeltaIntersectionFor()
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)
static constexpr double mm
Definition: SystemOfUnits.h:96
T sqr(const T &x)
Definition: templates.hh:128
#define G4ThreadLocal
Definition: tls.hh:77

References G4VIntegrationDriver::AccurateAdvance(), G4VIntersectionLocator::AdjustmentOfFoundIntersection(), G4ChordFinder::ApproxCurvePointS(), G4ChordFinder::ApproxCurvePointV(), CLHEP::Hep3Vector::dot(), FatalException, G4VIntersectionLocator::fiDeltaIntersection, fPreviousSafety, fPreviousSftOrigin, G4VIntersectionLocator::fVerboseLevel, G4cerr, G4cout, G4endl, G4Exception(), G4ThreadLocal, G4VIntersectionLocator::GetAdjustementOfFoundIntersection(), G4VIntersectionLocator::GetChordFinderFor(), G4FieldTrack::GetCurveLength(), G4VIntersectionLocator::GetDeltaIntersectionFor(), G4VIntersectionLocator::GetEpsilonStepFor(), G4ChordFinder::GetIntegrationDriver(), G4FieldTrack::GetMomentumDir(), G4FieldTrack::GetMomentumDirection(), G4VIntersectionLocator::GetNavigatorFor(), G4FieldTrack::GetPosition(), G4VIntersectionLocator::GetSurfaceNormal(), G4VIntersectionLocator::IntersectChord(), JustWarning, G4VIntersectionLocator::kCarTolerance, G4Navigator::LocateGlobalPointWithinVolume(), CLHEP::Hep3Vector::mag2(), max_depth, CLHEP::mm, G4VIntersectionLocator::printStatus(), ptrInterMedFT, G4VIntersectionLocator::ReEstimateEndpoint(), G4VIntersectionLocator::ReportTrialStep(), G4FieldTrack::SetPosition(), and sqr().

◆ GetAdjustementOfFoundIntersection()

G4bool G4VIntersectionLocator::GetAdjustementOfFoundIntersection ( )
inlineinherited

◆ GetCheckMode()

G4bool G4VIntersectionLocator::GetCheckMode ( )
inlineinherited

◆ GetChordFinderFor()

G4ChordFinder * G4VIntersectionLocator::GetChordFinderFor ( )
inlineinherited

◆ GetDeltaIntersectionFor()

G4double G4VIntersectionLocator::GetDeltaIntersectionFor ( )
inlineinherited

◆ GetEpsilonStepFor()

G4double G4VIntersectionLocator::GetEpsilonStepFor ( )
inlineinherited

◆ GetGlobalSurfaceNormal()

G4ThreeVector G4VIntersectionLocator::GetGlobalSurfaceNormal ( const G4ThreeVector CurrentE_Point,
G4bool validNormal 
)
protectedinherited

Definition at line 563 of file G4VIntersectionLocator.cc.

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}
static constexpr double perThousand
Definition: G4SIunits.hh:326
G4ThreeVector NetTranslation() const
G4RotationMatrix NetRotation() const
G4ThreeVector TransformAxis(const G4ThreeVector &axis) const
const G4AffineTransform GetLocalToGlobalTransform() const
G4ThreeVector GetLocalSurfaceNormal(const G4ThreeVector &CurrentE_Point, G4bool &validNormal)

References G4VIntersectionLocator::fHelpingNavigator, G4endl, G4Exception(), G4VIntersectionLocator::GetLocalSurfaceNormal(), G4Navigator::GetLocalToGlobalTransform(), JustWarning, CLHEP::Hep3Vector::mag2(), G4AffineTransform::NetRotation(), G4AffineTransform::NetTranslation(), perThousand, and G4AffineTransform::TransformAxis().

Referenced by G4VIntersectionLocator::AdjustmentOfFoundIntersection().

◆ GetLastSurfaceNormal()

G4ThreeVector G4VIntersectionLocator::GetLastSurfaceNormal ( const G4ThreeVector intersectPoint,
G4bool validNormal 
) const
privateinherited

Definition at line 602 of file G4VIntersectionLocator.cc.

605{
606 G4ThreeVector normalVec;
607 G4bool validNorm;
608 normalVec = fiNavigator->GetGlobalExitNormal( intersectPoint, &validNorm );
609 normalIsValid = validNorm;
610
611 return normalVec;
612}
virtual G4ThreeVector GetGlobalExitNormal(const G4ThreeVector &point, G4bool *valid)

References G4VIntersectionLocator::fiNavigator, and G4Navigator::GetGlobalExitNormal().

Referenced by G4VIntersectionLocator::GetSurfaceNormal().

◆ GetLocalSurfaceNormal()

G4ThreeVector G4VIntersectionLocator::GetLocalSurfaceNormal ( const G4ThreeVector CurrentE_Point,
G4bool validNormal 
)
privateinherited

Definition at line 384 of file G4VIntersectionLocator.cc.

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}
CLHEP::Hep3Vector G4ThreeVector
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4VSolid * GetSolid() const
const G4AffineTransform & GetTopTransform() const
G4TouchableHistory * CreateTouchableHistory() const
virtual G4VPhysicalVolume * LocateGlobalPointAndSetup(const G4ThreeVector &point, const G4ThreeVector *direction=nullptr, const G4bool pRelativeSearch=true, const G4bool ignoreDirection=true)
Definition: G4Navigator.cc:132
void SetWorldVolume(G4VPhysicalVolume *pWorld)
const G4NavigationHistory * GetHistory() const
G4TouchableHistory * fpTouchable
G4LogicalVolume * GetLogicalVolume() const
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
@ kSurface
Definition: geomdefs.hh:69
static constexpr double perThousand

References G4Navigator::CreateTouchableHistory(), G4VSolid::DistanceToOut(), G4VIntersectionLocator::fHelpingNavigator, G4VIntersectionLocator::fpTouchable, G4cerr, G4endl, G4TouchableHistory::GetHistory(), G4VPhysicalVolume::GetLogicalVolume(), G4VIntersectionLocator::GetNavigatorFor(), G4LogicalVolume::GetSolid(), G4NavigationHistory::GetTopTransform(), G4VSolid::Inside(), G4VIntersectionLocator::kCarTolerance, kSurface, G4Navigator::LocateGlobalPointAndSetup(), CLHEP::Hep3Vector::mag(), CLHEP::Hep3Vector::mag2(), CLHEP::perThousand, G4Navigator::SetWorldVolume(), G4VSolid::SurfaceNormal(), and G4AffineTransform::TransformPoint().

Referenced by G4VIntersectionLocator::GetGlobalSurfaceNormal().

◆ GetNavigatorFor()

G4Navigator * G4VIntersectionLocator::GetNavigatorFor ( )
inlineinherited

◆ GetSurfaceNormal()

G4ThreeVector G4VIntersectionLocator::GetSurfaceNormal ( const G4ThreeVector CurrentInt_Point,
G4bool validNormal 
)
protectedinherited

Definition at line 518 of file G4VIntersectionLocator.cc.

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}
double mag() const
G4ThreeVector GetLastSurfaceNormal(const G4ThreeVector &intersectPoint, G4bool &validNormal) const

References G4endl, G4Exception(), G4VIntersectionLocator::GetLastSurfaceNormal(), JustWarning, CLHEP::Hep3Vector::mag(), CLHEP::Hep3Vector::mag2(), and perThousand.

Referenced by EstimateIntersectionPoint(), G4MultiLevelLocator::EstimateIntersectionPoint(), and G4SimpleLocator::EstimateIntersectionPoint().

◆ GetVerboseFor()

G4int G4VIntersectionLocator::GetVerboseFor ( )
inlineinherited

◆ IntersectChord()

G4bool G4VIntersectionLocator::IntersectChord ( const G4ThreeVector StartPointA,
const G4ThreeVector EndPointB,
G4double NewSafety,
G4double PreviousSafety,
G4ThreeVector PreviousSftOrigin,
G4double LinearStepLength,
G4ThreeVector IntersectionPoint,
G4bool calledNavigator = nullptr 
)
inlineinherited

◆ LocateGlobalPointWithinVolumeAndCheck()

G4bool G4VIntersectionLocator::LocateGlobalPointWithinVolumeAndCheck ( const G4ThreeVector pos)
protectedinherited

Definition at line 674 of file G4VIntersectionLocator.cc.

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}
void CheckMode(G4bool mode)
virtual G4TouchableHistoryHandle CreateTouchableHistoryHandle() const
G4bool IsCheckModeActive() const
const G4AffineTransform & GetGlobalToLocalTransform() const
virtual G4int GetCopyNo() const =0
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const =0
EInside
Definition: geomdefs.hh:67
@ kInside
Definition: geomdefs.hh:70
G4bool transform(G4String &input, const G4String &type)

References G4Navigator::CheckMode(), G4Navigator::CreateTouchableHistoryHandle(), G4VSolid::DistanceToIn(), G4VIntersectionLocator::fCheckMode, G4endl, G4Exception(), G4VPhysicalVolume::GetCopyNo(), G4Navigator::GetGlobalToLocalTransform(), G4VIntersectionLocator::GetNavigatorFor(), G4VSolid::Inside(), G4Navigator::IsCheckModeActive(), JustWarning, kInside, kSurface, G4Navigator::LocateGlobalPointAndSetup(), G4Navigator::LocateGlobalPointWithinVolume(), and G4coutFormatters::anonymous_namespace{G4coutFormatters.cc}::transform().

Referenced by G4VIntersectionLocator::LocateGlobalPointWithinVolumeCheckAndReport().

◆ LocateGlobalPointWithinVolumeCheckAndReport()

void G4VIntersectionLocator::LocateGlobalPointWithinVolumeCheckAndReport ( const G4ThreeVector pos,
const G4String CodeLocationInfo,
G4int  CheckMode 
)
protectedinherited

Definition at line 736 of file G4VIntersectionLocator.cc.

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}
G4bool LocateGlobalPointWithinVolumeAndCheck(const G4ThreeVector &pos)

References G4endl, G4Exception(), G4VIntersectionLocator::GetCheckMode(), JustWarning, G4VIntersectionLocator::LocateGlobalPointWithinVolumeAndCheck(), and G4VIntersectionLocator::SetCheckMode().

◆ printStatus() [1/2]

void G4VIntersectionLocator::printStatus ( const G4FieldTrack startFT,
const G4FieldTrack currentFT,
G4double  requestStep,
G4double  safety,
G4int  stepNum 
)
inherited

◆ printStatus() [2/2]

void G4VIntersectionLocator::printStatus ( const G4FieldTrack startFT,
const G4FieldTrack currentFT,
G4double  requestStep,
G4double  safety,
G4int  stepNum,
std::ostream &  oss,
G4int  verboseLevel 
)
staticinherited

Definition at line 91 of file G4VIntersectionLocator.cc.

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}
double z() const
double x() const
double y() const

References G4endl, G4FieldTrack::GetCurveLength(), G4FieldTrack::GetMomentum(), G4FieldTrack::GetMomentumDir(), G4FieldTrack::GetPosition(), CLHEP::Hep3Vector::mag(), G4VIntersectionLocator::printStatus(), CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

◆ ReEstimateEndpoint()

G4FieldTrack G4VIntersectionLocator::ReEstimateEndpoint ( const G4FieldTrack CurrentStateA,
const G4FieldTrack EstimtdEndStateB,
G4double  linearDistSq,
G4double  curveDist 
)
protectedinherited

Definition at line 188 of file G4VIntersectionLocator.cc.

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}

References G4VIntersectionLocator::fVerboseLevel, G4cout, G4endl, G4Exception(), G4VIntersectionLocator::GetChordFinderFor(), G4FieldTrack::GetCurveLength(), G4VIntersectionLocator::GetEpsilonStepFor(), G4ChordFinder::GetIntegrationDriver(), G4FieldTrack::GetPosition(), JustWarning, G4VIntersectionLocator::kCarTolerance, and CLHEP::Hep3Vector::mag().

Referenced by G4VIntersectionLocator::CheckAndReEstimateEndpoint(), EstimateIntersectionPoint(), and G4SimpleLocator::EstimateIntersectionPoint().

◆ ReportImmediateHit()

void G4VIntersectionLocator::ReportImmediateHit ( const char *  MethodName,
const G4ThreeVector StartPosition,
const G4ThreeVector TrialPoint,
G4double  tolerance,
unsigned long int  numCalls 
)
protectedinherited

Definition at line 845 of file G4VIntersectionLocator.cc.

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()
void Register(T *inst)
Definition: G4AutoDelete.hh:65
#define DBL_MAX
Definition: templates.hh:62

References DBL_MAX, G4cout, G4endl, G4ThreadLocal, and G4AutoDelete::Register().

Referenced by G4MultiLevelLocator::EstimateIntersectionPoint().

◆ ReportProgress()

void G4VIntersectionLocator::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 
)
protectedinherited

Definition at line 813 of file G4VIntersectionLocator.cc.

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}

References G4endl, and G4VIntersectionLocator::printStatus().

Referenced by G4MultiLevelLocator::EstimateIntersectionPoint().

◆ ReportReversedPoints()

void G4VIntersectionLocator::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 
)
protectedinherited

Definition at line 761 of file G4VIntersectionLocator.cc.

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}
static constexpr double millimeter
Definition: SystemOfUnits.h:63

References G4endl, G4FieldTrack::GetCurveLength(), CLHEP::millimeter, and G4VIntersectionLocator::printStatus().

Referenced by G4MultiLevelLocator::EstimateIntersectionPoint().

◆ ReportTrialStep()

void G4VIntersectionLocator::ReportTrialStep ( G4int  step_no,
const G4ThreeVector ChordAB_v,
const G4ThreeVector ChordEF_v,
const G4ThreeVector NewMomentumDir,
const G4ThreeVector NormalAtEntry,
G4bool  validNormal 
)
protectedinherited

Definition at line 618 of file G4VIntersectionLocator.cc.

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}

References CLHEP::Hep3Vector::dot(), G4cout, G4endl, G4Exception(), JustWarning, CLHEP::Hep3Vector::mag(), CLHEP::Hep3Vector::mag2(), and perThousand.

Referenced by EstimateIntersectionPoint(), G4MultiLevelLocator::EstimateIntersectionPoint(), and G4SimpleLocator::EstimateIntersectionPoint().

◆ SetCheckMode()

void G4VIntersectionLocator::SetCheckMode ( G4bool  value)
inlineinherited

◆ SetChordFinderFor()

void G4VIntersectionLocator::SetChordFinderFor ( G4ChordFinder fCFinder)
inlineinherited

◆ SetDeltaIntersectionFor()

void G4VIntersectionLocator::SetDeltaIntersectionFor ( G4double  deltaIntersection)
inlineinherited

◆ SetEpsilonStepFor()

void G4VIntersectionLocator::SetEpsilonStepFor ( G4double  EpsilonStep)
inlineinherited

◆ SetNavigatorFor()

void G4VIntersectionLocator::SetNavigatorFor ( G4Navigator fNavigator)
inlineinherited

◆ SetSafetyParametersFor()

void G4VIntersectionLocator::SetSafetyParametersFor ( G4bool  UseSafety)
inlineinherited

◆ SetVerboseFor()

void G4VIntersectionLocator::SetVerboseFor ( G4int  fVerbose)
inlineinherited

Field Documentation

◆ fCheckMode

G4bool G4VIntersectionLocator::fCheckMode = false
protectedinherited

◆ fHelpingNavigator

G4Navigator* G4VIntersectionLocator::fHelpingNavigator
protectedinherited

◆ fiChordFinder

G4ChordFinder* G4VIntersectionLocator::fiChordFinder = nullptr
protectedinherited

Definition at line 261 of file G4VIntersectionLocator.hh.

◆ fiDeltaIntersection

G4double G4VIntersectionLocator::fiDeltaIntersection = -1.0
protectedinherited

◆ fiEpsilonStep

G4double G4VIntersectionLocator::fiEpsilonStep = -1.0
protectedinherited

◆ fiNavigator

G4Navigator* G4VIntersectionLocator::fiNavigator
protectedinherited

◆ fiUseSafety

G4bool G4VIntersectionLocator::fiUseSafety = false
protectedinherited

Definition at line 257 of file G4VIntersectionLocator.hh.

◆ fpTouchable

G4TouchableHistory* G4VIntersectionLocator::fpTouchable = nullptr
protectedinherited

◆ fUseNormalCorrection

G4bool G4VIntersectionLocator::fUseNormalCorrection = false
protectedinherited

◆ fVerboseLevel

G4int G4VIntersectionLocator::fVerboseLevel = 0
protectedinherited

◆ kCarTolerance

G4double G4VIntersectionLocator::kCarTolerance
protectedinherited

◆ max_depth

const G4int G4BrentLocator::max_depth = 4
staticprivate

Definition at line 71 of file G4BrentLocator.hh.

Referenced by EstimateIntersectionPoint(), G4BrentLocator(), and ~G4BrentLocator().

◆ ptrInterMedFT

G4FieldTrack* G4BrentLocator::ptrInterMedFT[max_depth+1]
private

Definition at line 72 of file G4BrentLocator.hh.

Referenced by EstimateIntersectionPoint(), G4BrentLocator(), and ~G4BrentLocator().


The documentation for this class was generated from the following files: