Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Member Functions
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
 
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)
 
G4VIntersectionLocatorGetIntersectionLocator ()
 
void SetIntersectionLocator (G4VIntersectionLocator *pLocator)
 
G4double GetDeltaIntersection () const
 
G4double GetDeltaOneStep () const
 
G4FieldManagerGetCurrentFieldManager ()
 
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)
 

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 59 of file G4PropagatorInField.cc.

References G4cout, G4endl, G4GeometryTolerance::GetInstance(), G4FieldManager::GetMaximumEpsilonStep(), G4GeometryTolerance::GetSurfaceTolerance(), G4INCL::Math::max(), python.hepunit::meter, python.hepunit::micrometer, python.hepunit::millimeter, and RefreshIntersectionLocator().

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

Definition at line 111 of file G4PropagatorInField.cc.

112 {
113  if(fAllocatedLocator)delete fIntersectionLocator;
114 }

Member Function Documentation

void G4PropagatorInField::ClearPropagatorState ( )

Definition at line 613 of file G4PropagatorInField.cc.

Referenced by G4ITTransportation::StartTracking(), G4CoupledTransportation::StartTracking(), G4MonopoleTransportation::StartTracking(), and G4Transportation::StartTracking().

614 {
615  // Goal: Clear all memory of previous steps, cached information
616 
617  fParticleIsLooping= false;
618  fNoZeroStep= 0;
619 
620  End_PointAndTangent= G4FieldTrack( G4ThreeVector(0.,0.,0.),
621  G4ThreeVector(0.,0.,0.),
622  0.0,0.0,0.0,0.0,0.0);
623  fFull_CurveLen_of_LastAttempt = -1;
624  fLast_ProposedStepLength = -1;
625 
626  fPreviousSftOrigin= G4ThreeVector(0.,0.,0.);
627  fPreviousSafety= 0.0;
628 }
CLHEP::Hep3Vector G4ThreeVector
G4double G4PropagatorInField::ComputeStep ( G4FieldTrack pFieldTrack,
G4double  pCurrentProposedStepLength,
G4double pNewSafety,
G4VPhysicalVolume pPhysVol = 0 
)

Definition at line 130 of file G4PropagatorInField.cc.

References G4ChordFinder::AdvanceChordLimited(), python.hepunit::cm, G4VCurvedTrajectoryFilter::CreateNewTrajectorySegment(), FatalException, FindAndSetFieldManager(), G4cout, G4endl, G4Exception(), GetChordFinder(), G4FieldTrack::GetCurveLength(), G4FieldManager::GetDeltaOneStep(), G4VPhysicalVolume::GetLogicalVolume(), G4FieldManager::GetMaximumEpsilonStep(), G4FieldManager::GetMinimumEpsilonStep(), G4FieldTrack::GetMomentumDir(), G4VPhysicalVolume::GetName(), G4FieldTrack::GetPosition(), G4Navigator::GetWorldVolume(), IntersectChord(), JustWarning, G4Navigator::LocateGlobalPointWithinVolume(), G4INCL::Math::min(), python.hepunit::mm, printStatus(), PrintStepLengthDiagnostic(), RefreshIntersectionLocator(), SetEpsilonStep(), and G4VCurvedTrajectoryFilter::TakeIntermediatePoint().

Referenced by G4Transportation::AlongStepGetPhysicalInteractionLength(), G4MonopoleTransportation::AlongStepGetPhysicalInteractionLength(), and G4PathFinder::DoNextCurvedStep().

