Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4PropagatorInField Class Reference

#include <G4PropagatorInField.hh>

Public Member Functions

 G4PropagatorInField (G4Navigator *theNavigator, G4FieldManager *detectorFieldMgr, G4VIntersectionLocator *vLocator=0)
 
 ~G4PropagatorInField ()
 
G4double ComputeStep (G4FieldTrack &pFieldTrack, G4double pCurrentProposedStepLength, G4double &pNewSafety, G4VPhysicalVolume *pPhysVol=0)
 
G4ThreeVector EndPosition () const
 
G4ThreeVector EndMomentumDir () const
 
G4bool IsParticleLooping () const
 
G4double GetEpsilonStep () const
 
void SetEpsilonStep (G4double newEps)
 
G4FieldManagerFindAndSetFieldManager (G4VPhysicalVolume *pCurrentPhysVol)
 
G4ChordFinderGetChordFinder ()
 
G4int SetVerboseLevel (G4int verbose)
 
G4int GetVerboseLevel () const
 
G4int Verbose () const
 
void SetVerboseTrace (G4bool enable)
 
G4bool GetVerboseTrace ()
 
G4int GetMaxLoopCount () const
 
void SetMaxLoopCount (G4int new_max)
 
void printStatus (const G4FieldTrack &startFT, const G4FieldTrack &currentFT, G4double requestStep, G4double safety, G4int step, G4VPhysicalVolume *startVolume)
 
G4FieldTrack GetEndState () const
 
G4double GetMinimumEpsilonStep () const
 
void SetMinimumEpsilonStep (G4double newEpsMin)
 
G4double GetMaximumEpsilonStep () const
 
void SetMaximumEpsilonStep (G4double newEpsMax)
 
void SetLargestAcceptableStep (G4double newBigDist)
 
G4double GetLargestAcceptableStep ()
 
void SetTrajectoryFilter (G4VCurvedTrajectoryFilter *filter)
 
std::vector< G4ThreeVector > * GimmeTrajectoryVectorAndForgetIt () const
 
void ClearPropagatorState ()
 
void SetDetectorFieldManager (G4FieldManager *newGlobalFieldManager)
 
void SetUseSafetyForOptimization (G4bool)
 
G4bool GetUseSafetyForOptimization ()
 
G4bool IntersectChord (const G4ThreeVector &StartPointA, const G4ThreeVector &EndPointB, G4double &NewSafety, G4double &LinearStepLength, G4ThreeVector &IntersectionPoint)
 
G4bool IsFirstStepInVolume ()
 
G4bool IsLastStepInVolume ()
 
void PrepareNewTrack ()
 
G4VIntersectionLocatorGetIntersectionLocator ()
 
void SetIntersectionLocator (G4VIntersectionLocator *pLocator)
 
G4double GetDeltaIntersection () const
 
G4double GetDeltaOneStep () const
 
G4FieldManagerGetCurrentFieldManager ()
 
G4EquationOfMotionGetCurrentEquationOfMotion ()
 
void SetNavigatorForPropagating (G4Navigator *SimpleOrMultiNavigator)
 
G4NavigatorGetNavigatorForPropagating ()
 
void SetThresholdNoZeroStep (G4int noAct, G4int noHarsh, G4int noAbandon)
 
G4int GetThresholdNoZeroSteps (G4int i)
 
G4double GetZeroStepThreshold ()
 
void SetZeroStepThreshold (G4double newLength)
 
void RefreshIntersectionLocator ()
 

Protected Member Functions

void PrintStepLengthDiagnostic (G4double currentProposedStepLength, G4double decreaseFactor, G4double stepTrial, const G4FieldTrack &aFieldTrack)
 
void ReportLoopingParticle (G4int count, double StepTaken, G4VPhysicalVolume *pPhysVol)
 
void ReportStuckParticle (G4int noZeroSteps, G4double proposedStep, G4double lastTriedStep, G4VPhysicalVolume *physVol)
 

Detailed Description

Definition at line 64 of file G4PropagatorInField.hh.

Constructor & Destructor Documentation

G4PropagatorInField::G4PropagatorInField ( G4Navigator theNavigator,
G4FieldManager detectorFieldMgr,
G4VIntersectionLocator vLocator = 0 
)

Definition at line 60 of file G4PropagatorInField.cc.

