Belle II Software  release-08-01-10
G4LongLivedNeutralTransportation Class Reference

Concrete class that does the geometrical transport. More...

#include <G4LongLivedNeutralTransportation.h>

Inheritance diagram for G4LongLivedNeutralTransportation:
Collaboration diagram for G4LongLivedNeutralTransportation:

Public Member Functions

 G4LongLivedNeutralTransportation (G4int verbosityLevel=1)
 Constructor. More...
 
 ~G4LongLivedNeutralTransportation ()
 Destructor.
 
G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
 G4VProcess::AlongStepGetPhysicalInteractionLength() implementation,. More...
 
G4VParticleChange * AlongStepDoIt (const G4Track &track, const G4Step &stepData)
 G4VProcess::AlongStepDoIt() implementation,. More...
 
G4VParticleChange * PostStepDoIt (const G4Track &track, const G4Step &stepData)
 G4VProcess::PostStepDoIt() implementation,. More...
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &, G4double previousStepSize, G4ForceCondition *pForceCond)
 G4VProcess::PostStepGetPhysicalInteractionLength() implementation. More...
 
G4bool FieldExertedForce ()
 References fFieldExertedForce.
 
G4PropagatorInField * GetPropagatorInField ()
 Access fFieldPropagator, the assistant class that Propagate in a Field.
 
void SetPropagatorInField (G4PropagatorInField *pFieldPropagator)
 Access/set the assistant class that Propagate in a Field.
 
G4double GetThresholdWarningEnergy () const
 Access fThreshold_Warning_Energy.
 
G4double GetThresholdImportantEnergy () const
 Access fThreshold_Important_Energy.
 
G4int GetThresholdTrials () const
 Access fThresholdTrials.
 
void SetThresholdWarningEnergy (G4double newEnWarn)
 Set fThreshold_Warning_Energy.
 
void SetThresholdImportantEnergy (G4double newEnImp)
 Set fThreshold_Important_Energy.
 
void SetThresholdTrials (G4int newMaxTrials)
 Set fThresholdTrials.
 
void SetHighLooperThresholds ()
 Get/Set parameters for killing loopers: Above 'important' energy a 'looping' particle in field will NOT be abandoned, except after fThresholdTrials attempts. More...
 
void SetLowLooperThresholds ()
 Shortcut method - old values (meant for HEP)
 
void ReportLooperThresholds ()
 Set low thresholds - for low-E applications.
 
G4double GetMaxEnergyKilled () const
 Print values of looper thresholds. More...
 
G4double GetSumEnergyKilled () const
 Access fSumEnergyKilled.
 
void ResetKilledStatistics (G4int report=1)
 Statistics for tracks killed (currently due to looping in field)
 
void EnableShortStepOptimisation (G4bool optimise=true)
 Whether short steps < safety will avoid to call Navigator (if field=0)
 
G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 No operation in AtRestGPIL.
 
G4VParticleChange * AtRestDoIt (const G4Track &, const G4Step &)
 No operation in AtRestDoIt.
 
void StartTracking (G4Track *aTrack)
 Reset state for new (potentially resumed) track.
 
virtual void ProcessDescription (std::ostream &outFile) const
 G4LongLivedNeutralTransportation::ProcessDescription() More...
 
void PrintStatistics (std::ostream &outStr) const
 returns current logging info of the algorithm
 

Static Public Member Functions

static void SetSilenceLooperWarnings (G4bool val)
 Do not warn about 'looping' particles.
 
static G4bool GetSilenceLooperWarnings ()
 Do not throw exception about 'looping' particles.
 

Protected Member Functions

G4bool DoesGlobalFieldExist ()
 Checks whether a field exists for the "global" field manager.
 

Private Attributes

G4Navigator * fLinearNavigator
 Propagator used to transport the particle.
 
G4PropagatorInField * fFieldPropagator
 Propagator used to transport the particle.
 
G4ThreeVector fTransportEndPosition = G4ThreeVector(0.0, 0.0, 0.0)
 The particle's state after this Step, Store for DoIt.
 
G4ThreeVector fTransportEndMomentumDir = G4ThreeVector(0.0, 0.0, 0.0)
 The particle's state after this Step, Store for DoIt.
 