135 {
136  // If CurrentProposedStepLength is too small for finding Chords
137  // then return with no action (for now - TODO: some action)
138  //
139  if(CurrentProposedStepLength<kCarTolerance)
140  {
141  return kInfinity;
142  }
143 
144  // Introducing smooth trajectory display (jacek 01/11/2002)
145  //
146  if (fpTrajectoryFilter)
147  {
148  fpTrajectoryFilter->CreateNewTrajectorySegment();
149  }
150 
151  // Parameters for adaptive Runge-Kutta integration
152 
153  G4double h_TrialStepSize; // 1st Step Size
154  G4double TruePathLength = CurrentProposedStepLength;
155  G4double StepTaken = 0.0;
156  G4double s_length_taken, epsilon ;
157  G4bool intersects;
158  G4bool first_substep = true;
159 
160  G4double NewSafety;
161  fParticleIsLooping = false;
162 
163  // If not yet done,
164  // Set the field manager to the local one if the volume has one,
165  // or to the global one if not
166  //
167  if( !fSetFieldMgr ) fCurrentFieldMgr= FindAndSetFieldManager( pPhysVol );
168  // For the next call, the field manager must again be set
169  fSetFieldMgr= false;
170 
171  // Values for Intersection Locator has to be updated on each call for the
172  // case that CurrentFieldManager has changed from the one of previous step
174 
175  G4FieldTrack CurrentState(pFieldTrack);
176  G4FieldTrack OriginalState = CurrentState;
177 
178  // If the Step length is "infinite", then an approximate-maximum Step
179  // length (used to calculate the relative accuracy) must be guessed.
180  //
181  if( CurrentProposedStepLength >= fLargestAcceptableStep )
182  {
183  G4ThreeVector StartPointA, VelocityUnit;
184  StartPointA = pFieldTrack.GetPosition();
185  VelocityUnit = pFieldTrack.GetMomentumDir();
186 
187  G4double trialProposedStep = 1.e2 * ( 10.0 * cm +
188  fNavigator->GetWorldVolume()->GetLogicalVolume()->
189  GetSolid()->DistanceToOut(StartPointA, VelocityUnit) );
190  CurrentProposedStepLength= std::min( trialProposedStep,
191  fLargestAcceptableStep );
192  }
193  epsilon = fCurrentFieldMgr->GetDeltaOneStep() / CurrentProposedStepLength;
194  // G4double raw_epsilon= epsilon;
195  G4double epsilonMin= fCurrentFieldMgr->GetMinimumEpsilonStep();
196  G4double epsilonMax= fCurrentFieldMgr->GetMaximumEpsilonStep();;
197  if( epsilon < epsilonMin ) epsilon = epsilonMin;
198  if( epsilon > epsilonMax ) epsilon = epsilonMax;
199  SetEpsilonStep( epsilon );
200 
201  // G4cout << "G4PiF: Epsilon of current step - raw= " << raw_epsilon
202  // << " final= " << epsilon << G4endl;
203 
204  // Shorten the proposed step in case of earlier problems (zero steps)
205  //
206  if( fNoZeroStep > fActionThreshold_NoZeroSteps )
207  {
208  G4double stepTrial;
209 
210  stepTrial= fFull_CurveLen_of_LastAttempt;
211  if( (stepTrial <= 0.0) && (fLast_ProposedStepLength > 0.0) )
212  stepTrial= fLast_ProposedStepLength;
213 
214  G4double decreaseFactor = 0.9; // Unused default
215  if( (fNoZeroStep < fSevereActionThreshold_NoZeroSteps)
216  && (stepTrial > 100.0*fZeroStepThreshold) )
217  {
218  // Attempt quick convergence
219  //
220  decreaseFactor= 0.25;
221  }
222  else
223  {
224  // We are in significant difficulties, probably at a boundary that
225  // is either geometrically sharp or between very different materials.
226  // Careful decreases to cope with tolerance are required.
227  //
228  if( stepTrial > 100.0*fZeroStepThreshold )
229  decreaseFactor = 0.35; // Try decreasing slower
230  else if( stepTrial > 30.0*fZeroStepThreshold )
231  decreaseFactor= 0.5; // Try yet slower decrease
232  else if( stepTrial > 10.0*fZeroStepThreshold )
233  decreaseFactor= 0.75; // Try even slower decreases
234  else
235  decreaseFactor= 0.9; // Try very slow decreases
236  }
237  stepTrial *= decreaseFactor;
238 
239 #ifdef G4DEBUG_FIELD
240  G4cout << " G4PropagatorInField::ComputeStep(): " << G4endl
241  << " Decreasing step - "
242  << " decreaseFactor= " << std::setw(8) << decreaseFactor
243  << " stepTrial = " << std::setw(18) << stepTrial << " "
244  << " fZeroStepThreshold = " << fZeroStepThreshold << G4endl;
245  PrintStepLengthDiagnostic(CurrentProposedStepLength, decreaseFactor,
246  stepTrial, pFieldTrack);
247 #endif
248  if( stepTrial == 0.0 ) // Change to make it < 0.1 * kCarTolerance ??
249  {
250  std::ostringstream message;
251  message << "Particle abandoned due to lack of progress in field."
252  << G4endl
253  << " Properties : " << pFieldTrack << G4endl
254  << " Attempting a zero step = " << stepTrial << G4endl
255  << " while attempting to progress after " << fNoZeroStep
256  << " trial steps. Will abandon step.";
257  G4Exception("G4PropagatorInField::ComputeStep()", "GeomNav1002",
258  JustWarning, message);
259  fParticleIsLooping= true;
260  return 0; // = stepTrial;
261  }
262  if( stepTrial < CurrentProposedStepLength )
263  CurrentProposedStepLength = stepTrial;
264  }
265  fLast_ProposedStepLength = CurrentProposedStepLength;
266 
267  G4int do_loop_count = 0;
268  do
269  {
270  G4FieldTrack SubStepStartState = CurrentState;
271  G4ThreeVector SubStartPoint = CurrentState.GetPosition();
272 
273  if( !first_substep) {
274  fNavigator->LocateGlobalPointWithinVolume( SubStartPoint );
275  }
276 
277  // How far to attempt to move the particle !
278  //
279  h_TrialStepSize = CurrentProposedStepLength - StepTaken;
280 
281  // Integrate as far as "chord miss" rule allows.
282  //
283  s_length_taken = GetChordFinder()->AdvanceChordLimited(
284  CurrentState, // Position & velocity
285  h_TrialStepSize,
286  fEpsilonStep,
287  fPreviousSftOrigin,
288  fPreviousSafety
289  );
290  // CurrentState is now updated with the final position and velocity.
291 
292  fFull_CurveLen_of_LastAttempt = s_length_taken;
293 
294  G4ThreeVector EndPointB = CurrentState.GetPosition();
295  G4ThreeVector InterSectionPointE;
296  G4double LinearStepLength;
297 
298  // Intersect chord AB with geometry
299  intersects= IntersectChord( SubStartPoint, EndPointB,
300  NewSafety, LinearStepLength,
301  InterSectionPointE );
302  // E <- Intersection Point of chord AB and either volume A's surface
303  // or a daughter volume's surface ..
304 
305  if( first_substep ) {
306  currentSafety = NewSafety;
307  } // Updating safety in other steps is potential future extention
308 
309  if( intersects )
310  {
311  G4FieldTrack IntersectPointVelct_G(CurrentState); // FT-Def-Construct
312 
313  // Find the intersection point of AB true path with the surface
314  // of vol(A), if it exists. Start with point E as first "estimate".
315  G4bool recalculatedEndPt= false;
316 
317  G4bool found_intersection = fIntersectionLocator->
318  EstimateIntersectionPoint( SubStepStartState, CurrentState,
319  InterSectionPointE, IntersectPointVelct_G,
320  recalculatedEndPt,fPreviousSafety,fPreviousSftOrigin);
321  intersects = intersects && found_intersection;
322  if( found_intersection ) {
323  End_PointAndTangent= IntersectPointVelct_G; // G is our EndPoint ...
324  StepTaken = TruePathLength = IntersectPointVelct_G.GetCurveLength()
325  - OriginalState.GetCurveLength();
326  } else {
327  // intersects= false; // "Minor" chords do not intersect
328  if( recalculatedEndPt ){
329  CurrentState= IntersectPointVelct_G;
330  }
331  }
332  }
333  if( !intersects )
334  {
335  StepTaken += s_length_taken;
336  // For smooth trajectory display (jacek 01/11/2002)
337  if (fpTrajectoryFilter) {
338  fpTrajectoryFilter->TakeIntermediatePoint(CurrentState.GetPosition());
339  }
340  }
341  first_substep = false;
342 
343 #ifdef G4DEBUG_FIELD
344  if( fNoZeroStep > fActionThreshold_NoZeroSteps ) {
345  printStatus( SubStepStartState, // or OriginalState,
346  CurrentState, CurrentProposedStepLength,
347  NewSafety, do_loop_count, pPhysVol );
348  }
349  if( (fVerboseLevel > 1) && (do_loop_count > fMax_loop_count-10 )) {
350  if( do_loop_count == fMax_loop_count-9 ){
351  G4cout << " G4PropagatorInField::ComputeStep(): " << G4endl
352  << " Difficult track - taking many sub steps." << G4endl;
353  }
354  printStatus( SubStepStartState, CurrentState, CurrentProposedStepLength,
355  NewSafety, do_loop_count, pPhysVol );
356  }
357 #endif
358 
359  do_loop_count++;
360 
361  } while( (!intersects )
362  && (StepTaken + kCarTolerance < CurrentProposedStepLength)
363  && ( do_loop_count < fMax_loop_count ) );
364 
365  if( do_loop_count >= fMax_loop_count )
366  {
367  fParticleIsLooping = true;
368 
369  if ( fVerboseLevel > 0 )
370  {
371  G4cout << " G4PropagateInField::ComputeStep(): " << G4endl
372  << " Killing looping particle "
373  // << " of " << energy << " energy "
374  << " after " << do_loop_count << " field substeps "
375  << " totaling " << StepTaken / mm << " mm " ;
376  if( pPhysVol )
377  G4cout << " in volume " << pPhysVol->GetName() ;
378  else
379  G4cout << " in unknown or null volume. " ;
380  G4cout << G4endl;
381  }
382  }
383 
384  if( !intersects )
385  {
386  // Chord AB or "minor chords" do not intersect
387  // B is the endpoint Step of the current Step.
388  //
389  End_PointAndTangent = CurrentState;
390  TruePathLength = StepTaken;
391  }
392 
393  // Set pFieldTrack to the return value
394  //
395  pFieldTrack = End_PointAndTangent;
396 
397 #ifdef G4VERBOSE
398  // Check that "s" is correct
399  //
400  if( std::fabs(OriginalState.GetCurveLength() + TruePathLength
401  - End_PointAndTangent.GetCurveLength()) > 3.e-4 * TruePathLength )
402  {
403  std::ostringstream message;
404  message << "Curve length mis-match between original state "
405  << "and proposed endpoint of propagation." << G4endl
406  << " The curve length of the endpoint should be: "
407  << OriginalState.GetCurveLength() + TruePathLength << G4endl
408  << " and it is instead: "
409  << End_PointAndTangent.GetCurveLength() << "." << G4endl
410  << " A difference of: "
411  << OriginalState.GetCurveLength() + TruePathLength
412  - End_PointAndTangent.GetCurveLength() << G4endl
413  << " Original state = " << OriginalState << G4endl
414  << " Proposed state = " << End_PointAndTangent;
415  G4Exception("G4PropagatorInField::ComputeStep()",
416  "GeomNav0003", FatalException, message);
417  }
418 #endif
419 
420  // In particular anomalous cases, we can get repeated zero steps
421  // In order to correct this efficiently, we identify these cases
422  // and only take corrective action when they occur.
423  //
424  if( ( (TruePathLength < fZeroStepThreshold)
425  && ( TruePathLength+kCarTolerance < CurrentProposedStepLength )
426  )
427  || ( TruePathLength < 0.5*kCarTolerance )
428  )
429  {
430  fNoZeroStep++;
431  }
432  else{
433  fNoZeroStep = 0;
434  }
435 
436  if( fNoZeroStep > fAbandonThreshold_NoZeroSteps )
437  {
438  fParticleIsLooping = true;
439  std::ostringstream message;
440  message << "Particle is stuck; it will be killed." << G4endl
441  << " Zero progress for " << fNoZeroStep << " attempted steps."
442  << G4endl
443  << " Proposed Step is " << CurrentProposedStepLength
444  << " but Step Taken is "<< fFull_CurveLen_of_LastAttempt << G4endl;
445  if( pPhysVol )
446  message << " in volume " << pPhysVol->GetName() ;
447  else
448  message << " in unknown or null volume. " ;
449  G4Exception("G4PropagatorInField::ComputeStep()",
450  "GeomNav1002", JustWarning, message);
451  fNoZeroStep = 0;
452  }
453 
454  return TruePathLength;
455 }
void SetEpsilonStep(G4double newEps)
G4double GetCurveLength() const
virtual void TakeIntermediatePoint(G4ThreeVector newPoint)=0
void PrintStepLengthDiagnostic(G4double currentProposedStepLength, G4double decreaseFactor, G4double stepTrial, const G4FieldTrack &aFieldTrack)
const G4ThreeVector & GetMomentumDir() const
G4double GetDeltaOneStep() const
int G4int
Definition: G4Types.hh:78
G4ThreeVector GetPosition() const
G4GLOB_DLL std::ostream G4cout
const G4String & GetName() const
bool G4bool
Definition: G4Types.hh:79
G4double GetMaximumEpsilonStep() const
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
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
virtual void LocateGlobalPointWithinVolume(const G4ThreeVector &position)
Definition: G4Navigator.cc:550
void printStatus(const G4FieldTrack &startFT, const G4FieldTrack &currentFT, G4double requestStep, G4double safety, G4int step, G4VPhysicalVolume *startVolume)
G4double GetMinimumEpsilonStep() const
G4ThreeVector G4PropagatorInField::EndMomentumDir ( ) const
inline
G4ThreeVector G4PropagatorInField::EndPosition ( ) const
inline
G4FieldManager * G4PropagatorInField::FindAndSetFieldManager ( G4VPhysicalVolume pCurrentPhysVol)