63  :
64  fMax_loop_count(1000),
65  fUseSafetyForOptimisation(true), // (false) is less sensitive to incorrect safety
66  fZeroStepThreshold( 0.0 ), // length of what is recognised as 'zero' step
67  fDetectorFieldMgr(detectorFieldMgr),
68  fpTrajectoryFilter( 0 ),
69  fNavigator(theNavigator),
70  fCurrentFieldMgr(detectorFieldMgr),
71  fSetFieldMgr(false),
72  End_PointAndTangent(G4ThreeVector(0.,0.,0.),
73  G4ThreeVector(0.,0.,0.),0.0,0.0,0.0,0.0,0.0),
74  fParticleIsLooping(false),
75  fNoZeroStep(0),
76  fVerboseLevel(0),
77  fVerbTracePiF(false),
78  fFirstStepInVolume(true),
79  fLastStepInVolume(true),
80  fNewTrack(true)
81 {
82  if(fDetectorFieldMgr) { fEpsilonStep = fDetectorFieldMgr->GetMaximumEpsilonStep();}
83  else { fEpsilonStep= 1.0e-5; }
84  fActionThreshold_NoZeroSteps = 2;
85  fSevereActionThreshold_NoZeroSteps = 10;
86  fAbandonThreshold_NoZeroSteps = 50;
87  fFull_CurveLen_of_LastAttempt = -1;
88  fLast_ProposedStepLength = -1;
89  fLargestAcceptableStep = 1000.0 * meter;
90 
91  fPreviousSftOrigin= G4ThreeVector(0.,0.,0.);
92  fPreviousSafety= 0.0;
94  fZeroStepThreshold= std::max( 1.0e5 * kCarTolerance, 1.0e-1 * micrometer );
95 
96 #ifdef G4DEBUG_FIELD
97  G4cout << " PiF: Zero Step Threshold set to "
98  << fZeroStepThreshold / millimeter
99  << " mm." << G4endl;
100  G4cout << " PiF: Value of kCarTolerance = "
101  << kCarTolerance / millimeter
102  << " mm. " << G4endl;
103  fVerboseLevel = 3;
104  fVerbTracePiF = true;
105 #endif
106 
107  // Defining Intersection Locator and his parameters
108  if (vLocator==0)
109  {
110  fIntersectionLocator= new G4MultiLevelLocator(theNavigator);
111  fAllocatedLocator= true;
112  }
113  else
114  {
115  fIntersectionLocator= vLocator;
116  fAllocatedLocator= false;
117  }
118  RefreshIntersectionLocator(); // Copy all relevant parameters
119 }
CLHEP::Hep3Vector G4ThreeVector
G4double GetSurfaceTolerance() const
static constexpr double meter
Definition: G4SIunits.hh:82
G4GLOB_DLL std::ostream G4cout
G4double GetMaximumEpsilonStep() const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
static constexpr double millimeter
Definition: G4SIunits.hh:86
#define G4endl
Definition: G4ios.hh:61
static constexpr double micrometer
Definition: G4SIunits.hh:100
static G4GeometryTolerance * GetInstance()

Here is the call graph for this function:

G4PropagatorInField::~G4PropagatorInField ( )

Definition at line 123 of file G4PropagatorInField.cc.

124 {
125  if(fAllocatedLocator) { delete fIntersectionLocator; }
126 }

Member Function Documentation

void G4PropagatorInField::ClearPropagatorState ( )

Definition at line 669 of file G4PropagatorInField.cc.

670 {
671  // Goal: Clear all memory of previous steps, cached information
672 
673  fParticleIsLooping= false;
674  fNoZeroStep= 0;
675 
676  End_PointAndTangent= G4FieldTrack( G4ThreeVector(0.,0.,0.),
677  G4ThreeVector(0.,0.,0.),
678  0.0,0.0,0.0,0.0,0.0);
679  fFull_CurveLen_of_LastAttempt = -1;
680  fLast_ProposedStepLength = -1;
681 
682  fPreviousSftOrigin= G4ThreeVector(0.,0.,0.);
683  fPreviousSafety= 0.0;
684 }
CLHEP::Hep3Vector G4ThreeVector

Here is the caller graph for this function:

G4double G4PropagatorInField::ComputeStep ( G4FieldTrack pFieldTrack,
G4double  pCurrentProposedStepLength,
G4double pNewSafety,
G4VPhysicalVolume pPhysVol = 0 
)

Definition at line 145 of file G4PropagatorInField.cc.