G4double fTransportEndKineticEnergy = 0.0
 The particle's state after this Step, Store for DoIt.
 
G4ThreeVector fTransportEndSpin = G4ThreeVector(0.0, 0.0, 0.0)
 The particle's state after this Step, Store for DoIt.
 
G4bool fMomentumChanged = true
 The particle's state after this Step, Store for DoIt.
 
G4bool fEndGlobalTimeComputed = false
 The particle's state after this Step, Store for DoIt.
 
G4double fCandidateEndGlobalTime = 0.0
 The particle's state after this Step, Store for DoIt.
 
G4bool fAnyFieldExists = false
 Flag for existing fields.
 
G4bool fNewTrack = true
 Flag from StartTracking.
 
G4bool fFirstStepInVolume = true
 Flag first step in a geom. More...
 
G4bool fLastStepInVolume = false
 Flag last step in a geom. More...
 
G4bool fGeometryLimitedStep = true
 Flag to determine whether a boundary was reached.
 
G4bool fFieldExertedForce = false
 During current step.
 
G4TouchableHandle fCurrentTouchableHandle
 Current touchable handle.
 
G4ThreeVector fPreviousSftOrigin
 Remember last safety origin.
 
G4double fPreviousSafety
 Remember last safety value.
 
G4ParticleChangeForTransport fParticleChange
 New ParticleChange.
 
G4double fEndPointDistance
 Endpoint distance.
 
G4double fThreshold_Warning_Energy = 1.0 * CLHEP::keV
 Warn above this energy about looping particle.
 
G4double fThreshold_Important_Energy = 1.0 * CLHEP::MeV
 Give a few trial above this E for looping particle.
 
G4int fThresholdTrials = 10
 Number of trials an important looper survives.
 
G4int fAbandonUnstableTrials = 0
 Number of trials after which to abandon.
 
G4int fNoLooperTrials = 0
 Counter for steps in which particle reports 'looping', if it is above 'Important' Energy.
 
G4double fSumEnergyKilled = 0.0
 Sum of abandoned looping tracks energies.
 
G4bool fShortStepOptimisation
 Whether to avoid calling G4Navigator for short step ( < safety) If using it, the safety estimate for endpoint will likely be smaller.
 
G4SafetyHelper * fpSafetyHelper
 To pass it the safety value obtained.
 

Static Private Attributes

static G4bool fUseMagneticMoment
 Flag take into account magnetic moment.
 
static G4bool fUseGravity
 Flag take into account gravity.
 
static G4bool fSilenceLooperWarnings = false
 Flag to Supress all 'looper' warnings.
 

Friends

class G4CoupledTransportation
 

Detailed Description

Concrete class that does the geometrical transport.

Definition at line 36 of file G4LongLivedNeutralTransportation.h.

Constructor & Destructor Documentation

◆ G4LongLivedNeutralTransportation()

G4LongLivedNeutralTransportation ( G4int  verbosityLevel = 1)

Constructor.

Parameters
verbosityLevel

Definition at line 33 of file G4LongLivedNeutralTransportation.cc.

