Geant4  10.02.p03
G4Transportation Class Reference

#include <G4Transportation.hh>

Inheritance diagram for G4Transportation:
Collaboration diagram for G4Transportation:

Public Member Functions

 G4Transportation (G4int verbosityLevel=1)
 
 ~G4Transportation ()
 
G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
 
G4VParticleChange * AlongStepDoIt (const G4Track &track, const G4Step &stepData)
 
G4VParticleChange * PostStepDoIt (const G4Track &track, const G4Step &stepData)
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &, G4double previousStepSize, G4ForceCondition *pForceCond)
 
G4bool FieldExertedForce ()
 
G4PropagatorInFieldGetPropagatorInField ()
 
void SetPropagatorInField (G4PropagatorInField *pFieldPropagator)
 
void SetVerboseLevel (G4int verboseLevel)
 
G4int GetVerboseLevel () const
 
G4double GetThresholdWarningEnergy () const
 
G4double GetThresholdImportantEnergy () const
 
G4int GetThresholdTrials () const
 
void SetThresholdWarningEnergy (G4double newEnWarn)
 
void SetThresholdImportantEnergy (G4double newEnImp)
 
void SetThresholdTrials (G4int newMaxTrials)
 
G4double GetMaxEnergyKilled () const
 
G4double GetSumEnergyKilled () const
 
void ResetKilledStatistics (G4int report=1)
 
void EnableShortStepOptimisation (G4bool optimise=true)
 
G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
G4VParticleChange * AtRestDoIt (const G4Track &, const G4Step &)
 
void StartTracking (G4Track *aTrack)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4int operator== (const G4VProcess &right) const
 
G4int operator!= (const G4VProcess &right) const
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Static Public Member Functions

static G4bool EnableUseMagneticMoment (G4bool useMoment=true)
 
- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 

Protected Member Functions

G4bool DoesGlobalFieldExist ()
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Private Attributes

G4NavigatorfLinearNavigator
 
G4PropagatorInFieldfFieldPropagator
 
G4ThreeVector fTransportEndPosition
 
G4ThreeVector fTransportEndMomentumDir
 
G4double fTransportEndKineticEnergy
 
G4ThreeVector fTransportEndSpin
 
G4bool fMomentumChanged
 
G4bool fEndGlobalTimeComputed
 
G4double fCandidateEndGlobalTime
 
G4bool fParticleIsLooping
 
G4bool fNewTrack
 
G4bool fFirstStepInVolume
 
G4bool fLastStepInVolume
 
G4bool fGeometryLimitedStep
 
G4bool fFieldExertedForce
 
G4TouchableHandle fCurrentTouchableHandle
 
G4ThreeVector fPreviousSftOrigin
 
G4double fPreviousSafety
 
G4ParticleChangeForTransport fParticleChange
 
G4double fEndPointDistance
 
G4double fThreshold_Warning_Energy
 
G4double fThreshold_Important_Energy
 
G4int fThresholdTrials
 
G4int fNoLooperTrials
 
G4double fSumEnergyKilled
 
G4double fMaxEnergyKilled
 
G4bool fShortStepOptimisation
 
G4SafetyHelperfpSafetyHelper
 
G4int fVerboseLevel
 

Static Private Attributes

static G4bool fUseMagneticMoment =false
 

Friends

class G4CoupledTransportation
 

Additional Inherited Members

- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChange * pParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 

Detailed Description

Definition at line 59 of file G4Transportation.hh.

Constructor & Destructor Documentation

◆ G4Transportation()

G4Transportation::G4Transportation ( G4int  verbosityLevel = 1)

Definition at line 83 of file G4Transportation.cc.