150 {
151  // If CurrentProposedStepLength is too small for finding Chords
152  // then return with no action (for now - TODO: some action)
153  //
154  if(CurrentProposedStepLength<kCarTolerance)
155  {
156  return kInfinity;
157  }
158 
159  // Introducing smooth trajectory display (jacek 01/11/2002)
160  //
161  if (fpTrajectoryFilter)
162  {
163  fpTrajectoryFilter->CreateNewTrajectorySegment();
164  }
165 
166  fFirstStepInVolume = fNewTrack ? true : fLastStepInVolume;
167  fLastStepInVolume= false;
168  fNewTrack= false;
169 
170  if( fVerboseLevel > 2 )
171  {
172  G4cout << "G4PropagatorInField::ComputeStep() called" << G4endl;
173  G4cout << " Starting FT: " << pFieldTrack;
174  G4cout << " Requested length = " << CurrentProposedStepLength << G4endl;
175  G4cout << " PhysVol = ";
176  if( pPhysVol )
177  G4cout << pPhysVol->GetName() << G4endl;
178  else
179  G4cout << " N/A ";
180  G4cout << G4endl;
181  }
182 
183  // Parameters for adaptive Runge-Kutta integration
184 
185  G4double h_TrialStepSize; // 1st Step Size
186  G4double TruePathLength = CurrentProposedStepLength;
187  G4double StepTaken = 0.0;
188  G4double s_length_taken, epsilon ;
189  G4bool intersects;
190  G4bool first_substep = true;
191 
192  G4double NewSafety;
193  fParticleIsLooping = false;
194 
195  // If not yet done,
196  // Set the field manager to the local one if the volume has one,
197  // or to the global one if not
198  //
199  if( !fSetFieldMgr ) fCurrentFieldMgr= FindAndSetFieldManager( pPhysVol );
200  // For the next call, the field manager must again be set
201  fSetFieldMgr= false;
202 
203  G4FieldTrack CurrentState(pFieldTrack);
204  G4FieldTrack OriginalState = CurrentState;
205 
206  // If the Step length is "infinite", then an approximate-maximum Step
207  // length (used to calculate the relative accuracy) must be guessed.
208  //
209  if( CurrentProposedStepLength >= fLargestAcceptableStep )
210  {
211  G4ThreeVector StartPointA, VelocityUnit;
212  StartPointA = pFieldTrack.GetPosition();
213  VelocityUnit = pFieldTrack.GetMomentumDir();
214 
215  G4double trialProposedStep = 1.e2 * ( 10.0 * cm +
216  fNavigator->GetWorldVolume()->GetLogicalVolume()->
217  GetSolid()->DistanceToOut(StartPointA, VelocityUnit) );
218  CurrentProposedStepLength= std::min( trialProposedStep,
219  fLargestAcceptableStep );
220  }
221  epsilon = fCurrentFieldMgr->GetDeltaOneStep() / CurrentProposedStepLength;
222  // G4double raw_epsilon= epsilon;
223  G4double epsilonMin= fCurrentFieldMgr->GetMinimumEpsilonStep();
224  G4double epsilonMax= fCurrentFieldMgr->GetMaximumEpsilonStep();
225  if( epsilon < epsilonMin ) epsilon = epsilonMin;
226  if( epsilon > epsilonMax ) epsilon = epsilonMax;
227  SetEpsilonStep( epsilon );
228 
229 
230  // Values for Intersection Locator has to be updated on each call for the
231  // case that CurrentFieldManager has changed from the one of previous step
233 
234  // G4cout << "G4PiF: Epsilon of current step - raw= " << raw_epsilon
235  // << " final= " << epsilon << G4endl;
236 
237  // Shorten the proposed step in case of earlier problems (zero steps)
238  //
239  if( fNoZeroStep > fActionThreshold_NoZeroSteps )
240  {
241  G4double stepTrial;
242 
243  stepTrial= fFull_CurveLen_of_LastAttempt;
244  if( (stepTrial <= 0.0) && (fLast_ProposedStepLength > 0.0) )
245  stepTrial= fLast_ProposedStepLength;
246 
247  G4double decreaseFactor = 0.9; // Unused default
248  if( (fNoZeroStep < fSevereActionThreshold_NoZeroSteps)
249  && (stepTrial > 100.0*fZeroStepThreshold) )
250  {
251  // Attempt quick convergence
252  //
253  decreaseFactor= 0.25;
254  }
255  else
256  {
257  // We are in significant difficulties, probably at a boundary that
258  // is either geometrically sharp or between very different materials.
259  // Careful decreases to cope with tolerance are required.
260  //
261  if( stepTrial > 100.0*fZeroStepThreshold )
262  decreaseFactor = 0.35; // Try decreasing slower
263  else if( stepTrial > 30.0*fZeroStepThreshold )
264  decreaseFactor= 0.5; // Try yet slower decrease
265  else if( stepTrial > 10.0*fZeroStepThreshold )
266  decreaseFactor= 0.75; // Try even slower decreases
267  else
268  decreaseFactor= 0.9; // Try very slow decreases
269  }
270  stepTrial *= decreaseFactor;
271 
272 #ifdef G4DEBUG_FIELD
273  G4cerr << " G4PropagatorInField::ComputeStep(): " << G4endl
274  << " Decreasing step - in volume " << pPhysVol;
275  if( pPhysVol )
276  G4cerr << " with name " << pPhysVol->GetName();
277  G4cerr << G4endl;
278  PrintStepLengthDiagnostic(CurrentProposedStepLength, decreaseFactor,
279  stepTrial, pFieldTrack);
280 #endif
281  if( stepTrial == 0.0 ) // Change to make it < 0.1 * kCarTolerance ??
282  {
283  std::ostringstream message;
284  message << "Particle abandoned due to lack of progress in field."
285  << G4endl
286  << " Properties : " << pFieldTrack << G4endl
287  << " Attempting a zero step = " << stepTrial << G4endl
288  << " while attempting to progress after " << fNoZeroStep
289  << " trial steps. Will abandon step.";
290  G4Exception("G4PropagatorInField::ComputeStep()", "GeomNav1002",
291  JustWarning, message);
292  fParticleIsLooping= true;
293  return 0; // = stepTrial;
294  }
295  if( stepTrial < CurrentProposedStepLength )
296  CurrentProposedStepLength = stepTrial;
297  }
298  fLast_ProposedStepLength = CurrentProposedStepLength;
299 
300  G4int do_loop_count = 0;
301  do // Loop checking, 07.10.2016, J.Apostolakis
302  {
303  G4FieldTrack SubStepStartState = CurrentState;
304  G4ThreeVector SubStartPoint = CurrentState.GetPosition();
305 
306  if(!first_substep)
307  {
308  if( fVerboseLevel > 4 )
309  {
310  G4cout << " PiF: Calling Nav/Locate Global Point within-Volume "
311  << G4endl;
312  }
313  fNavigator->LocateGlobalPointWithinVolume( SubStartPoint );
314  }
315 
316  // How far to attempt to move the particle !
317  //
318  h_TrialStepSize = CurrentProposedStepLength - StepTaken;
319 
320  // Integrate as far as "chord miss" rule allows.
321  //
322  s_length_taken = GetChordFinder()->AdvanceChordLimited(
323  CurrentState, // Position & velocity
324  h_TrialStepSize,
325  fEpsilonStep,
326  fPreviousSftOrigin,
327  fPreviousSafety
328  );
329  // CurrentState is now updated with the final position and velocity.
330 
331  fFull_CurveLen_of_LastAttempt = s_length_taken;
332 
333  G4ThreeVector EndPointB = CurrentState.GetPosition();
334  G4ThreeVector InterSectionPointE;
335  G4double LinearStepLength;
336 
337  // Intersect chord AB with geometry
338  intersects= IntersectChord( SubStartPoint, EndPointB,
339  NewSafety, LinearStepLength,
340  InterSectionPointE );
341  // E <- Intersection Point of chord AB and either volume A's surface
342  // or a daughter volume's surface ..
343 
344  if( first_substep )
345  {
346  currentSafety = NewSafety;
347  } // Updating safety in other steps is potential future extention
348 
349  if( intersects )
350  {
351  G4FieldTrack IntersectPointVelct_G(CurrentState); // FT-Def-Construct
352 
353  // Find the intersection point of AB true path with the surface
354  // of vol(A), if it exists. Start with point E as first "estimate".
355  G4bool recalculatedEndPt= false;
356 
357  G4bool found_intersection = fIntersectionLocator->
358  EstimateIntersectionPoint( SubStepStartState, CurrentState,
359  InterSectionPointE, IntersectPointVelct_G,
360  recalculatedEndPt, fPreviousSafety,
361  fPreviousSftOrigin);
362  intersects = found_intersection;
363  if( found_intersection )
364  {
365  End_PointAndTangent= IntersectPointVelct_G; // G is our EndPoint ...
366  StepTaken = TruePathLength = IntersectPointVelct_G.GetCurveLength()
367  - OriginalState.GetCurveLength();
368  }
369  else
370  {
371  // Either "minor" chords do not intersect
372  // or else stopped (due to too many steps)
373  //
374  if( recalculatedEndPt )
375  {
376  G4double endAchieved = IntersectPointVelct_G.GetCurveLength();
377  G4double endExpected = CurrentState.GetCurveLength();
378 
379  // Detect failure - due to too many steps
380  G4bool shortEnd = endAchieved
381  < (endExpected*(1.0-CLHEP::perMillion));
382 
383  G4double stepAchieved = endAchieved
384  - SubStepStartState.GetCurveLength();
385 
386  // Update remaining state - must work for 'full' step or
387  // abandonned intersection
388  //
389  CurrentState= IntersectPointVelct_G;
390  s_length_taken = stepAchieved;
391  if( shortEnd )
392  {
393  fParticleIsLooping = true;
394  }
395  }
396  }
397  }
398  if( !intersects )
399  {
400  StepTaken += s_length_taken;
401  // For smooth trajectory display (jacek 01/11/2002)
402  if (fpTrajectoryFilter) {
403  fpTrajectoryFilter->TakeIntermediatePoint(CurrentState.GetPosition());
404  }
405  }
406  first_substep = false;
407 
408 #ifdef G4DEBUG_FIELD
409  if( fNoZeroStep > fActionThreshold_NoZeroSteps )
410  {
411  printStatus( SubStepStartState, // or OriginalState,
412  CurrentState, CurrentProposedStepLength,
413  NewSafety, do_loop_count, pPhysVol );
414  }
415  if( (fVerboseLevel > 1) && (do_loop_count > fMax_loop_count-10 ))
416  {
417  if( do_loop_count == fMax_loop_count-9 )
418  {
419  G4cout << " G4PropagatorInField::ComputeStep(): " << G4endl
420  << " Difficult track - taking many sub steps." << G4endl;
421  }
422  printStatus( SubStepStartState, CurrentState, CurrentProposedStepLength,
423  NewSafety, do_loop_count, pPhysVol );
424  }
425 #endif
426 
427  do_loop_count++;
428 
429  } while( (!intersects )
430  && (!fParticleIsLooping)
431  && (StepTaken + kCarTolerance < CurrentProposedStepLength)
432  && ( do_loop_count < fMax_loop_count ) );
433 
434  if( do_loop_count >= fMax_loop_count )
435  {
436  fParticleIsLooping = true;
437  }
438  if ( fParticleIsLooping && (fVerboseLevel > 0) )
439  {
440  ReportLoopingParticle( do_loop_count, StepTaken, pPhysVol );
441  }
442 
443  if( !intersects )
444  {
445  // Chord AB or "minor chords" do not intersect
446  // B is the endpoint Step of the current Step.
447  //
448  End_PointAndTangent = CurrentState;
449  TruePathLength = StepTaken; // Original code
450  // Tried the following to avoid potential issue with round-off error
451  // - but has issues... Suppressing this change JA 2015/05/02
452  // TruePathLength = CurrentProposedStepLength;
453  }
454  fLastStepInVolume = intersects;
455 
456  // Set pFieldTrack to the return value
457  //
458  pFieldTrack = End_PointAndTangent;
459 
460 #ifdef G4VERBOSE
461  // Check that "s" is correct
462  //
463  if( std::fabs(OriginalState.GetCurveLength() + TruePathLength
464  - End_PointAndTangent.GetCurveLength()) > 3.e-4 * TruePathLength )
465  {
466  std::ostringstream message;
467  message << "Curve length mis-match between original state "
468  << "and proposed endpoint of propagation." << G4endl
469  << " The curve length of the endpoint should be: "
470  << OriginalState.GetCurveLength() + TruePathLength << G4endl
471  << " and it is instead: "
472  << End_PointAndTangent.GetCurveLength() << "." << G4endl
473  << " A difference of: "
474  << OriginalState.GetCurveLength() + TruePathLength
475  - End_PointAndTangent.GetCurveLength() << G4endl
476  << " Original state = " << OriginalState << G4endl
477  << " Proposed state = " << End_PointAndTangent;
478  G4Exception("G4PropagatorInField::ComputeStep()",
479  "GeomNav0003", FatalException, message);
480  }
481 #endif
482 
483  if( TruePathLength+kCarTolerance >= CurrentProposedStepLength )
484  {
485  fNoZeroStep = 0;
486  }
487  else
488  {
489  // In particular anomalous cases, we can get repeated zero steps
490  // We identify these cases and take corrective action when they occur.
491  //
492  if( TruePathLength < std::max( fZeroStepThreshold, 0.5*kCarTolerance ) )
493  {
494  fNoZeroStep++;
495  }
496  else{
497  fNoZeroStep = 0;
498  }
499  }
500  if( fNoZeroStep > fAbandonThreshold_NoZeroSteps )
501  {
502  fParticleIsLooping = true;
503  ReportStuckParticle( fNoZeroStep, CurrentProposedStepLength, fFull_CurveLen_of_LastAttempt,
504  pPhysVol );
505  fNoZeroStep = 0;
506  }
507 
508  return TruePathLength;
509 }
void SetEpsilonStep(G4double newEps)
static const G4double kInfinity
Definition: geomdefs.hh:42
G4double GetCurveLength() const
virtual void TakeIntermediatePoint(G4ThreeVector newPoint)=0
void PrintStepLengthDiagnostic(G4double currentProposedStepLength, G4double decreaseFactor, G4double stepTrial, const G4FieldTrack &aFieldTrack)
G4double GetDeltaOneStep() const
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
const G4String & GetName() const
void ReportLoopingParticle(G4int count, double StepTaken, G4VPhysicalVolume *pPhysVol)
bool G4bool
Definition: G4Types.hh:79
G4double GetMaximumEpsilonStep() const
static constexpr double cm
Definition: G4SIunits.hh:119
static constexpr double perMillion
G4FieldManager * FindAndSetFieldManager(G4VPhysicalVolume *pCurrentPhysVol)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4double AdvanceChordLimited(G4FieldTrack &yCurrent, G4double stepInitial, G4double epsStep_Relative, const G4ThreeVector latestSafetyOrigin, G4double lasestSafetyRadius)
G4LogicalVolume * GetLogicalVolume() const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
void ReportStuckParticle(G4int noZeroSteps, G4double proposedStep, G4double lastTriedStep, G4VPhysicalVolume *physVol)
G4ChordFinder * GetChordFinder()
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4bool IntersectChord(const G4ThreeVector &StartPointA, const G4ThreeVector &EndPointB, G4double &NewSafety, G4double &LinearStepLength, G4ThreeVector &IntersectionPoint)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4VPhysicalVolume * GetWorldVolume() const
double epsilon(double density, double temperature)
virtual void LocateGlobalPointWithinVolume(const G4ThreeVector &position)
Definition: G4Navigator.cc:582
G4GLOB_DLL std::ostream G4cerr
void printStatus(const G4FieldTrack &startFT, const G4FieldTrack &currentFT, G4double requestStep, G4double safety, G4int step, G4VPhysicalVolume *startVolume)
G4double GetMinimumEpsilonStep() const