Definition at line 631 of file G4PropagatorInField.cc.

References G4Region::GetFieldManager(), G4LogicalVolume::GetFieldManager(), G4VPhysicalVolume::GetLogicalVolume(), and G4LogicalVolume::GetRegion().

Referenced by G4Transportation::AlongStepGetPhysicalInteractionLength(), G4CoupledTransportation::AlongStepGetPhysicalInteractionLength(), G4MonopoleTransportation::AlongStepGetPhysicalInteractionLength(), G4ITTransportation::AlongStepGetPhysicalInteractionLength(), G4PathFinder::ComputeStep(), ComputeStep(), G4SynchrotronRadiation::GetMeanFreePath(), G4SynchrotronRadiationInMat::GetMeanFreePath(), G4SynchrotronRadiationInMat::GetPhotonEnergy(), G4SynchrotronRadiation::PostStepDoIt(), and G4SynchrotronRadiationInMat::PostStepDoIt().

632 {
633  G4FieldManager* currentFieldMgr;
634 
635  currentFieldMgr = fDetectorFieldMgr;
636  if( pCurrentPhysicalVolume)
637  {
638  G4FieldManager *pRegionFieldMgr= 0, *localFieldMgr = 0;
639  G4LogicalVolume* pLogicalVol= pCurrentPhysicalVolume->GetLogicalVolume();
640 
641  if( pLogicalVol ) {
642  // Value for Region, if any, Overrides
643  G4Region* pRegion= pLogicalVol->GetRegion();
644  if( pRegion ) {
645  pRegionFieldMgr= pRegion->GetFieldManager();
646  if( pRegionFieldMgr )
647  currentFieldMgr= pRegionFieldMgr;
648  }
649 
650  // 'Local' Value from logical volume, if any, Overrides
651  localFieldMgr= pLogicalVol->GetFieldManager();
652  if ( localFieldMgr )
653  currentFieldMgr = localFieldMgr;
654  }
655  }
656  fCurrentFieldMgr= currentFieldMgr;
657 
658  // Flag that field manager has been set.
659  fSetFieldMgr= true;
660 
661  return currentFieldMgr;
662 }
G4Region * GetRegion() const
G4FieldManager * GetFieldManager() const
G4FieldManager * GetFieldManager() const
G4ChordFinder* G4PropagatorInField::GetChordFinder ( )
inline
G4FieldManager* G4PropagatorInField::GetCurrentFieldManager ( )
inline
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
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
G4double G4PropagatorInField::GetZeroStepThreshold ( )
inline
std::vector< G4ThreeVector > * G4PropagatorInField::GimmeTrajectoryVectorAndForgetIt ( ) const