84  : G4VProcess( G4String("Transportation"), fTransportation ),
85  fTransportEndPosition( 0.0, 0.0, 0.0 ),
86  fTransportEndMomentumDir( 0.0, 0.0, 0.0 ),
88  fTransportEndSpin( 0.0, 0.0, 0.0 ),
89  fMomentumChanged(true),
90  fEndGlobalTimeComputed(false),
92  fParticleIsLooping( false ),
93  fNewTrack( true ),
94  fFirstStepInVolume( true ),
95  fLastStepInVolume( false ),
97  fFieldExertedForce( false ),
98  fPreviousSftOrigin( 0.,0.,0. ),
99  fPreviousSafety( 0.0 ),
100  // fParticleChange(),
101  fEndPointDistance( -1.0 ),
102  fThreshold_Warning_Energy( 100 * MeV ),
104  fThresholdTrials( 10 ),
105  fNoLooperTrials( 0 ),
106  fSumEnergyKilled( 0.0 ), fMaxEnergyKilled( 0.0 ),
107  fShortStepOptimisation( false ), // Old default: true (=fast short steps)
108  fVerboseLevel( verbosity )
109 {
110  // set Process Sub Type
111  SetProcessSubType(static_cast<G4int>(TRANSPORTATION));
112  pParticleChange= &fParticleChange; // Required to conform to G4VProcess
113 
114  G4TransportationManager* transportMgr ;
115 
117 
118  fLinearNavigator = transportMgr->GetNavigatorForTracking() ;
119 
120  fFieldPropagator = transportMgr->GetPropagatorInField() ;
121 
122  fpSafetyHelper = transportMgr->GetSafetyHelper(); // New
123 
124  // Cannot determine whether a field exists here, as it would
125  // depend on the relative order of creating the detector's
126  // field and this process. That order is not guaranted.
127  // Instead later the method DoesGlobalFieldExist() is called
128 
129  static G4ThreadLocal G4TouchableHandle* pNullTouchableHandle = 0;
130  if ( !pNullTouchableHandle) { pNullTouchableHandle = new G4TouchableHandle; }
131  fCurrentTouchableHandle = *pNullTouchableHandle;
132  // Points to (G4VTouchable*) 0
133 
134 #ifdef G4VERBOSE
135  if( fVerboseLevel > 0)
136  {
137  G4cout << " G4Transportation constructor> set fShortStepOptimisation to ";
138  if ( fShortStepOptimisation ) G4cout << "true" << G4endl;
139  else G4cout << "false" << G4endl;
140  }
141 #endif
142 }
static const double MeV
Definition: G4SIunits.hh:211
G4PropagatorInField * fFieldPropagator
G4VProcess(const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
Definition: G4VProcess.cc:52
G4double fTransportEndKineticEnergy
G4ThreeVector fTransportEndPosition
#define G4ThreadLocal
Definition: tls.hh:89
G4SafetyHelper * fpSafetyHelper
G4ReferenceCountedHandle< G4VTouchable > G4TouchableHandle
G4GLOB_DLL std::ostream G4cout
G4double fThreshold_Important_Energy
G4SafetyHelper * GetSafetyHelper() const
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
G4Navigator * GetNavigatorForTracking() const
G4ThreeVector fPreviousSftOrigin
G4TouchableHandle fCurrentTouchableHandle
static G4TransportationManager * GetTransportationManager()
G4double fThreshold_Warning_Energy
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
#define G4endl
Definition: G4ios.hh:61
G4ThreeVector fTransportEndSpin
G4Navigator * fLinearNavigator
G4PropagatorInField * GetPropagatorInField() const
G4double fCandidateEndGlobalTime
G4ParticleChangeForTransport fParticleChange
G4ThreeVector fTransportEndMomentumDir
Here is the call graph for this function:

◆ ~G4Transportation()

G4Transportation::~G4Transportation ( )

Definition at line 146 of file G4Transportation.cc.

147 {
148  if( (fVerboseLevel > 0) && (fSumEnergyKilled > 0.0 ) )
149  {
150  G4cout << " G4Transportation: Statistics for looping particles " << G4endl;
151  G4cout << " Sum of energy of loopers killed: " << fSumEnergyKilled << G4endl;
152  G4cout << " Max energy of loopers killed: " << fMaxEnergyKilled << G4endl;
153  }
154 }
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
Here is the call graph for this function:

Member Function Documentation

◆ AlongStepDoIt()

G4VParticleChange * G4Transportation::AlongStepDoIt ( const G4Track &  track,
const G4Step &  stepData 
)
virtual

Implements G4VProcess.

Definition at line 525 of file G4Transportation.cc.

527 {
528  static G4ThreadLocal G4int noCalls=0;
529  noCalls++;
530 
531  fParticleChange.Initialize(track) ;
532 
533  // Code for specific process
534  //
535  fParticleChange.ProposePosition(fTransportEndPosition) ;
536  fParticleChange.ProposeMomentumDirection(fTransportEndMomentumDir) ;
538  fParticleChange.SetMomentumChanged(fMomentumChanged) ;
539 
540  fParticleChange.ProposePolarization(fTransportEndSpin);
541 
542  G4double deltaTime = 0.0 ;
543 
544  // Calculate Lab Time of Flight (ONLY if field Equations used it!)
545  // G4double endTime = fCandidateEndGlobalTime;
546  // G4double delta_time = endTime - startTime;
547 
548  G4double startTime = track.GetGlobalTime() ;
549 
551  {
552  // The time was not integrated .. make the best estimate possible
553  //
554  G4double initialVelocity = stepData.GetPreStepPoint()->GetVelocity();
555  G4double stepLength = track.GetStepLength();
556 
557  deltaTime= 0.0; // in case initialVelocity = 0
558  if ( initialVelocity > 0.0 ) { deltaTime = stepLength/initialVelocity; }
559 
560  fCandidateEndGlobalTime = startTime + deltaTime ;
561  fParticleChange.ProposeLocalTime( track.GetLocalTime() + deltaTime) ;
562  }
563  else
564  {
565  deltaTime = fCandidateEndGlobalTime - startTime ;
566  fParticleChange.ProposeGlobalTime( fCandidateEndGlobalTime ) ;
567  }
568 
569 
570  // Now Correct by Lorentz factor to get delta "proper" Time
571 
572  G4double restMass = track.GetDynamicParticle()->GetMass() ;
573  G4double deltaProperTime = deltaTime*( restMass/track.GetTotalEnergy() ) ;
574 
575  fParticleChange.ProposeProperTime(track.GetProperTime() + deltaProperTime) ;
576  //fParticleChange. ProposeTrueStepLength( track.GetStepLength() ) ;
577 
578  // If the particle is caught looping or is stuck (in very difficult
579  // boundaries) in a magnetic field (doing many steps)
580  // THEN this kills it ...
581  //
582  if ( fParticleIsLooping )
583  {
585 
586  if( (endEnergy < fThreshold_Important_Energy)
588  {
589  // Kill the looping particle
590  //
591  fParticleChange.ProposeTrackStatus( fStopAndKill ) ;
592 
593  // 'Bare' statistics
594  fSumEnergyKilled += endEnergy;
595  if( endEnergy > fMaxEnergyKilled) { fMaxEnergyKilled= endEnergy; }
596 
597 #ifdef G4VERBOSE
598  if( (fVerboseLevel > 1) &&
599  ( endEnergy > fThreshold_Warning_Energy ) )
600  {
601  G4cout << " G4Transportation is killing track that is looping or stuck "
602  << G4endl
603  << " This track has " << track.GetKineticEnergy() / MeV
604  << " MeV energy." << G4endl;
605  G4cout << " Number of trials = " << fNoLooperTrials
606  << " No of calls to AlongStepDoIt = " << noCalls
607  << G4endl;
608  }
609 #endif
610  fNoLooperTrials=0;
611  }
612  else
613  {
614  fNoLooperTrials ++;
615 #ifdef G4VERBOSE
616  if( (fVerboseLevel > 2) )
617  {
618  G4cout << " G4Transportation::AlongStepDoIt(): Particle looping - "
619  << " Number of trials = " << fNoLooperTrials
620  << " No of calls to = " << noCalls
621  << G4endl;
622  }
623 #endif
624  }
625  }
626  else
627  {
628  fNoLooperTrials=0;
629  }
630 
631  // Another (sometimes better way) is to use a user-limit maximum Step size
632  // to alleviate this problem ..
633 
634  // Introduce smooth curved trajectories to particle-change
635  //
636  fParticleChange.SetPointerToVectorOfAuxiliaryPoints
638 
639  return &fParticleChange ;
640 }
static const double MeV
Definition: G4SIunits.hh:211
G4PropagatorInField * fFieldPropagator
G4double fTransportEndKineticEnergy
G4ThreeVector fTransportEndPosition
#define G4ThreadLocal
Definition: tls.hh:89
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
G4double fThreshold_Important_Energy
std::vector< G4ThreeVector > * GimmeTrajectoryVectorAndForgetIt() const
G4double fThreshold_Warning_Energy
#define G4endl
Definition: G4ios.hh:61
G4ThreeVector fTransportEndSpin
double G4double
Definition: G4Types.hh:76
G4double fCandidateEndGlobalTime
G4ParticleChangeForTransport fParticleChange
G4ThreeVector fTransportEndMomentumDir
Here is the call graph for this function:

◆ AlongStepGetPhysicalInteractionLength()

G4double G4Transportation::AlongStepGetPhysicalInteractionLength ( const G4Track &  track,
G4double  previousStepSize,
G4double  currentMinimumStep,
G4double currentSafety,
G4GPILSelection *  selection 
)
virtual

Implements G4VProcess.

Definition at line 164 of file G4Transportation.cc.

169 {
170  G4double geometryStepLength= -1.0, newSafety= -1.0;
171  fParticleIsLooping = false ;
172 
173  // Initial actions moved to StartTrack()
174  // --------------------------------------
175  // Note: in case another process changes touchable handle
176  // it will be necessary to add here (for all steps)
177  // fCurrentTouchableHandle = aTrack->GetTouchableHandle();
178 
179  // GPILSelection is set to defaule value of CandidateForSelection
180  // It is a return value
181  //
182  *selection = CandidateForSelection ;
183 
185  // G4cout << " Transport::AlongStep GPIL: 1st-step= " << fFirstStepInVolume << " newTrack= " << fNewTrack << " fLastStep-in-Vol= " << fLastStepInVolume << G4endl;
186  fLastStepInVolume= false;
187  fNewTrack = false;
188 
189  fParticleChange.ProposeFirstStepInVolume(fFirstStepInVolume);
190 
191  // Get initial Energy/Momentum of the track
192  //
193  const G4DynamicParticle* pParticle = track.GetDynamicParticle() ;
194  const G4ParticleDefinition* pParticleDef = pParticle->GetDefinition() ;
195  G4ThreeVector startMomentumDir = pParticle->GetMomentumDirection() ;
196  G4ThreeVector startPosition = track.GetPosition() ;
197 
198  // G4double theTime = track.GetGlobalTime() ;
199 
200  // The Step Point safety can be limited by other geometries and/or the
201  // assumptions of any process - it's not always the geometrical safety.
202  // We calculate the starting point's isotropic safety here.
203  //
204  G4ThreeVector OriginShift = startPosition - fPreviousSftOrigin ;
205  G4double MagSqShift = OriginShift.mag2() ;
206  if( MagSqShift >= sqr(fPreviousSafety) )
207  {
208  currentSafety = 0.0 ;
209  }
210  else
211  {
212  currentSafety = fPreviousSafety - std::sqrt(MagSqShift) ;
213  }
214 
215  // Is the particle charged or has it a magnetic moment?
216  //
217  G4double particleCharge = pParticle->GetCharge() ;
218  G4double magneticMoment = pParticle->GetMagneticMoment() ;
219  G4double restMass = pParticle->GetMass() ;
220 
221  fGeometryLimitedStep = false ;
222  // fEndGlobalTimeComputed = false ;
223 
224  // There is no need to locate the current volume. It is Done elsewhere:
225  // On track construction
226  // By the tracking, after all AlongStepDoIts, in "Relocation"
227 
228  // Check if the particle has a force, EM or gravitational, exerted on it
229  //
230  G4FieldManager* fieldMgr=0;
231  G4bool fieldExertsForce = false ;
232 
233  G4bool gravityOn = false;
234  G4bool fieldExists= false; // Field is not 0 (null pointer)
235 
236  fieldMgr = fFieldPropagator->FindAndSetFieldManager( track.GetVolume() );
237  if( fieldMgr != 0 )
238  {
239  // Message the field Manager, to configure it for this track
240  fieldMgr->ConfigureForTrack( &track );
241  // Is here to allow a transition from no-field pointer
242  // to finite field (non-zero pointer).
243 
244  // If the field manager has no field ptr, the field is zero
245  // by definition ( = there is no field ! )
246  const G4Field* ptrField= fieldMgr->GetDetectorField();
247  fieldExists = (ptrField!=0) ;
248  if( fieldExists )
249  {
250  gravityOn= ptrField->IsGravityActive();
251 
252  if( (particleCharge != 0.0)
253  || (fUseMagneticMoment && (magneticMoment != 0.0) )
254  || (gravityOn && (restMass != 0.0) )
255  )
256  {
257  fieldExertsForce = fieldExists;
258  }
259  }
260  }
261  // G4cout << " G4Transport: field exerts force= " << fieldExertsForce
262  // << " fieldMgr= " << fieldMgr << G4endl;
263  fFieldExertedForce = fieldExertsForce;
264 
265  if( !fieldExertsForce )
266  {
267  G4double linearStepLength ;
268  if( fShortStepOptimisation && (currentMinimumStep <= currentSafety) )
269  {
270  // The Step is guaranteed to be taken
271  //
272  geometryStepLength = currentMinimumStep ;
273  fGeometryLimitedStep = false ;
274  }
275  else
276  {
277  // Find whether the straight path intersects a volume
278  //
279  linearStepLength = fLinearNavigator->ComputeStep( startPosition,
280  startMomentumDir,
281  currentMinimumStep,
282  newSafety) ;
283  // Remember last safety origin & value.
284  //
285  fPreviousSftOrigin = startPosition ;
286  fPreviousSafety = newSafety ;
287  fpSafetyHelper->SetCurrentSafety( newSafety, startPosition);
288 
289  currentSafety = newSafety ;
290 
291  fGeometryLimitedStep= (linearStepLength <= currentMinimumStep);
293  {
294  // The geometry limits the Step size (an intersection was found.)
295  geometryStepLength = linearStepLength ;
296  }
297  else
298  {
299  // The full Step is taken.
300  geometryStepLength = currentMinimumStep ;
301  }
302  }
303  fEndPointDistance = geometryStepLength ;
304 
305  // Calculate final position
306  //
307  fTransportEndPosition = startPosition+geometryStepLength*startMomentumDir ;
308 
309  // Momentum direction, energy and polarisation are unchanged by transport
310  //
311  fTransportEndMomentumDir = startMomentumDir ;
312  fTransportEndKineticEnergy = track.GetKineticEnergy() ;
313  fTransportEndSpin = track.GetPolarization();
314  fParticleIsLooping = false ;
315  fMomentumChanged = false ;
316  fEndGlobalTimeComputed = false ;
317  }
318  else // A field exerts force
319  {
320  G4double momentumMagnitude = pParticle->GetTotalMomentum() ;
321  G4ThreeVector EndUnitMomentum ;
322  G4double lengthAlongCurve ;
323 
324  G4ChargeState chargeState(particleCharge, // The charge can change (dynamic)
325  magneticMoment,
326  pParticleDef->GetPDGSpin() );
327  // For insurance, could set it again
328  // chargeState.SetPDGSpin(pParticleDef->GetPDGSpin() ); // Provisionally in same object
329 
330  G4EquationOfMotion* equationOfMotion =
332  ->GetEquationOfMotion();
333 
334 // equationOfMotion->SetChargeMomentumMass( particleCharge,
335  equationOfMotion->SetChargeMomentumMass( chargeState,
336  momentumMagnitude,
337  restMass);
338 
339  G4FieldTrack aFieldTrack = G4FieldTrack( startPosition,
340  track.GetGlobalTime(), // Lab.
341  // track.GetProperTime(), // Particle rest frame
342  track.GetMomentumDirection(),
343  track.GetKineticEnergy(),
344  restMass,
345  particleCharge,
346  track.GetPolarization(),
347  pParticleDef->GetPDGMagneticMoment(),
348  0.0, // Length along track
349  pParticleDef->GetPDGSpin()
350  ) ;
351 
352  if( currentMinimumStep > 0 )
353  {
354  // Do the Transport in the field (non recti-linear)
355  //
356  lengthAlongCurve = fFieldPropagator->ComputeStep( aFieldTrack,
357  currentMinimumStep,
358  currentSafety,
359  track.GetVolume() ) ;
360 
362  // It is possible that step was reduced in PropagatorInField due to previous zero steps
363  // To cope with case that reduced step is taken in full, we must rely on PiF to obtain this
364  // value.
365 
366  geometryStepLength = std::min( lengthAlongCurve, currentMinimumStep );
367 
368  // Remember last safety origin & value.
369  //
370  fPreviousSftOrigin = startPosition ;
371  fPreviousSafety = currentSafety ;
372  fpSafetyHelper->SetCurrentSafety( currentSafety, startPosition);
373  }
374  else
375  {
376  geometryStepLength = lengthAlongCurve= 0.0 ;
377  fGeometryLimitedStep = false ;
378  }
379 
380  // Get the End-Position and End-Momentum (Dir-ection)
381  //
382  fTransportEndPosition = aFieldTrack.GetPosition() ;
383 
384  // Momentum: Magnitude and direction can be changed too now ...
385  //
386  fMomentumChanged = true ;
387  fTransportEndMomentumDir = aFieldTrack.GetMomentumDir() ;
388 
390 
392  {
393  // If the field can change energy, then the time must be integrated
394  // - so this should have been updated
395  //
397  fEndGlobalTimeComputed = true;
398 
399  // was ( fCandidateEndGlobalTime != track.GetGlobalTime() );
400  // a cleaner way is to have FieldTrack knowing whether time is updated.
401  }
402  else
403  {
404  // The energy should be unchanged by field transport,
405  // - so the time changed will be calculated elsewhere
406  //
407  fEndGlobalTimeComputed = false;
408 
409  // Check that the integration preserved the energy
410  // - and if not correct this!
411  G4double startEnergy= track.GetKineticEnergy();
413 
414  static G4ThreadLocal G4int no_inexact_steps=0, no_large_ediff;
415  G4double absEdiff = std::fabs(startEnergy- endEnergy);
416  if( absEdiff > perMillion * endEnergy )
417  {
418  no_inexact_steps++;
419  // Possible statistics keeping here ...
420  }
421  if( fVerboseLevel > 1 )
422  {
423  if( std::fabs(startEnergy- endEnergy) > perThousand * endEnergy )
424  {
425  static G4ThreadLocal G4int no_warnings= 0, warnModulo=1,
426  moduloFactor= 10;
427  no_large_ediff ++;
428  if( (no_large_ediff% warnModulo) == 0 )
429  {
430  no_warnings++;
431  G4cout << "WARNING - G4Transportation::AlongStepGetPIL() "
432  << " Energy change in Step is above 1^-3 relative value. " << G4endl
433  << " Relative change in 'tracking' step = "
434  << std::setw(15) << (endEnergy-startEnergy)/startEnergy << G4endl
435  << " Starting E= " << std::setw(12) << startEnergy / MeV << " MeV " << G4endl
436  << " Ending E= " << std::setw(12) << endEnergy / MeV << " MeV " << G4endl;
437  G4cout << " Energy has been corrected -- however, review"
438  << " field propagation parameters for accuracy." << G4endl;
439  if( (fVerboseLevel > 2 ) || (no_warnings<4) || (no_large_ediff == warnModulo * moduloFactor) )
440  {
441  G4cout << " These include EpsilonStepMax(/Min) in G4FieldManager "
442  << " which determine fractional error per step for integrated quantities. " << G4endl
443  << " Note also the influence of the permitted number of integration steps."
444  << G4endl;
445  }
446  G4cerr << "ERROR - G4Transportation::AlongStepGetPIL()" << G4endl
447  << " Bad 'endpoint'. Energy change detected"
448  << " and corrected. "
449  << " Has occurred already "
450  << no_large_ediff << " times." << G4endl;
451  if( no_large_ediff == warnModulo * moduloFactor )
452  {
453  warnModulo *= moduloFactor;
454  }
455  }
456  }
457  } // end of if (fVerboseLevel)
458 
459  // Correct the energy for fields that conserve it
460  // This - hides the integration error
461  // - but gives a better physical answer
462  fTransportEndKineticEnergy= track.GetKineticEnergy();
463  }
464 
465  fTransportEndSpin = aFieldTrack.GetSpin();
467  fEndPointDistance = (fTransportEndPosition - startPosition).mag() ;
468  }
469 
470  // If we are asked to go a step length of 0, and we are on a boundary
471  // then a boundary will also limit the step -> we must flag this.
472  //
473  if( currentMinimumStep == 0.0 )
474  {
475  if( currentSafety == 0.0 ) { fGeometryLimitedStep = true; }
476  }
477 
478  // Update the safety starting from the end-point,
479  // if it will become negative at the end-point.
480  //
481  if( currentSafety < fEndPointDistance )
482  {
483  if( particleCharge != 0.0 )
484  {
485  G4double endSafety =
487  currentSafety = endSafety ;
488  fPreviousSftOrigin = fTransportEndPosition ;
489  fPreviousSafety = currentSafety ;
491 
492  // Because the Stepping Manager assumes it is from the start point,
493  // add the StepLength
494  //
495  currentSafety += fEndPointDistance ;
496 
497 #ifdef G4DEBUG_TRANSPORT
498  G4cout.precision(12) ;
499  G4cout << "***G4Transportation::AlongStepGPIL ** " << G4endl ;
500  G4cout << " Called Navigator->ComputeSafety at " << fTransportEndPosition
501  << " and it returned safety= " << endSafety << G4endl ;
502  G4cout << " Adding endpoint distance " << fEndPointDistance
503  << " to obtain pseudo-safety= " << currentSafety << G4endl ;
504  }
505  else
506  {
507  G4cout << "***G4Transportation::AlongStepGPIL ** " << G4endl ;
508  G4cout << " Avoiding call to ComputeSafety : " << G4endl;
509  G4cout << " charge = " << particleCharge << G4endl;
510  G4cout << " mag moment = " << magneticMoment << G4endl;
511 #endif
512  }
513  }
514 
515  fParticleChange.ProposeTrueStepLength(geometryStepLength) ;
516 
517  return geometryStepLength ;
518 }
G4double GetMagneticMoment() const
G4double GetMass() const
static const double MeV
Definition: G4SIunits.hh:211
G4PropagatorInField * fFieldPropagator
virtual void SetChargeMomentumMass(G4ChargeState particleCharge, G4double MomentumXc, G4double MassXc2)=0
G4bool IsParticleLooping() const
G4double fTransportEndKineticEnergy
virtual G4double ComputeStep(const G4ThreeVector &pGlobalPoint, const G4ThreeVector &pDirection, const G4double pCurrentProposedStepLength, G4double &pNewSafety)
Definition: G4Navigator.cc:747
G4bool DoesFieldChangeEnergy() const
G4double GetTotalMomentum() const
static const double perThousand
Definition: G4SIunits.hh:330
G4ThreeVector fTransportEndPosition
#define G4ThreadLocal
Definition: tls.hh:89
G4SafetyHelper * fpSafetyHelper
int G4int
Definition: G4Types.hh:78
double mag2() const
virtual void ConfigureForTrack(const G4Track *)
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
G4bool IsGravityActive() const
Definition: G4Field.hh:98
G4double GetPDGMagneticMoment() const
const G4MagIntegratorStepper * GetStepper() const
G4ThreeVector fPreviousSftOrigin
G4FieldManager * FindAndSetFieldManager(G4VPhysicalVolume *pCurrentPhysVol)
static const double perMillion
Definition: G4SIunits.hh:331
G4double ComputeStep(G4FieldTrack &pFieldTrack, G4double pCurrentProposedStepLength, G4double &pNewSafety, G4VPhysicalVolume *pPhysVol=0)
G4double GetCharge() const
const G4ThreeVector & GetMomentumDirection() const
G4ChordFinder * GetChordFinder()
G4double GetLabTimeOfFlight() const
G4FieldManager * GetCurrentFieldManager()
const G4ThreeVector & GetMomentumDir() const
G4ThreeVector GetSpin() const
void SetCurrentSafety(G4double val, const G4ThreeVector &pos)
G4ParticleDefinition * GetDefinition() const
#define G4endl
Definition: G4ios.hh:61
G4ThreeVector fTransportEndSpin
G4double GetKineticEnergy() const
T sqr(const T &x)
Definition: templates.hh:145
G4Navigator * fLinearNavigator
G4ThreeVector GetPosition() const
double G4double
Definition: G4Types.hh:76
virtual G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4double pProposedMaxLength=DBL_MAX, const G4bool keepState=true)
G4MagInt_Driver * GetIntegrationDriver()
G4double fCandidateEndGlobalTime
G4ParticleChangeForTransport fParticleChange
const G4Field * GetDetectorField() const
G4ThreeVector fTransportEndMomentumDir
G4GLOB_DLL std::ostream G4cerr
static G4bool fUseMagneticMoment
Here is the call graph for this function:
Here is the caller graph for this function:

◆ AtRestDoIt()

G4VParticleChange* G4Transportation::AtRestDoIt ( const G4Track &  ,
const G4Step &   
)
inlinevirtual

Implements G4VProcess.

Definition at line 138 of file G4Transportation.hh.

141  {return 0;};
Here is the call graph for this function:

◆ AtRestGetPhysicalInteractionLength()

G4double G4Transportation::AtRestGetPhysicalInteractionLength ( const G4Track &  ,
G4ForceCondition *   
)
inlinevirtual

Implements G4VProcess.

Definition at line 132 of file G4Transportation.hh.

135  { return -1.0; };

◆ DoesGlobalFieldExist()

G4bool G4Transportation::DoesGlobalFieldExist ( )
protected
Here is the caller graph for this function:

◆ EnableShortStepOptimisation()

void G4Transportation::EnableShortStepOptimisation ( G4bool  optimise = true)
inline
Here is the caller graph for this function:

◆ EnableUseMagneticMoment()

G4bool G4Transportation::EnableUseMagneticMoment ( G4bool  useMoment = true)
static

Definition at line 821 of file G4Transportation.cc.

822 {
823  G4bool lastValue= fUseMagneticMoment;
824  fUseMagneticMoment= useMoment;
826  return lastValue;
827 }
bool G4bool
Definition: G4Types.hh:79
static G4bool fUseMagneticMoment
Here is the caller graph for this function:

◆ FieldExertedForce()

G4bool G4Transportation::FieldExertedForce ( )
inline

Definition at line 95 of file G4Transportation.hh.

95 { return fFieldExertedForce; }
Here is the call graph for this function:

◆ GetMaxEnergyKilled()

G4double G4Transportation::GetMaxEnergyKilled ( ) const
inline
Here is the caller graph for this function:

◆ GetPropagatorInField()

G4PropagatorInField* G4Transportation::GetPropagatorInField ( )
Here is the caller graph for this function:

◆ GetSumEnergyKilled()

G4double G4Transportation::GetSumEnergyKilled ( ) const
inline
Here is the caller graph for this function:

◆ GetThresholdImportantEnergy()

G4double G4Transportation::GetThresholdImportantEnergy ( ) const
inline
Here is the caller graph for this function:

◆ GetThresholdTrials()

G4int G4Transportation::GetThresholdTrials ( ) const
inline
Here is the caller graph for this function:

◆ GetThresholdWarningEnergy()

G4double G4Transportation::GetThresholdWarningEnergy ( ) const
inline
Here is the caller graph for this function:

◆ GetVerboseLevel()

G4int G4Transportation::GetVerboseLevel ( ) const
inline
Here is the caller graph for this function:

◆ PostStepDoIt()

G4VParticleChange * G4Transportation::PostStepDoIt ( const G4Track &  track,
const G4Step &  stepData 
)
virtual

Implements G4VProcess.

Definition at line 661 of file G4Transportation.cc.

663 {
664  G4TouchableHandle retCurrentTouchable ; // The one to return
665  G4bool isLastStep= false;
666 
667  // Initialize ParticleChange (by setting all its members equal
668  // to corresponding members in G4Track)
669  // fParticleChange.Initialize(track) ; // To initialise TouchableChange
670 
671  fParticleChange.ProposeTrackStatus(track.GetTrackStatus()) ;
672 
673  // If the Step was determined by the volume boundary,
674  // logically relocate the particle
675 
677  {
678  // fCurrentTouchable will now become the previous touchable,
679  // and what was the previous will be freed.
680  // (Needed because the preStepPoint can point to the previous touchable)
681 
684  LocateGlobalPointAndUpdateTouchableHandle( track.GetPosition(),
685  track.GetMomentumDirection(),
687  true ) ;
688  // Check whether the particle is out of the world volume
689  // If so it has exited and must be killed.
690  //
691  if( fCurrentTouchableHandle->GetVolume() == 0 )
692  {
693  fParticleChange.ProposeTrackStatus( fStopAndKill ) ;
694  }
695  retCurrentTouchable = fCurrentTouchableHandle ;
696  fParticleChange.SetTouchableHandle( fCurrentTouchableHandle ) ;
697 
698  // Update the Step flag which identifies the Last Step in a volume
699  if( !fFieldExertedForce )
700  isLastStep = fLinearNavigator->ExitedMotherVolume()
702  else
703  isLastStep = fFieldPropagator->IsLastStepInVolume();
704  }
705  else // fGeometryLimitedStep is false
706  {
707  // This serves only to move the Navigator's location
708  //
710 
711  // The value of the track's current Touchable is retained.
712  // (and it must be correct because we must use it below to
713  // overwrite the (unset) one in particle change)
714  // It must be fCurrentTouchable too ??
715  //
716  fParticleChange.SetTouchableHandle( track.GetTouchableHandle() ) ;
717  retCurrentTouchable = track.GetTouchableHandle() ;
718 
719  isLastStep= false;
720  } // endif ( fGeometryLimitedStep )
721  fLastStepInVolume= isLastStep;
722 
723  fParticleChange.ProposeFirstStepInVolume(fFirstStepInVolume);
724  fParticleChange.ProposeLastStepInVolume(isLastStep);
725 
726  const G4VPhysicalVolume* pNewVol = retCurrentTouchable->GetVolume() ;
727  const G4Material* pNewMaterial = 0 ;
728  const G4VSensitiveDetector* pNewSensitiveDetector = 0 ;
729 
730  if( pNewVol != 0 )
731  {
732  pNewMaterial= pNewVol->GetLogicalVolume()->GetMaterial();
733  pNewSensitiveDetector= pNewVol->GetLogicalVolume()->GetSensitiveDetector();
734  }
735 
736  // ( <const_cast> pNewMaterial ) ;
737  // ( <const_cast> pNewSensitiveDetector) ;
738 
739  fParticleChange.SetMaterialInTouchable( (G4Material *) pNewMaterial ) ;
740  fParticleChange.SetSensitiveDetectorInTouchable( (G4VSensitiveDetector *) pNewSensitiveDetector ) ;
741 
742  const G4MaterialCutsCouple* pNewMaterialCutsCouple = 0;
743  if( pNewVol != 0 )
744  {
745  pNewMaterialCutsCouple=pNewVol->GetLogicalVolume()->GetMaterialCutsCouple();
746  }
747 
748  if( pNewVol!=0 && pNewMaterialCutsCouple!=0 && pNewMaterialCutsCouple->GetMaterial()!=pNewMaterial )
749  {
750  // for parametrized volume
751  //
752  pNewMaterialCutsCouple =
754  ->GetMaterialCutsCouple(pNewMaterial,
755  pNewMaterialCutsCouple->GetProductionCuts());
756  }
757  fParticleChange.SetMaterialCutsCoupleInTouchable( pNewMaterialCutsCouple );
758 
759  // temporarily until Get/Set Material of ParticleChange,
760  // and StepPoint can be made const.
761  // Set the touchable in ParticleChange
762  // this must always be done because the particle change always
763  // uses this value to overwrite the current touchable pointer.
764  //
765  fParticleChange.SetTouchableHandle(retCurrentTouchable) ;
766 
767  return &fParticleChange ;
768 }
G4bool EnteredDaughterVolume() const
G4bool ExitedMotherVolume() const
G4PropagatorInField * fFieldPropagator
const G4Material * GetMaterial() const
G4Material * GetMaterial() const
G4ProductionCuts * GetProductionCuts() const
void SetGeometricallyLimitedStep()
bool G4bool
Definition: G4Types.hh:79
G4TouchableHandle fCurrentTouchableHandle
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:44
static G4ProductionCutsTable * GetProductionCutsTable()
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
G4VSensitiveDetector * GetSensitiveDetector() const
G4Navigator * fLinearNavigator
G4LogicalVolume * GetLogicalVolume() const
G4ParticleChangeForTransport fParticleChange
virtual void LocateGlobalPointWithinVolume(const G4ThreeVector &position)
Definition: G4Navigator.cc:582
Here is the call graph for this function:

◆ PostStepGetPhysicalInteractionLength()

G4double G4Transportation::PostStepGetPhysicalInteractionLength ( const G4Track &  ,
G4double  previousStepSize,
G4ForceCondition *  pForceCond 
)
virtual

Implements G4VProcess.

Definition at line 649 of file G4Transportation.cc.

652 {
653  fFieldExertedForce = false; // Not known
654  *pForceCond = Forced ;
655  return DBL_MAX ; // was kInfinity ; but convention now is DBL_MAX
656 }
#define DBL_MAX
Definition: templates.hh:83
Here is the caller graph for this function:

◆ ResetKilledStatistics()

void G4Transportation::ResetKilledStatistics ( G4int  report = 1)
inline
Here is the caller graph for this function:

◆ SetPropagatorInField()

void G4Transportation::SetPropagatorInField ( G4PropagatorInField pFieldPropagator)
Here is the caller graph for this function:

◆ SetThresholdImportantEnergy()

void G4Transportation::SetThresholdImportantEnergy ( G4double  newEnImp)
inline
Here is the caller graph for this function:

◆ SetThresholdTrials()

void G4Transportation::SetThresholdTrials ( G4int  newMaxTrials)
inline
Here is the caller graph for this function:

◆ SetThresholdWarningEnergy()

void G4Transportation::SetThresholdWarningEnergy ( G4double  newEnWarn)
inline
Here is the caller graph for this function:

◆ SetVerboseLevel()

void G4Transportation::SetVerboseLevel ( G4int  verboseLevel)
inline
Here is the caller graph for this function:

◆ StartTracking()

void G4Transportation::StartTracking ( G4Track *  aTrack)
virtual

Reimplemented from G4VProcess.

Definition at line 774 of file G4Transportation.cc.

775 {
777  fNewTrack= true;
778  fFirstStepInVolume= true;
779  fLastStepInVolume= false;
780 
781  // The actions here are those that were taken in AlongStepGPIL
782  // when track.GetCurrentStepNumber()==1
783 
784  // reset safety value and center
785  //
786  fPreviousSafety = 0.0 ;
787  fPreviousSftOrigin = G4ThreeVector(0.,0.,0.) ;
788 
789  // reset looping counter -- for motion in field
790  fNoLooperTrials= 0;
791  // Must clear this state .. else it depends on last track's value
792  // --> a better solution would set this from state of suspended track TODO ?
793  // Was if( aTrack->GetCurrentStepNumber()==1 ) { .. }
794 
795  // ChordFinder reset internal state
796  //
797  if( DoesGlobalFieldExist() )
798  {
800  // Resets all state of field propagator class (ONLY)
801  // including safety values (in case of overlaps and to wipe for first track).
802 
803  // G4ChordFinder* chordF= fFieldPropagator->GetChordFinder();
804  // if( chordF ) chordF->ResetStepEstimate();
805  }
806 
807  // Make sure to clear the chord finders of all fields (ie managers)
808  //
810  fieldMgrStore->ClearAllChordFindersState();
811 
812  // Update the current touchable handle (from the track's)
813  //
814  fCurrentTouchableHandle = aTrack->GetTouchableHandle();
815 
816  // Inform field propagator of new track
818 }
G4PropagatorInField * fFieldPropagator
CLHEP::Hep3Vector G4ThreeVector
virtual void StartTracking(G4Track *)
Definition: G4VProcess.cc:101
G4bool DoesGlobalFieldExist()
G4ThreeVector fPreviousSftOrigin
G4TouchableHandle fCurrentTouchableHandle
static G4FieldManagerStore * GetInstance()
Here is the call graph for this function:
Here is the caller graph for this function:

Friends And Related Function Documentation

◆ G4CoupledTransportation

friend class G4CoupledTransportation
friend

Definition at line 219 of file G4Transportation.hh.

Member Data Documentation

◆ fCandidateEndGlobalTime

G4double G4Transportation::fCandidateEndGlobalTime
private

Definition at line 167 of file G4Transportation.hh.

◆ fCurrentTouchableHandle

G4TouchableHandle G4Transportation::fCurrentTouchableHandle
private

Definition at line 179 of file G4Transportation.hh.

◆ fEndGlobalTimeComputed

G4bool G4Transportation::fEndGlobalTimeComputed
private

Definition at line 166 of file G4Transportation.hh.

◆ fEndPointDistance

G4double G4Transportation::fEndPointDistance
private

Definition at line 188 of file G4Transportation.hh.

◆ fFieldExertedForce

G4bool G4Transportation::fFieldExertedForce
private

Definition at line 177 of file G4Transportation.hh.

◆ fFieldPropagator

G4PropagatorInField* G4Transportation::fFieldPropagator
private

Definition at line 155 of file G4Transportation.hh.

◆ fFirstStepInVolume

G4bool G4Transportation::fFirstStepInVolume
private

Definition at line 172 of file G4Transportation.hh.

◆ fGeometryLimitedStep

G4bool G4Transportation::fGeometryLimitedStep
private

Definition at line 175 of file G4Transportation.hh.

◆ fLastStepInVolume

G4bool G4Transportation::fLastStepInVolume
private

Definition at line 173 of file G4Transportation.hh.

◆ fLinearNavigator

G4Navigator* G4Transportation::fLinearNavigator
private

Definition at line 154 of file G4Transportation.hh.

◆ fMaxEnergyKilled

G4double G4Transportation::fMaxEnergyKilled
private

Definition at line 203 of file G4Transportation.hh.

◆ fMomentumChanged

G4bool G4Transportation::fMomentumChanged
private

Definition at line 165 of file G4Transportation.hh.

◆ fNewTrack

G4bool G4Transportation::fNewTrack
private

Definition at line 171 of file G4Transportation.hh.

◆ fNoLooperTrials

G4int G4Transportation::fNoLooperTrials
private

Definition at line 200 of file G4Transportation.hh.

◆ fParticleChange

G4ParticleChangeForTransport G4Transportation::fParticleChange
private

Definition at line 185 of file G4Transportation.hh.

◆ fParticleIsLooping

G4bool G4Transportation::fParticleIsLooping
private

Definition at line 170 of file G4Transportation.hh.

◆ fPreviousSafety

G4double G4Transportation::fPreviousSafety
private

Definition at line 182 of file G4Transportation.hh.

◆ fPreviousSftOrigin

G4ThreeVector G4Transportation::fPreviousSftOrigin
private

Definition at line 181 of file G4Transportation.hh.

◆ fpSafetyHelper

G4SafetyHelper* G4Transportation::fpSafetyHelper
private

Definition at line 209 of file G4Transportation.hh.

◆ fShortStepOptimisation

G4bool G4Transportation::fShortStepOptimisation
private

Definition at line 207 of file G4Transportation.hh.

◆ fSumEnergyKilled

G4double G4Transportation::fSumEnergyKilled
private

Definition at line 202 of file G4Transportation.hh.

◆ fThreshold_Important_Energy

G4double G4Transportation::fThreshold_Important_Energy
private

Definition at line 193 of file G4Transportation.hh.

◆ fThreshold_Warning_Energy

G4double G4Transportation::fThreshold_Warning_Energy
private

Definition at line 192 of file G4Transportation.hh.

◆ fThresholdTrials

G4int G4Transportation::fThresholdTrials
private

Definition at line 194 of file G4Transportation.hh.

◆ fTransportEndKineticEnergy

G4double G4Transportation::fTransportEndKineticEnergy
private

Definition at line 163 of file G4Transportation.hh.

◆ fTransportEndMomentumDir

G4ThreeVector G4Transportation::fTransportEndMomentumDir
private

Definition at line 162 of file G4Transportation.hh.

◆ fTransportEndPosition

G4ThreeVector G4Transportation::fTransportEndPosition
private

Definition at line 161 of file G4Transportation.hh.

◆ fTransportEndSpin

G4ThreeVector G4Transportation::fTransportEndSpin
private

Definition at line 164 of file G4Transportation.hh.

◆ fUseMagneticMoment

G4bool G4Transportation::fUseMagneticMoment =false
staticprivate

Definition at line 220 of file G4Transportation.hh.

◆ fVerboseLevel

G4int G4Transportation::fVerboseLevel
private

Definition at line 212 of file G4Transportation.hh.


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