Here is the call graph for this function:

Here is the caller graph for this function:

G4ThreeVector G4PropagatorInField::EndMomentumDir ( ) const
inline
G4ThreeVector G4PropagatorInField::EndPosition ( ) const
inline
G4FieldManager * G4PropagatorInField::FindAndSetFieldManager ( G4VPhysicalVolume pCurrentPhysVol)

Definition at line 687 of file G4PropagatorInField.cc.

688 {
689  G4FieldManager* currentFieldMgr;
690 
691  currentFieldMgr = fDetectorFieldMgr;
692  if( pCurrentPhysicalVolume)
693  {
694  G4FieldManager *pRegionFieldMgr= 0, *localFieldMgr = 0;
695  G4LogicalVolume* pLogicalVol= pCurrentPhysicalVolume->GetLogicalVolume();
696 
697  if( pLogicalVol ) {
698  // Value for Region, if any, Overrides
699  G4Region* pRegion= pLogicalVol->GetRegion();
700  if( pRegion ) {
701  pRegionFieldMgr= pRegion->GetFieldManager();
702  if( pRegionFieldMgr )
703  currentFieldMgr= pRegionFieldMgr;
704  }
705 
706  // 'Local' Value from logical volume, if any, Overrides
707  localFieldMgr= pLogicalVol->GetFieldManager();
708  if ( localFieldMgr )
709  currentFieldMgr = localFieldMgr;
710  }
711  }
712  fCurrentFieldMgr= currentFieldMgr;
713 
714  // Flag that field manager has been set
715  //
716  fSetFieldMgr= true;
717 
718  return currentFieldMgr;
719 }
G4Region * GetRegion() const
G4FieldManager * GetFieldManager() const
G4FieldManager * GetFieldManager() const

