Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PropagatorInField.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 G4PropagatorInField Implementation
27//
28// Author: John Apostolakis (CERN), 25 October 1996
29// ---------------------------------------------------------------------------
30
31#include <iomanip>
32
34#include "G4ios.hh"
35#include "G4SystemOfUnits.hh"
36#include "G4ThreeVector.hh"
37#include "G4Material.hh"
38#include "G4VPhysicalVolume.hh"
39#include "G4Navigator.hh"
42#include "G4ChordFinder.hh"
44
45
46// ---------------------------------------------------------------------------
47// Constructors and destructor
48//
50 G4FieldManager* detectorFieldMgr,
51 G4VIntersectionLocator* vLocator )
52 : fDetectorFieldMgr(detectorFieldMgr),
53 fNavigator(theNavigator),
54 fCurrentFieldMgr(detectorFieldMgr),
55 End_PointAndTangent(G4ThreeVector(0.,0.,0.),
56 G4ThreeVector(0.,0.,0.),0.0,0.0,0.0,0.0,0.0)
57{
58 fEpsilonStep = (fDetectorFieldMgr != nullptr)
59 ? fDetectorFieldMgr->GetMaximumEpsilonStep() : 1.0e-5;
60
61
62 fPreviousSftOrigin = G4ThreeVector(0.,0.,0.);
64 fZeroStepThreshold = std::max( 1.0e5 * kCarTolerance, 1.0e-1 * micrometer );
65
66 fLargestAcceptableStep = 100.0 * meter; // Reduced from 1000.0 * meter
67 fMaxStepSizeMultiplier= 0.1 ; // 0.1 in git (larger for tests.) // Reduced from 100;
68 fMinBigDistance= 100. * CLHEP::mm;
69#ifdef G4DEBUG_FIELD
70 G4cout << " PiF: Zero Step Threshold set to "
71 << fZeroStepThreshold / millimeter
72 << " mm." << G4endl;
73 G4cout << " PiF: Value of kCarTolerance = "
74 << kCarTolerance / millimeter
75 << " mm. " << G4endl;
76 fVerboseLevel = 2;
77 fVerbTracePiF = true;
78#endif
79
80 // Defining Intersection Locator and his parameters
81 if ( vLocator == nullptr )
82 {
83 fIntersectionLocator = new G4MultiLevelLocator(theNavigator);
84 fAllocatedLocator = true;
85 }
86 else
87 {
88 fIntersectionLocator = vLocator;
89 fAllocatedLocator = false;
90 }
91 RefreshIntersectionLocator(); // Copy all relevant parameters
92}
93
94// ---------------------------------------------------------------------------
95//
97{
98 if(fAllocatedLocator) { delete fIntersectionLocator; }
99}
100
101// ---------------------------------------------------------------------------
102// Update the IntersectionLocator with current parameters
103//
105{
106 fIntersectionLocator->SetEpsilonStepFor(fEpsilonStep);
107 fIntersectionLocator->SetDeltaIntersectionFor(fCurrentFieldMgr->GetDeltaIntersection());
108 fIntersectionLocator->SetChordFinderFor(GetChordFinder());
109 fIntersectionLocator->SetSafetyParametersFor( fUseSafetyForOptimisation);
110}
111
112// ---------------------------------------------------------------------------
113// Compute the next geometric Step
114//
116 G4FieldTrack& pFieldTrack,
117 G4double CurrentProposedStepLength,
118 G4double& currentSafety, // IN/OUT
119 G4VPhysicalVolume* pPhysVol,
120 G4bool canRelaxDeltaChord)
121{
122 GetChordFinder()->OnComputeStep(&pFieldTrack);
123 const G4double deltaChord = GetChordFinder()->GetDeltaChord();
124
125 // If CurrentProposedStepLength is too small for finding Chords
126 // then return with no action (for now - TODO: some action)
127 //
128 const char* methodName = "G4PropagatorInField::ComputeStep";
129 if (CurrentProposedStepLength<kCarTolerance)
130 {
131 return kInfinity;
132 }
133
134 // Introducing smooth trajectory display (jacek 01/11/2002)
135 //
136 if (fpTrajectoryFilter != nullptr)
137 {
138 fpTrajectoryFilter->CreateNewTrajectorySegment();
139 }
140
141 fFirstStepInVolume = fNewTrack ? true : fLastStepInVolume;
142 fLastStepInVolume = false;
143 fNewTrack = false;
144
145 if( fVerboseLevel > 2 )
146 {
147 G4cout << methodName << " called" << G4endl;
148 G4cout << " Starting FT: " << pFieldTrack;
149 G4cout << " Requested length = " << CurrentProposedStepLength << G4endl;
150 G4cout << " PhysVol = ";
151 if( pPhysVol != nullptr )
152 {
153 G4cout << pPhysVol->GetName() << G4endl;
154 }
155 else
156 {
157 G4cout << " N/A ";
158 }
159 G4cout << G4endl;
160 }
161
162 // Parameters for adaptive Runge-Kutta integration
163
164 G4double h_TrialStepSize; // 1st Step Size
165 G4double TruePathLength = CurrentProposedStepLength;
166 G4double StepTaken = 0.0;
167 G4double s_length_taken, epsilon;
168 G4bool intersects;
169 G4bool first_substep = true;
170
171 G4double NewSafety;
172 fParticleIsLooping = false;
173
174 // If not yet done,
175 // Set the field manager to the local one if the volume has one,
176 // or to the global one if not
177 //
178 if( !fSetFieldMgr )
179 {
180 fCurrentFieldMgr = FindAndSetFieldManager( pPhysVol );
181 }
182 fSetFieldMgr = false; // For next call, the field manager must be set again
183
184 G4FieldTrack CurrentState(pFieldTrack);
185 G4FieldTrack OriginalState = CurrentState;
186
187 // If the Step length is "infinite", then an approximate-maximum Step
188 // length (used to calculate the relative accuracy) must be guessed
189 //
190 if( CurrentProposedStepLength >= fLargestAcceptableStep )
191 {
192 G4ThreeVector StartPointA, VelocityUnit;
193 StartPointA = pFieldTrack.GetPosition();
194 VelocityUnit = pFieldTrack.GetMomentumDir();
195
196 G4double trialProposedStep = fMaxStepSizeMultiplier * ( fMinBigDistance +
197 fNavigator->GetWorldVolume()->GetLogicalVolume()->
198 GetSolid()->DistanceToOut(StartPointA, VelocityUnit) );
199 CurrentProposedStepLength = std::min( trialProposedStep,
200 fLargestAcceptableStep );
201 }
202 epsilon = fCurrentFieldMgr->GetDeltaOneStep() / CurrentProposedStepLength;
203 G4double epsilonMin= fCurrentFieldMgr->GetMinimumEpsilonStep();
204 G4double epsilonMax= fCurrentFieldMgr->GetMaximumEpsilonStep();
205 if( epsilon < epsilonMin ) { epsilon = epsilonMin; }
206 if( epsilon > epsilonMax ) { epsilon = epsilonMax; }
208
209 // Values for Intersection Locator has to be updated on each call for the
210 // case that CurrentFieldManager has changed from the one of previous step
211 //
213
214 // Shorten the proposed step in case of earlier problems (zero steps)
215 //
216 if( fNoZeroStep > fActionThreshold_NoZeroSteps )
217 {
218 G4double stepTrial;
219
220 stepTrial = fFull_CurveLen_of_LastAttempt;
221 if( (stepTrial <= 0.0) && (fLast_ProposedStepLength > 0.0) )
222 {
223 stepTrial = fLast_ProposedStepLength;
224 }
225
226 G4double decreaseFactor = 0.9; // Unused default
227 if( (fNoZeroStep < fSevereActionThreshold_NoZeroSteps)
228 && (stepTrial > 100.0*fZeroStepThreshold) )
229 {
230 // Attempt quick convergence
231 //
232 decreaseFactor= 0.25;
233 }
234 else
235 {
236 // We are in significant difficulties, probably at a boundary that
237 // is either geometrically sharp or between very different materials.
238 // Careful decreases to cope with tolerance are required
239 //
240 if( stepTrial > 100.0*fZeroStepThreshold ) {
241 decreaseFactor = 0.35; // Try decreasing slower
242 } else if( stepTrial > 30.0*fZeroStepThreshold ) {
243 decreaseFactor= 0.5; // Try yet slower decrease
244 } else if( stepTrial > 10.0*fZeroStepThreshold ) {
245 decreaseFactor= 0.75; // Try even slower decreases
246 } else {
247 decreaseFactor= 0.9; // Try very slow decreases
248 }
249 }
250 stepTrial *= decreaseFactor;
251
252#ifdef G4DEBUG_FIELD
253 if( fVerboseLevel > 2
254 || (fNoZeroStep >= fSevereActionThreshold_NoZeroSteps) )
255 {
256 G4cerr << " " << methodName
257 << " Decreasing step after " << fNoZeroStep << " zero steps "
258 << " - in volume " << pPhysVol;
259 if( pPhysVol )
260 G4cerr << " with name " << pPhysVol->GetName();
261 else
262 G4cerr << " i.e. *unknown* volume.";
263 G4cerr << G4endl;
264 PrintStepLengthDiagnostic(CurrentProposedStepLength, decreaseFactor,
265 stepTrial, pFieldTrack);
266 }
267#endif
268 if( stepTrial == 0.0 ) // Change to make it < 0.1 * kCarTolerance ??
269 {
270 std::ostringstream message;
271 message << "Particle abandoned due to lack of progress in field."
272 << G4endl
273 << " Properties : " << pFieldTrack << G4endl
274 << " Attempting a zero step = " << stepTrial << G4endl
275 << " while attempting to progress after " << fNoZeroStep
276 << " trial steps. Will abandon step.";
277 G4Exception(methodName, "GeomNav1002", JustWarning, message);
278 fParticleIsLooping = true;
279 return 0; // = stepTrial;
280 }
281 if( stepTrial < CurrentProposedStepLength )
282 {
283 CurrentProposedStepLength = stepTrial;
284 }
285 }
286 fLast_ProposedStepLength = CurrentProposedStepLength;
287
288 G4int do_loop_count = 0;
289 do // Loop checking, 07.10.2016, JA
290 {
291 G4FieldTrack SubStepStartState = CurrentState;
292 G4ThreeVector SubStartPoint = CurrentState.GetPosition();
293
294 if(!first_substep)
295 {
296 if( fVerboseLevel > 4 )
297 {
298 G4cout << " PiF: Calling Nav/Locate Global Point within-Volume "
299 << G4endl;
300 }
301 fNavigator->LocateGlobalPointWithinVolume( SubStartPoint );
302 }
303
304 // How far to attempt to move the particle !
305 //
306 h_TrialStepSize = CurrentProposedStepLength - StepTaken;
307
308 if (canRelaxDeltaChord &&
309 fIncreaseChordDistanceThreshold > 0 &&
310 do_loop_count > fIncreaseChordDistanceThreshold &&
311 do_loop_count % fIncreaseChordDistanceThreshold == 0)
312 {
314 GetChordFinder()->GetDeltaChord() * 2.0
315 );
316 }
317
318 // Integrate as far as "chord miss" rule allows.
319 //
320 s_length_taken = GetChordFinder()->AdvanceChordLimited(
321 CurrentState, // Position & velocity
322 h_TrialStepSize,
323 fEpsilonStep,
324 fPreviousSftOrigin,
325 fPreviousSafety );
326 // CurrentState is now updated with the final position and velocity
327
328 fFull_CurveLen_of_LastAttempt = s_length_taken;
329
330 G4ThreeVector EndPointB = CurrentState.GetPosition();
331 G4ThreeVector InterSectionPointE;
332 G4double LinearStepLength;
333
334 // Intersect chord AB with geometry
335 //
336 intersects= IntersectChord( SubStartPoint, EndPointB,
337 NewSafety, LinearStepLength,
338 InterSectionPointE );
339 // E <- Intersection Point of chord AB and either volume A's surface
340 // or a daughter volume's surface ..
341
342 if( first_substep )
343 {
344 currentSafety = NewSafety;
345 } // Updating safety in other steps is potential future extention
346
347 if( intersects )
348 {
349 G4FieldTrack IntersectPointVelct_G(CurrentState); // FT-Def-Construct
350
351 // Find the intersection point of AB true path with the surface
352 // of vol(A), if it exists. Start with point E as first "estimate".
353 G4bool recalculatedEndPt = false;
354
355 G4bool found_intersection = fIntersectionLocator->
356 EstimateIntersectionPoint( SubStepStartState, CurrentState,
357 InterSectionPointE, IntersectPointVelct_G,
358 recalculatedEndPt, fPreviousSafety,
359 fPreviousSftOrigin);
360 intersects = found_intersection;
361 if( found_intersection )
362 {
363 End_PointAndTangent= IntersectPointVelct_G; // G is our EndPoint ...
364 StepTaken = TruePathLength = IntersectPointVelct_G.GetCurveLength()
365 - OriginalState.GetCurveLength();
366 }
367 else
368 {
369 // Either "minor" chords do not intersect
370 // or else stopped (due to too many steps)
371 //
372 if( recalculatedEndPt )
373 {
374 G4double endAchieved = IntersectPointVelct_G.GetCurveLength();
375 G4double endExpected = CurrentState.GetCurveLength();
376
377 // Detect failure - due to too many steps
378 G4bool shortEnd = endAchieved
379 < (endExpected*(1.0-CLHEP::perMillion));
380
381 G4double stepAchieved = endAchieved
382 - SubStepStartState.GetCurveLength();
383
384 // Update remaining state - must work for 'full' step or
385 // abandonned intersection
386 //
387 CurrentState = IntersectPointVelct_G;
388 s_length_taken = stepAchieved;
389 if( shortEnd )
390 {
391 fParticleIsLooping = true;
392 }
393 }
394 }
395 }
396 if( !intersects )
397 {
398 StepTaken += s_length_taken;
399
400 if (fpTrajectoryFilter != nullptr) // For smooth trajectory display (jacek 1/11/2002)
401 {
402 fpTrajectoryFilter->TakeIntermediatePoint(CurrentState.GetPosition());
403 }
404 }
405 first_substep = false;
406
407#ifdef G4DEBUG_FIELD
408 if( fNoZeroStep > fActionThreshold_NoZeroSteps )
409 {
410 if( fNoZeroStep > fSevereActionThreshold_NoZeroSteps )
411 G4cout << " Above 'Severe Action' threshold -- for Zero steps. ";
412 else
413 G4cout << " Above 'action' threshold -- for Zero steps. ";
414 G4cout << " Number of zero steps = " << fNoZeroStep << G4endl;
415 printStatus( SubStepStartState, // or OriginalState,
416 CurrentState, CurrentProposedStepLength,
417 NewSafety, do_loop_count, pPhysVol );
418 }
419 if( (fVerboseLevel > 1) && (do_loop_count > fMax_loop_count-10 ))
420 {
421 if( do_loop_count == fMax_loop_count-9 )
422 {
423 G4cout << " G4PropagatorInField::ComputeStep(): " << G4endl
424 << " Difficult track - taking many sub steps." << G4endl;
425 printStatus( SubStepStartState, SubStepStartState, CurrentProposedStepLength,
426 NewSafety, 0, pPhysVol );
427 }
428 printStatus( SubStepStartState, CurrentState, CurrentProposedStepLength,
429 NewSafety, do_loop_count, pPhysVol );
430 }
431#endif
432
433 ++do_loop_count;
434
435 } while( (!intersects )
436 && (!fParticleIsLooping)
437 && (StepTaken + kCarTolerance < CurrentProposedStepLength)
438 && ( do_loop_count < fMax_loop_count ) );
439
440 if( do_loop_count >= fMax_loop_count
441 && (StepTaken + kCarTolerance < CurrentProposedStepLength) )
442 {
443 fParticleIsLooping = true;
444 }
445 if ( ( fParticleIsLooping ) && (fVerboseLevel > 0) )
446 {
447 ReportLoopingParticle( do_loop_count, StepTaken,
448 CurrentProposedStepLength, methodName,
449 CurrentState.GetMomentum(), pPhysVol );
450 }
451
452 if( !intersects )
453 {
454 // Chord AB or "minor chords" do not intersect
455 // B is the endpoint Step of the current Step.
456 //
457 End_PointAndTangent = CurrentState;
458 TruePathLength = StepTaken; // Original code
459
460 // Tried the following to avoid potential issue with round-off error
461 // - but has issues... Suppressing this change JA 2015/05/02
462 // TruePathLength = CurrentProposedStepLength;
463 }
464 fLastStepInVolume = intersects;
465
466 // Set pFieldTrack to the return value
467 //
468 pFieldTrack = End_PointAndTangent;
469
470#ifdef G4VERBOSE
471 // Check that "s" is correct
472 //
473 if( std::fabs(OriginalState.GetCurveLength() + TruePathLength
474 - End_PointAndTangent.GetCurveLength()) > 3.e-4 * TruePathLength )
475 {
476 std::ostringstream message;
477 message << "Curve length mis-match between original state "
478 << "and proposed endpoint of propagation." << G4endl
479 << " The curve length of the endpoint should be: "
480 << OriginalState.GetCurveLength() + TruePathLength << G4endl
481 << " and it is instead: "
482 << End_PointAndTangent.GetCurveLength() << "." << G4endl
483 << " A difference of: "
484 << OriginalState.GetCurveLength() + TruePathLength
485 - End_PointAndTangent.GetCurveLength() << G4endl
486 << " Original state = " << OriginalState << G4endl
487 << " Proposed state = " << End_PointAndTangent;
488 G4Exception(methodName, "GeomNav0003", FatalException, message);
489 }
490#endif
491
492 if( TruePathLength+kCarTolerance >= CurrentProposedStepLength )
493 {
494 fNoZeroStep = 0;
495 }
496 else
497 {
498 // In particular anomalous cases, we can get repeated zero steps
499 // We identify these cases and take corrective action when they occur.
500 //
501 if( TruePathLength < std::max( fZeroStepThreshold, 0.5*kCarTolerance ) )
502 {
503 ++fNoZeroStep;
504 }
505 else
506 {
507 fNoZeroStep = 0;
508 }
509 }
510 if( fNoZeroStep > fAbandonThreshold_NoZeroSteps )
511 {
512 fParticleIsLooping = true;
513 ReportStuckParticle( fNoZeroStep, CurrentProposedStepLength,
514 fFull_CurveLen_of_LastAttempt, pPhysVol );
515 fNoZeroStep = 0;
516 }
517
518 GetChordFinder()->SetDeltaChord(deltaChord);
519 return TruePathLength;
520}
521
522// ---------------------------------------------------------------------------
523// Dumps status of propagator
524//
525void
527 const G4FieldTrack& CurrentFT,
528 G4double requestStep,
529 G4double safety,
530 G4int stepNo,
531 G4VPhysicalVolume* startVolume)
532{
533 const G4int verboseLevel = fVerboseLevel;
534 const G4ThreeVector StartPosition = StartFT.GetPosition();
535 const G4ThreeVector StartUnitVelocity = StartFT.GetMomentumDir();
536 const G4ThreeVector CurrentPosition = CurrentFT.GetPosition();
537 const G4ThreeVector CurrentUnitVelocity = CurrentFT.GetMomentumDir();
538
539 G4double step_len = CurrentFT.GetCurveLength() - StartFT.GetCurveLength();
540
541 G4long oldprec; // cout/cerr precision settings
542
543 if( ((stepNo == 0) && (verboseLevel <3)) || (verboseLevel >= 3) )
544 {
545 oldprec = G4cout.precision(4);
546 G4cout << std::setw( 5) << "Step#"
547 << std::setw(10) << " s " << " "
548 << std::setw(10) << "X(mm)" << " "
549 << std::setw(10) << "Y(mm)" << " "
550 << std::setw(10) << "Z(mm)" << " "
551 << std::setw( 7) << " N_x " << " "
552 << std::setw( 7) << " N_y " << " "
553 << std::setw( 7) << " N_z " << " " ;
554 G4cout << std::setw( 7) << " Delta|N|" << " "
555 << std::setw( 9) << "StepLen" << " "
556 << std::setw(12) << "StartSafety" << " "
557 << std::setw( 9) << "PhsStep" << " ";
558 if( startVolume != nullptr )
559 { G4cout << std::setw(18) << "NextVolume" << " "; }
560 G4cout.precision(oldprec);
561 G4cout << G4endl;
562 }
563 if((stepNo == 0) && (verboseLevel <=3))
564 {
565 // Recurse to print the start values
566 //
567 printStatus( StartFT, StartFT, -1.0, safety, -1, startVolume);
568 }
569 if( verboseLevel <= 3 )
570 {
571 if( stepNo >= 0)
572 { G4cout << std::setw( 4) << stepNo << " "; }
573 else
574 { G4cout << std::setw( 5) << "Start" ; }
575 oldprec = G4cout.precision(8);
576 G4cout << std::setw(10) << CurrentFT.GetCurveLength() << " ";
577 G4cout.precision(8);
578 G4cout << std::setw(10) << CurrentPosition.x() << " "
579 << std::setw(10) << CurrentPosition.y() << " "
580 << std::setw(10) << CurrentPosition.z() << " ";
581 G4cout.precision(4);
582 G4cout << std::setw( 7) << CurrentUnitVelocity.x() << " "
583 << std::setw( 7) << CurrentUnitVelocity.y() << " "
584 << std::setw( 7) << CurrentUnitVelocity.z() << " ";
585 G4cout.precision(3);
586 G4cout << std::setw( 7)
587 << CurrentFT.GetMomentum().mag()-StartFT.GetMomentum().mag() << " ";
588 G4cout << std::setw( 9) << step_len << " ";
589 G4cout << std::setw(12) << safety << " ";
590 if( requestStep != -1.0 )
591 { G4cout << std::setw( 9) << requestStep << " "; }
592 else
593 { G4cout << std::setw( 9) << "Init/NotKnown" << " "; }
594 if( startVolume != nullptr)
595 { G4cout << std::setw(12) << startVolume->GetName() << " "; }
596 G4cout.precision(oldprec);
597 G4cout << G4endl;
598 }
599 else // if( verboseLevel > 3 )
600 {
601 // Multi-line output
602
603 G4cout << "Step taken was " << step_len
604 << " out of PhysicalStep = " << requestStep << G4endl;
605 G4cout << "Final safety is: " << safety << G4endl;
606 G4cout << "Chord length = " << (CurrentPosition-StartPosition).mag()
607 << G4endl;
608 G4cout << G4endl;
609 }
610}
611
612// ---------------------------------------------------------------------------
613// Prints Step diagnostics
614//
615void
617 G4double CurrentProposedStepLength,
618 G4double decreaseFactor,
619 G4double stepTrial,
620 const G4FieldTrack& )
621{
622 G4long iprec= G4cout.precision(8);
623 G4cout << " " << std::setw(12) << " PiF: NoZeroStep "
624 << " " << std::setw(20) << " CurrentProposed len "
625 << " " << std::setw(18) << " Full_curvelen_last"
626 << " " << std::setw(18) << " last proposed len "
627 << " " << std::setw(18) << " decrease factor "
628 << " " << std::setw(15) << " step trial "
629 << G4endl;
630
631 G4cout << " " << std::setw(10) << fNoZeroStep << " "
632 << " " << std::setw(20) << CurrentProposedStepLength
633 << " " << std::setw(18) << fFull_CurveLen_of_LastAttempt
634 << " " << std::setw(18) << fLast_ProposedStepLength
635 << " " << std::setw(18) << decreaseFactor
636 << " " << std::setw(15) << stepTrial
637 << G4endl;
638 G4cout.precision( iprec );
639}
640
641// Access the points which have passed through the filter. The
642// points are stored as ThreeVectors for the initial impelmentation
643// only (jacek 30/10/2002)
644// Responsibility for deleting the points lies with
645// SmoothTrajectoryPoint, which is the points' final
646// destination. The points pointer is set to NULL, to ensure that
647// the points are not re-used in subsequent steps, therefore THIS
648// METHOD MUST BE CALLED EXACTLY ONCE PER STEP. (jacek 08/11/2002)
649
650std::vector<G4ThreeVector>*
652{
653 // NB, GimmeThePointsAndForgetThem really forgets them, so it can
654 // only be called (exactly) once for each step.
655
656 if (fpTrajectoryFilter != nullptr)
657 {
658 return fpTrajectoryFilter->GimmeThePointsAndForgetThem();
659 }
660 return nullptr;
661}
662
663// ---------------------------------------------------------------------------
664//
665void
667{
668 fpTrajectoryFilter = filter;
669}
670
671// ---------------------------------------------------------------------------
672//
674{
675 // Goal: Clear all memory of previous steps, cached information
676
677 fParticleIsLooping = false;
678 fNoZeroStep = 0;
679
680 fSetFieldMgr = false; // Has field-manager been set for the current step?
681 fEpsilonStep= 1.0e-5; // Relative accuracy of current Step
682
683 End_PointAndTangent= G4FieldTrack( G4ThreeVector(0.,0.,0.),
684 G4ThreeVector(0.,0.,0.),
685 0.0,0.0,0.0,0.0,0.0);
686 fFull_CurveLen_of_LastAttempt = -1;
687 fLast_ProposedStepLength = -1;
688
689 fPreviousSftOrigin= G4ThreeVector(0.,0.,0.);
690 fPreviousSafety= 0.0;
691
692 fNewTrack = true;
693}
694
695// ---------------------------------------------------------------------------
696//
698FindAndSetFieldManager( G4VPhysicalVolume* pCurrentPhysicalVolume )
699{
700 G4FieldManager* currentFieldMgr;
701
702 currentFieldMgr = fDetectorFieldMgr;
703 if( pCurrentPhysicalVolume != nullptr )
704 {
705 G4FieldManager *pRegionFieldMgr = nullptr, *localFieldMgr = nullptr;
706 G4LogicalVolume* pLogicalVol = pCurrentPhysicalVolume->GetLogicalVolume();
707
708 if( pLogicalVol != nullptr )
709 {
710 // Value for Region, if any, overrides
711 //
712 G4Region* pRegion = pLogicalVol->GetRegion();
713 if( pRegion != nullptr )
714 {
715 pRegionFieldMgr = pRegion->GetFieldManager();
716 if( pRegionFieldMgr != nullptr )
717 {
718 currentFieldMgr= pRegionFieldMgr;
719 }
720 }
721
722 // 'Local' Value from logical volume, if any, overrides
723 //
724 localFieldMgr = pLogicalVol->GetFieldManager();
725 if ( localFieldMgr != nullptr )
726 {
727 currentFieldMgr = localFieldMgr;
728 }
729 }
730 }
731 fCurrentFieldMgr = currentFieldMgr;
732
733 // Flag that field manager has been set
734 //
735 fSetFieldMgr = true;
736
737 return currentFieldMgr;
738}
739
740// ---------------------------------------------------------------------------
741//
743{
744 G4int oldval = fVerboseLevel;
745 fVerboseLevel = level;
746
747 // Forward the verbose level 'reduced' to ChordFinder,
748 // MagIntegratorDriver ... ?
749 //
750 auto integrDriver = GetChordFinder()->GetIntegrationDriver();
751 integrDriver->SetVerboseLevel( fVerboseLevel - 2 );
752 G4cout << "Set Driver verbosity to " << fVerboseLevel - 2 << G4endl;
753
754 return oldval;
755}
756
757// ---------------------------------------------------------------------------
758//
760 G4double StepTaken,
761 G4double StepRequested,
762 const char* methodName,
763 const G4ThreeVector& momentumVec,
764 G4VPhysicalVolume* pPhysVol )
765{
766 std::ostringstream message;
767 G4double fraction = StepTaken / StepRequested;
768 message << " Unfinished integration of track (likely looping particle) "
769 << " of momentum " << momentumVec << " ( magnitude = "
770 << momentumVec.mag() << " ) " << G4endl
771 << " after " << count << " field substeps "
772 << " totaling " << std::setprecision(12) << StepTaken / mm << " mm "
773 << " out of requested step " << std::setprecision(12)
774 << StepRequested / mm << " mm ";
775 message << " a fraction of ";
776 G4int prec = 4;
777 if( fraction > 0.99 )
778 {
779 prec = 7;
780 }
781 else
782 {
783 if (fraction > 0.97 ) { prec = 5; }
784 }
785 message << std::setprecision(prec)
786 << 100. * StepTaken / StepRequested << " % " << G4endl ;
787 if( pPhysVol != nullptr )
788 {
789 message << " in volume " << pPhysVol->GetName() ;
790 auto material = pPhysVol->GetLogicalVolume()->GetMaterial();
791 if( material != nullptr )
792 {
793 message << " with material " << material->GetName()
794 << " ( density = "
795 << material->GetDensity() / ( g/(cm*cm*cm) ) << " g / cm^3 ) ";
796 }
797 }
798 else
799 {
800 message << " in unknown (null) volume. " ;
801 }
802 G4Exception(methodName, "GeomNav1002", JustWarning, message);
803}
804
805// ---------------------------------------------------------------------------
806//
808 G4double proposedStep,
809 G4double lastTriedStep,
810 G4VPhysicalVolume* physVol )
811{
812 std::ostringstream message;
813 message << "Particle is stuck; it will be killed." << G4endl
814 << " Zero progress for " << noZeroSteps << " attempted steps."
815 << G4endl
816 << " Proposed Step is " << proposedStep
817 << " but Step Taken is "<< lastTriedStep << G4endl;
818 if( physVol != nullptr )
819 {
820 message << " in volume " << physVol->GetName() ;
821 }
822 else
823 {
824 message << " in unknown or null volume. " ;
825 }
826 G4Exception("G4PropagatorInField::ComputeStep()",
827 "GeomNav1002", JustWarning, message);
828}
829
830// ------------------------------------------------------------------------
831
832// ----------------------------------------------
833// Methods to alter Parameters
834// ----------------------------------------------
835
836// Was a data member (of an object) -- now moved to class member
838{
839 return fLargestAcceptableStep;
840}
841
842// ------------------------------------------------------------------------
843//
845{
846 if( fLargestAcceptableStep>0.0 )
847 {
848 fLargestAcceptableStep = newBigDist;
849 }
850}
851
852// ---------------------------------------------------------------------------
853
855{
856 return fMaxStepSizeMultiplier;
857}
858
859// ---------------------------------------------------------------------------
860
862{
863 fMaxStepSizeMultiplier=vm;
864}
865
866// ---------------------------------------------------------------------------
867
869{
870 return fMinBigDistance;
871}
872
873// ---------------------------------------------------------------------------
874
876{
877 fMinBigDistance= val;
878}
G4double epsilon(G4double density, G4double temperature)
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
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:67
G4GLOB_DLL std::ostream G4cout
double z() const
double x() const
double y() const
double mag() const
G4double GetDeltaChord() const
G4VIntegrationDriver * GetIntegrationDriver()
void OnComputeStep(const G4FieldTrack *track)
void SetDeltaChord(G4double newval)
G4double AdvanceChordLimited(G4FieldTrack &yCurrent, G4double stepInitial, G4double epsStep_Relative, const G4ThreeVector &latestSafetyOrigin, G4double lasestSafetyRadius)
G4FieldManager is a manager (store) for a pointer to the Field subclass that describes the field of a...
G4FieldTrack defines a data structure bringing together a magnetic track's state (position,...
const G4ThreeVector & GetMomentumDir() const
G4double GetCurveLength() const
G4ThreeVector GetPosition() const
G4ThreeVector GetMomentum() const
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()
G4LogicalVolume represents a leaf node or unpositioned subtree in the geometry hierarchy....
G4Region * GetRegion() const
G4Material * GetMaterial() const
G4FieldManager * GetFieldManager() const
G4MultiLevelLocator implements the calculation of the intersection point with a boundary when G4Propa...
G4Navigator is a class for use by the tracking management, able to obtain/calculate dynamic tracking ...
void printStatus(const G4FieldTrack &startFT, const G4FieldTrack &currentFT, G4double requestStep, G4double safety, G4int step, G4VPhysicalVolume *startVolume)
G4FieldManager * FindAndSetFieldManager(G4VPhysicalVolume *pCurrentPhysVol)
G4int SetVerboseLevel(G4int verbose)
void PrintStepLengthDiagnostic(G4double currentProposedStepLength, G4double decreaseFactor, G4double stepTrial, const G4FieldTrack &aFieldTrack)
G4PropagatorInField(G4Navigator *theNavigator, G4FieldManager *detectorFieldMgr, G4VIntersectionLocator *vLocator=nullptr)
void ReportStuckParticle(G4int noZeroSteps, G4double proposedStep, G4double lastTriedStep, G4VPhysicalVolume *physVol)
void SetMinBigDistance(G4double val)
G4ChordFinder * GetChordFinder()
void SetMaxStepSizeMultiplier(G4double vm)
G4bool IntersectChord(const G4ThreeVector &StartPointA, const G4ThreeVector &EndPointB, G4double &NewSafety, G4double &LinearStepLength, G4ThreeVector &IntersectionPoint)
std::vector< G4ThreeVector > * GimmeTrajectoryVectorAndForgetIt() const
void ReportLoopingParticle(G4int count, G4double StepTaken, G4double stepRequest, const char *methodName, const G4ThreeVector &momentumVec, G4VPhysicalVolume *physVol)
void SetTrajectoryFilter(G4VCurvedTrajectoryFilter *filter)
G4double ComputeStep(G4FieldTrack &pFieldTrack, G4double pCurrentProposedStepLength, G4double &pNewSafety, G4VPhysicalVolume *pPhysVol=nullptr, G4bool canRelaxDeltaChord=false)
void SetLargestAcceptableStep(G4double newBigDist)
void SetEpsilonStep(G4double newEps)
G4Region defines a region or a group of regions in the detector geometry setup, sharing properties as...
Definition G4Region.hh:90
G4FieldManager * GetFieldManager() const
G4VCurvedTrajectoryFilter defines a filter for deciding which intermediate points on a curved traject...
virtual void SetVerboseLevel(G4int level)=0
G4VIntersectionLocator is a base class for the calculation of the intersection point with a boundary ...
G4VPhysicalVolume is an abstract base class for the representation of a positioned volume....
G4LogicalVolume * GetLogicalVolume() const
const G4String & GetName() const