Definition at line 592 of file G4PropagatorInField.cc.

References G4VCurvedTrajectoryFilter::GimmeThePointsAndForgetThem().

Referenced by G4Transportation::AlongStepDoIt(), G4MonopoleTransportation::AlongStepDoIt(), and G4ITTransportation::AlongStepDoIt().

593 {
594  // NB, GimmeThePointsAndForgetThem really forgets them, so it can
595  // only be called (exactly) once for each step.
596 
597  if (fpTrajectoryFilter)
598  {
599  return fpTrajectoryFilter->GimmeThePointsAndForgetThem();
600  }
601  else
602  {
603  return 0;
604  }
605 }
std::vector< G4ThreeVector > * GimmeThePointsAndForgetThem()
G4bool G4PropagatorInField::IntersectChord ( const G4ThreeVector StartPointA,
const G4ThreeVector EndPointB,
G4double NewSafety,
G4double LinearStepLength,
G4ThreeVector IntersectionPoint 
)
inline

Referenced by ComputeStep().

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

Definition at line 462 of file G4PropagatorInField.cc.

References G4cout, G4endl, G4FieldTrack::GetCurveLength(), G4FieldTrack::GetMomentum(), G4FieldTrack::GetMomentumDir(), G4VPhysicalVolume::GetName(), G4FieldTrack::GetPosition(), CLHEP::Hep3Vector::mag(), CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