Here is the call graph for this function:

Here is the caller graph for this function:

G4ChordFinder* G4PropagatorInField::GetChordFinder ( )
inline

Here is the caller graph for this function:

G4EquationOfMotion* G4PropagatorInField::GetCurrentEquationOfMotion ( )
inline
G4FieldManager* G4PropagatorInField::GetCurrentFieldManager ( )
inline

Here is the caller graph for this function:

G4double G4PropagatorInField::GetDeltaIntersection ( ) const
inline
G4double G4PropagatorInField::GetDeltaOneStep ( ) const
inline
G4FieldTrack G4PropagatorInField::GetEndState ( ) const
inline
G4double G4PropagatorInField::GetEpsilonStep ( ) const
inline
G4VIntersectionLocator* G4PropagatorInField::GetIntersectionLocator ( )
inline
G4double G4PropagatorInField::GetLargestAcceptableStep ( )
inline

Here is the caller graph for this function:

G4double G4PropagatorInField::GetMaximumEpsilonStep ( ) const
inline
G4int G4PropagatorInField::GetMaxLoopCount ( ) const
inline
G4double G4PropagatorInField::GetMinimumEpsilonStep ( ) const
inline
G4Navigator* G4PropagatorInField::GetNavigatorForPropagating ( )
inline
G4int G4PropagatorInField::GetThresholdNoZeroSteps ( G4int  i)
inline
G4bool G4PropagatorInField::GetUseSafetyForOptimization ( )
inline
G4int G4PropagatorInField::GetVerboseLevel ( ) const
inline
G4bool G4PropagatorInField::GetVerboseTrace ( )
inline