34  : G4VProcess(G4String("LongLivedNeutralTransportation"), fTransportation),
35  fFieldExertedForce(false),
36  fPreviousSftOrigin(0., 0., 0.),
37  fPreviousSafety(0.0),
38  fEndPointDistance(-1.0),
39  fShortStepOptimisation(false) // Old default: true (=fast short steps)
40 {
41  SetProcessSubType(static_cast<G4int>(TRANSPORTATION));
42  pParticleChange = &fParticleChange; // Required to conform to G4VProcess
43  SetVerboseLevel(verbosity);
44 
45  G4TransportationManager* transportMgr ;
46 
47  transportMgr = G4TransportationManager::GetTransportationManager() ;
48 
49  fLinearNavigator = transportMgr->GetNavigatorForTracking() ;
50 
51  fFieldPropagator = transportMgr->GetPropagatorInField() ;
52 
53 
54 
55 
57  // Use the old defaults: Warning = 100 MeV, Important = 250 MeV, No Trials = 10;
58 
59  // Cannot determine whether a field exists here, as it would
60  // depend on the relative order of creating the detector's
61  // field and this process. That order is not guaranteed.
63  // This value must be updated using DoesAnyFieldExist() at least at the
64  // start of each Run -- for now this is at the Start of every Track. TODO
65 
66  static G4ThreadLocal G4TouchableHandle* pNullTouchableHandle = 0;
67  pNullTouchableHandle = new G4TouchableHandle;
68  fCurrentTouchableHandle = *pNullTouchableHandle;
69  // Points to (G4VTouchable*) 0
70 
71 
72 #ifdef G4VERBOSE
73  if (verboseLevel > 0) {
74  G4cout << " G4LongLivedNeutralTransportation constructor> set fShortStepOptimisation to ";
75  if (fShortStepOptimisation) { G4cout << "true" << G4endl; }
76  else { G4cout << "false" << G4endl; }
77  }
78 #endif
79 }
G4TouchableHandle fCurrentTouchableHandle
Current touchable handle.
void SetHighLooperThresholds()
Get/Set parameters for killing loopers: Above 'important' energy a 'looping' particle in field will N...
G4ThreeVector fPreviousSftOrigin
Remember last safety origin.
G4Navigator * fLinearNavigator
Propagator used to transport the particle.
G4bool fShortStepOptimisation
Whether to avoid calling G4Navigator for short step ( < safety) If using it, the safety estimate for ...
G4ParticleChangeForTransport fParticleChange
New ParticleChange.
G4bool DoesGlobalFieldExist()
Checks whether a field exists for the "global" field manager.
G4PropagatorInField * fFieldPropagator
Propagator used to transport the particle.

Member Function Documentation

◆ AlongStepDoIt()

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

G4VProcess::AlongStepDoIt() implementation,.

Parameters
trackPropagating particle track reference
stepDataCurrent step reference

Definition at line 259 of file G4LongLivedNeutralTransportation.cc.

◆ AlongStepGetPhysicalInteractionLength()

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

G4VProcess::AlongStepGetPhysicalInteractionLength() implementation,.

Parameters
trackPropagating particle track reference
previousStepSizeThis argument of base function is ignored
currentMinimumStepCurrent minimum step size
currentSafetyReference to current step safety
selectionPointer for return value of GPILSelection, which is set to default value of CandidateForSelection
Returns
Next geometry step length

Definition at line 108 of file G4LongLivedNeutralTransportation.cc.

◆ GetMaxEnergyKilled()

G4double GetMaxEnergyKilled ( ) const
inline

Print values of looper thresholds.

Access fMaxEnergyKilled

◆ PostStepDoIt()

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

G4VProcess::PostStepDoIt() implementation,.

Parameters
track
stepData

Definition at line 346 of file G4LongLivedNeutralTransportation.cc.

◆ PostStepGetPhysicalInteractionLength()

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

G4VProcess::PostStepGetPhysicalInteractionLength() implementation.

Parameters
trackThis argument of base function is ignored
previousStepSizeThis argument of base function is ignored
pForceCondForce condition by default

Forces the PostStepDoIt action to be called, but does not limit the step

Definition at line 333 of file G4LongLivedNeutralTransportation.cc.

◆ ProcessDescription()

void ProcessDescription ( std::ostream &  outFile) const
virtual

G4LongLivedNeutralTransportation::ProcessDescription()

@outfile Description of process

Definition at line 558 of file G4LongLivedNeutralTransportation.cc.

◆ SetHighLooperThresholds()

void SetHighLooperThresholds ( )

Get/Set parameters for killing loopers: Above 'important' energy a 'looping' particle in field will NOT be abandoned, except after fThresholdTrials attempts.

Below Warning energy, no verbosity for looping particles is issued

Definition at line 523 of file G4LongLivedNeutralTransportation.cc.

Member Data Documentation

◆ fFirstStepInVolume

G4bool fFirstStepInVolume = true
private

Flag first step in a geom.

volume

Definition at line 224 of file G4LongLivedNeutralTransportation.h.

◆ fLastStepInVolume

G4bool fLastStepInVolume = false
private

Flag last step in a geom.

volume

Definition at line 225 of file G4LongLivedNeutralTransportation.h.


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