Referenced by ComputeStep().

468 {
469  const G4int verboseLevel=fVerboseLevel;
470  const G4ThreeVector StartPosition = StartFT.GetPosition();
471  const G4ThreeVector StartUnitVelocity = StartFT.GetMomentumDir();
472  const G4ThreeVector CurrentPosition = CurrentFT.GetPosition();
473  const G4ThreeVector CurrentUnitVelocity = CurrentFT.GetMomentumDir();
474 
475  G4double step_len = CurrentFT.GetCurveLength() - StartFT.GetCurveLength();
476 
477  G4int oldprec; // cout/cerr precision settings
478 
479  if( ((stepNo == 0) && (verboseLevel <3)) || (verboseLevel >= 3) )
480  {
481  oldprec = G4cout.precision(4);
482  G4cout << std::setw( 6) << " "
483  << std::setw( 25) << " Current Position and Direction" << " "
484  << G4endl;
485  G4cout << std::setw( 5) << "Step#"
486  << std::setw(10) << " s " << " "
487  << std::setw(10) << "X(mm)" << " "
488  << std::setw(10) << "Y(mm)" << " "
489  << std::setw(10) << "Z(mm)" << " "
490  << std::setw( 7) << " N_x " << " "
491  << std::setw( 7) << " N_y " << " "
492  << std::setw( 7) << " N_z " << " " ;
493  G4cout << std::setw( 7) << " Delta|N|" << " "
494  << std::setw( 9) << "StepLen" << " "
495  << std::setw(12) << "StartSafety" << " "
496  << std::setw( 9) << "PhsStep" << " ";
497  if( startVolume )
498  { G4cout << std::setw(18) << "NextVolume" << " "; }
499  G4cout.precision(oldprec);
500  G4cout << G4endl;
501  }
502  if((stepNo == 0) && (verboseLevel <=3))
503  {
504  // Recurse to print the start values
505  //
506  printStatus( StartFT, StartFT, -1.0, safety, -1, startVolume);
507  }
508  if( verboseLevel <= 3 )
509  {
510  if( stepNo >= 0)
511  { G4cout << std::setw( 4) << stepNo << " "; }
512  else
513  { G4cout << std::setw( 5) << "Start" ; }
514  oldprec = G4cout.precision(8);
515  G4cout << std::setw(10) << CurrentFT.GetCurveLength() << " ";
516  G4cout.precision(8);
517  G4cout << std::setw(10) << CurrentPosition.x() << " "
518  << std::setw(10) << CurrentPosition.y() << " "
519  << std::setw(10) << CurrentPosition.z() << " ";
520  G4cout.precision(4);
521  G4cout << std::setw( 7) << CurrentUnitVelocity.x() << " "
522  << std::setw( 7) << CurrentUnitVelocity.y() << " "
523  << std::setw( 7) << CurrentUnitVelocity.z() << " ";
524  G4cout.precision(3);
525  G4cout << std::setw( 7)
526  << CurrentFT.GetMomentum().mag()-StartFT.GetMomentum().mag() << " ";
527  G4cout << std::setw( 9) << step_len << " ";
528  G4cout << std::setw(12) << safety << " ";
529  if( requestStep != -1.0 )
530  { G4cout << std::setw( 9) << requestStep << " "; }
531  else
532  { G4cout << std::setw( 9) << "Init/NotKnown" << " "; }
533  if( startVolume != 0)
534  { G4cout << std::setw(12) << startVolume->GetName() << " "; }
535  G4cout.precision(oldprec);
536  G4cout << G4endl;
537  }
538  else // if( verboseLevel > 3 )
539  {
540  // Multi-line output
541 
542  G4cout << "Step taken was " << step_len
543  << " out of PhysicalStep = " << requestStep << G4endl;
544  G4cout << "Final safety is: " << safety << G4endl;
545  G4cout << "Chord length = " << (CurrentPosition-StartPosition).mag()
546  << G4endl;
547  G4cout << G4endl;
548  }
549 }
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)
void G4PropagatorInField::PrintStepLengthDiagnostic ( G4double  currentProposedStepLength,
G4double  decreaseFactor,
G4double  stepTrial,
const G4FieldTrack aFieldTrack 
)
protected