Definition at line 103 of file G4PropagatorInField.hh.

103 { return fVerbTracePiF; }
G4double G4PropagatorInField::GetZeroStepThreshold ( )
inline
std::vector< G4ThreeVector > * G4PropagatorInField::GimmeTrajectoryVectorAndForgetIt ( ) const

Definition at line 646 of file G4PropagatorInField.cc.

647 {
648  // NB, GimmeThePointsAndForgetThem really forgets them, so it can
649  // only be called (exactly) once for each step.
650 
651  if (fpTrajectoryFilter)
652  {
653  return fpTrajectoryFilter->GimmeThePointsAndForgetThem();
654  }
655  else
656  {
657  return 0;
658  }
659 }
std::vector< G4ThreeVector > * GimmeThePointsAndForgetThem()

Here is the call graph for this function:

Here is the caller graph for this function:

G4bool G4PropagatorInField::IntersectChord ( const G4ThreeVector StartPointA,
const G4ThreeVector EndPointB,
G4double NewSafety,
G4double LinearStepLength,
G4ThreeVector IntersectionPoint 
)
inline

Here is the caller graph for this function:

G4bool G4PropagatorInField::IsFirstStepInVolume ( )
inline

Definition at line 162 of file G4PropagatorInField.hh.

162 { return fFirstStepInVolume; }
G4bool G4PropagatorInField::IsLastStepInVolume ( )
inline

Definition at line 163 of file G4PropagatorInField.hh.

163 { return fLastStepInVolume; }

Here is the caller graph for this function:

G4bool G4PropagatorInField::IsParticleLooping ( ) const
inline

Here is the caller graph for this function:

void G4PropagatorInField::PrepareNewTrack ( )
inline

Definition at line 164 of file G4PropagatorInField.hh.

164 { fNewTrack = true; fFirstStepInVolume=false; fLastStepInVolume=false; }

Here is the caller graph for this function:

void G4PropagatorInField::printStatus ( const G4FieldTrack startFT,
const G4FieldTrack currentFT,
G4double  requestStep,
G4double  safety,
G4int  step,
G4VPhysicalVolume startVolume 
)

Definition at line 516 of file G4PropagatorInField.cc.

522 {
523  const G4int verboseLevel=fVerboseLevel;
524  const G4ThreeVector StartPosition = StartFT.GetPosition();
525  const G4ThreeVector StartUnitVelocity = StartFT.GetMomentumDir();
526  const G4ThreeVector CurrentPosition = CurrentFT.GetPosition();
527  const G4ThreeVector CurrentUnitVelocity = CurrentFT.GetMomentumDir();
528 
529  G4double step_len = CurrentFT.GetCurveLength() - StartFT.GetCurveLength();
530 
531  G4int oldprec; // cout/cerr precision settings
532 
533  if( ((stepNo == 0) && (verboseLevel <3)) || (verboseLevel >= 3) )
534  {
535  oldprec = G4cout.precision(4);
536  G4cout << std::setw( 6) << " "
537  << std::setw( 25) << " Current Position and Direction" << " "
538  << G4endl;
539  G4cout << std::setw( 5) << "Step#"
540  << std::setw(10) << " s " << " "
541  << std::setw(10) << "X(mm)" << " "
542  << std::setw(10) << "Y(mm)" << " "
543  << std::setw(10) << "Z(mm)" << " "
544  << std::setw( 7) << " N_x " << " "
545  << std::setw( 7) << " N_y " << " "
546  << std::setw( 7) << " N_z " << " " ;
547  G4cout << std::setw( 7) << " Delta|N|" << " "
548  << std::setw( 9) << "StepLen" << " "
549  << std::setw(12) << "StartSafety" << " "
550  << std::setw( 9) << "PhsStep" << " ";
551  if( startVolume )
552  { G4cout << std::setw(18) << "NextVolume" << " "; }
553  G4cout.precision(oldprec);
554  G4cout << G4endl;
555  }
556  if((stepNo == 0) && (verboseLevel <=3))
557  {
558  // Recurse to print the start values
559  //
560  printStatus( StartFT, StartFT, -1.0, safety, -1, startVolume);
561  }
562  if( verboseLevel <= 3 )
563  {
564  if( stepNo >= 0)
565  { G4cout << std::setw( 4) << stepNo << " "; }
566  else
567  { G4cout << std::setw( 5) << "Start" ; }
568  oldprec = G4cout.precision(8);
569  G4cout << std::setw(10) << CurrentFT.GetCurveLength() << " ";
570  G4cout.precision(8);
571  G4cout << std::setw(10) << CurrentPosition.x() << " "
572  << std::setw(10) << CurrentPosition.y() << " "
573  << std::setw(10) << CurrentPosition.z() << " ";
574  G4cout.precision(4);
575  G4cout << std::setw( 7) << CurrentUnitVelocity.x() << " "
576  << std::setw( 7) << CurrentUnitVelocity.y() << " "
577  << std::setw( 7) << CurrentUnitVelocity.z() << " ";
578  G4cout.precision(3);
579  G4cout << std::setw( 7)
580  << CurrentFT.GetMomentum().mag()-StartFT.GetMomentum().mag() << " ";
581  G4cout << std::setw( 9) << step_len << " ";
582  G4cout << std::setw(12) << safety << " ";
583  if( requestStep != -1.0 )
584  { G4cout << std::setw( 9) << requestStep << " "; }
585  else
586  { G4cout << std::setw( 9) << "Init/NotKnown" << " "; }
587  if( startVolume != 0)
588  { G4cout << std::setw(12) << startVolume->GetName() << " "; }
589  G4cout.precision(oldprec);
590  G4cout << G4endl;
591  }
592  else // if( verboseLevel > 3 )
593  {
594  // Multi-line output
595 
596  G4cout << "Step taken was " << step_len
597  << " out of PhysicalStep = " << requestStep << G4endl;
598  G4cout << "Final safety is: " << safety << G4endl;
599  G4cout << "Chord length = " << (CurrentPosition-StartPosition).mag()
600  << G4endl;
601  G4cout << G4endl;
602  }
603 }
double x() const
int G4int
Definition: G4Types.hh:78
double z() const
G4GLOB_DLL std::ostream G4cout
const G4String & GetName() const
double y() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void printStatus(const G4FieldTrack &startFT, const G4FieldTrack &currentFT, G4double requestStep, G4double safety, G4int step, G4VPhysicalVolume *startVolume)

Here is the call graph for this function:

Here is the caller graph for this function:

void G4PropagatorInField::PrintStepLengthDiagnostic ( G4double  currentProposedStepLength,
G4double  decreaseFactor,
G4double  stepTrial,
const G4FieldTrack aFieldTrack 
)
protected

Definition at line 610 of file G4PropagatorInField.cc.

615 {
616  G4int iprec= G4cout.precision(8);
617  G4cout << " " << std::setw(12) << " PiF: NoZeroStep "
618  << " " << std::setw(20) << " CurrentProposed len "
619  << " " << std::setw(18) << " Full_curvelen_last"
620  << " " << std::setw(18) << " last proposed len "
621  << " " << std::setw(18) << " decrease factor "
622  << " " << std::setw(15) << " step trial "
623  << G4endl;
624 
625  G4cout << " " << std::setw(10) << fNoZeroStep << " "
626  << " " << std::setw(20) << CurrentProposedStepLength
627  << " " << std::setw(18) << fFull_CurveLen_of_LastAttempt
628  << " " << std::setw(18) << fLast_ProposedStepLength
629  << " " << std::setw(18) << decreaseFactor
630  << " " << std::setw(15) << stepTrial
631  << G4endl;
632  G4cout.precision( iprec );
633 
634 }
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61

Here is the caller graph for this function:

void G4PropagatorInField::RefreshIntersectionLocator ( )

Definition at line 132 of file G4PropagatorInField.cc.