Definition at line 556 of file G4PropagatorInField.cc.

References G4cout, and G4endl.

Referenced by ComputeStep().

561 {
562  G4int iprec= G4cout.precision(8);
563  G4cout << " " << std::setw(12) << " PiF: NoZeroStep "
564  << " " << std::setw(20) << " CurrentProposed len "
565  << " " << std::setw(18) << " Full_curvelen_last"
566  << " " << std::setw(18) << " last proposed len "
567  << " " << std::setw(18) << " decrease factor "
568  << " " << std::setw(15) << " step trial "
569  << G4endl;
570 
571  G4cout << " " << std::setw(10) << fNoZeroStep << " "
572  << " " << std::setw(20) << CurrentProposedStepLength
573  << " " << std::setw(18) << fFull_CurveLen_of_LastAttempt
574  << " " << std::setw(18) << fLast_ProposedStepLength
575  << " " << std::setw(18) << decreaseFactor
576  << " " << std::setw(15) << stepTrial
577  << G4endl;
578  G4cout.precision( iprec );
579 
580 }
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
void G4PropagatorInField::RefreshIntersectionLocator ( )

Definition at line 118 of file G4PropagatorInField.cc.

References GetChordFinder(), G4FieldManager::GetDeltaIntersection(), G4VIntersectionLocator::SetChordFinderFor(), G4VIntersectionLocator::SetDeltaIntersectionFor(), G4VIntersectionLocator::SetEpsilonStepFor(), and G4VIntersectionLocator::SetSafetyParametersFor().