133 {
134  fIntersectionLocator->SetEpsilonStepFor(fEpsilonStep);
135  fIntersectionLocator->SetDeltaIntersectionFor(fCurrentFieldMgr->GetDeltaIntersection());
136  fIntersectionLocator->SetChordFinderFor(GetChordFinder());
137  fIntersectionLocator->SetSafetyParametersFor( fUseSafetyForOptimisation);
138 }
void SetChordFinderFor(G4ChordFinder *fCFinder)
G4double GetDeltaIntersection() const
G4ChordFinder * GetChordFinder()
void SetSafetyParametersFor(G4bool UseSafety)
void SetEpsilonStepFor(G4double EpsilonStep)
void SetDeltaIntersectionFor(G4double deltaIntersection)

Here is the call graph for this function:

Here is the caller graph for this function:

void G4PropagatorInField::ReportLoopingParticle ( G4int  count,
double  StepTaken,
G4VPhysicalVolume pPhysVol 
)
protected

Definition at line 736 of file G4PropagatorInField.cc.

739 {
740  std::ostringstream message;
741  message << " Killing looping particle "
742  << " after " << count << " field substeps "
743  << " totaling " << StepTaken / mm << " mm " ;
744  if( pPhysVol )
745  {
746  message << " in *volume* " << pPhysVol->GetName() ;
747  }
748  else
749  {
750  message << " in unknown or null volume. " ;
751  }
752  G4Exception("G4PropagatorInField::ComputeStep()", "GeomNav1002",
753  JustWarning, message);
754 }
static constexpr double mm
Definition: G4SIunits.hh:115
const G4String & GetName() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

Here is the call graph for this function:

Here is the caller graph for this function:

void G4PropagatorInField::ReportStuckParticle ( G4int  noZeroSteps,
G4double  proposedStep,
G4double  lastTriedStep,
G4VPhysicalVolume physVol 
)
protected

Definition at line 756 of file G4PropagatorInField.cc.

760 {
761  std::ostringstream message;
762  message << "Particle is stuck; it will be killed." << G4endl
763  << " Zero progress for " << noZeroSteps << " attempted steps."
764  << G4endl
765  << " Proposed Step is " << proposedStep
766  << " but Step Taken is "<< lastTriedStep << G4endl;
767  if( physVol )
768  message << " in volume " << physVol->GetName() ;
769  else
770  message << " in unknown or null volume. " ;
771  G4Exception("G4PropagatorInField::ComputeStep()",
772  "GeomNav1002", JustWarning, message);
773 }
const G4String & GetName() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

Here is the caller graph for this function:

void G4PropagatorInField::SetDetectorFieldManager ( G4FieldManager newGlobalFieldManager)
inline
void G4PropagatorInField::SetEpsilonStep ( G4double  newEps)
inline

Here is the caller graph for this function:

void G4PropagatorInField::SetIntersectionLocator ( G4VIntersectionLocator pLocator)
inline
void G4PropagatorInField::SetLargestAcceptableStep ( G4double  newBigDist)
inline

Here is the caller graph for this function:

void G4PropagatorInField::SetMaximumEpsilonStep ( G4double  newEpsMax)
inline

Here is the caller graph for this function:

void G4PropagatorInField::SetMaxLoopCount ( G4int  new_max)
inline
void G4PropagatorInField::SetMinimumEpsilonStep ( G4double  newEpsMin)
inline

Here is the caller graph for this function:

void G4PropagatorInField::SetNavigatorForPropagating ( G4Navigator SimpleOrMultiNavigator)
inline

Here is the caller graph for this function:

void G4PropagatorInField::SetThresholdNoZeroStep ( G4int  noAct,
G4int  noHarsh,
G4int  noAbandon 
)
inline
void G4PropagatorInField::SetTrajectoryFilter ( G4VCurvedTrajectoryFilter filter)

Definition at line 664 of file G4PropagatorInField.cc.

665 {
666  fpTrajectoryFilter = filter;
667 }

Here is the caller graph for this function:

void G4PropagatorInField::SetUseSafetyForOptimization ( G4bool  )
inline
G4int G4PropagatorInField::SetVerboseLevel ( G4int  verbose)

Definition at line 721 of file G4PropagatorInField.cc.

722 {
723  G4int oldval= fVerboseLevel;
724  fVerboseLevel= level;
725 
726  // Forward the verbose level 'reduced' to ChordFinder,
727  // MagIntegratorDriver ... ?
728  //
730  integrDriver->SetVerboseLevel( fVerboseLevel - 2 );
731  G4cout << "Set Driver verbosity to " << fVerboseLevel - 2 << G4endl;
732 
733  return oldval;
734 }
void SetVerboseLevel(G4int newLevel)
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
G4ChordFinder * GetChordFinder()
#define G4endl
Definition: G4ios.hh:61
G4MagInt_Driver * GetIntegrationDriver()

Here is the call graph for this function:

void G4PropagatorInField::SetVerboseTrace ( G4bool  enable)
inline

Definition at line 102 of file G4PropagatorInField.hh.

102 { fVerbTracePiF = enable; }
void G4PropagatorInField::SetZeroStepThreshold ( G4double  newLength)
inline
G4int G4PropagatorInField::Verbose ( ) const
inline

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