Referenced by ComputeStep(), and G4PropagatorInField().

119 {
120  fIntersectionLocator->SetEpsilonStepFor(fEpsilonStep);
121  fIntersectionLocator->SetDeltaIntersectionFor(fCurrentFieldMgr->GetDeltaIntersection());
122  fIntersectionLocator->SetChordFinderFor(GetChordFinder());
123  fIntersectionLocator->SetSafetyParametersFor( fUseSafetyForOptimisation);
124 }
void SetChordFinderFor(G4ChordFinder *fCFinder)
G4double GetDeltaIntersection() const
G4ChordFinder * GetChordFinder()
void SetSafetyParametersFor(G4bool UseSafety)
void SetEpsilonStepFor(G4double EpsilonStep)
void SetDeltaIntersectionFor(G4double deltaIntersection)
void G4PropagatorInField::SetDetectorFieldManager ( G4FieldManager newGlobalFieldManager)
inline
void G4PropagatorInField::SetEpsilonStep ( G4double  newEps)
inline

Referenced by ComputeStep().

void G4PropagatorInField::SetIntersectionLocator ( G4VIntersectionLocator pLocator)
inline
void G4PropagatorInField::SetLargestAcceptableStep ( G4double  newBigDist)
inline
void G4PropagatorInField::SetMaximumEpsilonStep ( G4double  newEpsMax)
inline
void G4PropagatorInField::SetMaxLoopCount ( G4int  new_max)
inline
void G4PropagatorInField::SetMinimumEpsilonStep ( G4double  newEpsMin)
inline
void G4PropagatorInField::SetNavigatorForPropagating ( G4Navigator SimpleOrMultiNavigator)
inline
void G4PropagatorInField::SetThresholdNoZeroStep ( G4int  noAct,
G4int  noHarsh,
G4int  noAbandon 
)
inline
void G4PropagatorInField::SetTrajectoryFilter ( G4VCurvedTrajectoryFilter filter)

Definition at line 608 of file G4PropagatorInField.cc.

References filter.

Referenced by G4TrackingMessenger::SetNewValue().

609 {
610  fpTrajectoryFilter = filter;
611 }
pid_t filter
Definition: tracer.cxx:30
void G4PropagatorInField::SetUseSafetyForOptimization ( G4bool  )
inline
G4int G4PropagatorInField::SetVerboseLevel ( G4int  verbose)

Definition at line 664 of file G4PropagatorInField.cc.

References G4cout, G4endl, GetChordFinder(), G4ChordFinder::GetIntegrationDriver(), and G4MagInt_Driver::SetVerboseLevel().

665 {
666  G4int oldval= fVerboseLevel;
667  fVerboseLevel= level;
668 
669  // Forward the verbose level 'reduced' to ChordFinder,
670  // MagIntegratorDriver ... ?
671  //
673  integrDriver->SetVerboseLevel( fVerboseLevel - 2 );
674  G4cout << "Set Driver verbosity to " << fVerboseLevel - 2 << G4endl;
675 
676  return oldval;
677 }
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()
void G4PropagatorInField::SetZeroStepThreshold ( G4double  newLength)
inline
G4int G4PropagatorInField::Verbose ( ) const
inline

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