Belle II Software development
G4LongLivedNeutralTransportation Class Reference

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

#include <G4LongLivedNeutralTransportation.h>

Inheritance diagram for G4LongLivedNeutralTransportation:

Public Member Functions

 G4LongLivedNeutralTransportation (G4int verbosityLevel=1)
 Constructor.
 
 ~G4LongLivedNeutralTransportation ()
 Destructor.
 
G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
 G4VProcess::AlongStepGetPhysicalInteractionLength() implementation,.
 
G4VParticleChange * AlongStepDoIt (const G4Track &track, const G4Step &stepData)
 G4VProcess::AlongStepDoIt() implementation,.
 
G4VParticleChange * PostStepDoIt (const G4Track &track, const G4Step &stepData)
 G4VProcess::PostStepDoIt() implementation,.
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *pForceCond)
 G4VProcess::PostStepGetPhysicalInteractionLength() implementation.
 
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.
 
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.
 
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 &outStream) const
 G4LongLivedNeutralTransportation::ProcessDescription()
 
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.
 
G4bool fLastStepInVolume = false
 Flag last step in a geom.
 
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),
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.

◆ ~G4LongLivedNeutralTransportation()

Destructor.

Definition at line 83 of file G4LongLivedNeutralTransportation.cc.

84{
85 if (fSumEnergyKilled > 0.0) {
86 PrintStatistics(G4cout);
87 }
88
89}
void PrintStatistics(std::ostream &outStr) const
returns current logging info of the algorithm
G4double fSumEnergyKilled
Sum of abandoned looping tracks energies.

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.

261{
262 // static G4ThreadLocal G4long noCallsASDI = 0;
263 // noCallsASDI++;
264
265 fParticleChange.Initialize(track) ;
266
267 // Code for specific process
268 //
269 fParticleChange.ProposePosition(fTransportEndPosition) ;
270 fParticleChange.ProposeMomentumDirection(fTransportEndMomentumDir) ;
272 fParticleChange.SetMomentumChanged(fMomentumChanged) ;
273
274 fParticleChange.ProposePolarization(fTransportEndSpin);
275
276 G4double deltaTime = 0.0 ;
277
278 // Calculate Lab Time of Flight (ONLY if field Equations used it!)
279 // G4double endTime = fCandidateEndGlobalTime;
280 // G4double delta_time = endTime - startTime;
281
282 G4double startTime = track.GetGlobalTime() ;
283
284 const G4DynamicParticle* pParticle = track.GetDynamicParticle() ;
285 const G4PrimaryParticle* primaryParticle = pParticle->GetPrimaryParticle();
286
287
289 // The time was not integrated .. make the best estimate possible
290 //
291
292 // use c*beta from primary particle, llp will not be slowed down anyways
293 // default was: stepData.GetPreStepPoint()->GetVelocity();, but assumes wrong values from particle definition
294 B2DEBUG(0, "Velocity calculated from Step data: " << stepData.GetPreStepPoint()->GetVelocity());
295
296 G4double initialVelocity = c_light * primaryParticle->GetTotalMomentum() / primaryParticle->GetTotalEnergy();
297 G4double stepLength = track.GetStepLength();
298
299 deltaTime = 0.0; // in case initialVelocity = 0
300 if (initialVelocity > 0.0) { deltaTime = stepLength / initialVelocity; }
301
302 fCandidateEndGlobalTime = startTime + deltaTime ;
303 fParticleChange.ProposeLocalTime(track.GetLocalTime() + deltaTime) ;
304 } else {
305 deltaTime = fCandidateEndGlobalTime - startTime ;
306 fParticleChange.ProposeGlobalTime(fCandidateEndGlobalTime) ;
307 }
308
309 // Now Correct by Lorentz factor to get delta "proper" Time
310 G4double restMass = primaryParticle->GetMass() ;
311 G4double deltaProperTime = deltaTime * (restMass / primaryParticle->GetTotalEnergy()) ;
312
313 fParticleChange.ProposeProperTime(track.GetProperTime() + deltaProperTime) ;
314 //fParticleChange.ProposeTrueStepLength( track.GetStepLength() ) ;
315
316 // Another (sometimes better way) is to use a user-limit maximum Step size
317 // to alleviate this problem ..
318
319 // Introduce smooth curved trajectories to particle-change
320 //
321 fParticleChange.SetPointerToVectorOfAuxiliaryPoints
322 (fFieldPropagator->GimmeTrajectoryVectorAndForgetIt());
323
324 return &fParticleChange ;
325}
G4double fCandidateEndGlobalTime
The particle's state after this Step, Store for DoIt.
G4ThreeVector fTransportEndMomentumDir
The particle's state after this Step, Store for DoIt.
G4ThreeVector fTransportEndSpin
The particle's state after this Step, Store for DoIt.
G4ThreeVector fTransportEndPosition
The particle's state after this Step, Store for DoIt.
G4bool fEndGlobalTimeComputed
The particle's state after this Step, Store for DoIt.
G4bool fMomentumChanged
The particle's state after this Step, Store for DoIt.
G4double fTransportEndKineticEnergy
The particle's state after this Step, Store for DoIt.

◆ 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.

114{
115 G4double geometryStepLength = -1.0, newSafety = -1.0;
116
117 // Initial actions moved to StartTrack()
118 // --------------------------------------
119 // Note: in case another process changes touchable handle
120 // it will be necessary to add here (for all steps)
121 // fCurrentTouchableHandle = aTrack->GetTouchableHandle();
122
123 // GPILSelection is set to default value of CandidateForSelection
124 // It is a return value
125 //
126 *selection = CandidateForSelection ;
127
129 fLastStepInVolume = false;
130 fNewTrack = false;
131
132 fParticleChange.ProposeFirstStepInVolume(fFirstStepInVolume);
133
134 // Get initial Energy/Momentum of the track
135 //
136 const G4DynamicParticle* pParticle = track.GetDynamicParticle() ;
137 const G4PrimaryParticle* primaryParticle = pParticle->GetPrimaryParticle();
138 G4ThreeVector startMomentumDir = primaryParticle->GetMomentumDirection() ;
139 G4ThreeVector startPosition = track.GetPosition() ;
140
141 // The Step Point safety can be limited by other geometries and/or the
142 // assumptions of any process - it's not always the geometrical safety.
143 // We calculate the starting point's isotropic safety here.
144 //
145 G4ThreeVector OriginShift = startPosition - fPreviousSftOrigin ;
146 G4double MagSqShift = OriginShift.mag2() ;
147 if (MagSqShift >= sqr(fPreviousSafety)) {
148 currentSafety = 0.0 ;
149 } else {
150 currentSafety = fPreviousSafety - std::sqrt(MagSqShift) ;
151 }
152
153 // Is the particle charged or has it a magnetic moment?
154 //
155 G4double particleCharge = 0 ;
156
157 fGeometryLimitedStep = false ;
158
159 // There is no need to locate the current volume. It is Done elsewhere:
160 // On track construction
161 // By the tracking, after all AlongStepDoIts, in "Relocation"
162
163 G4double linearStepLength ;
164 if (fShortStepOptimisation && (currentMinimumStep <= currentSafety)) {
165 // The Step is guaranteed to be taken
166 //
167 geometryStepLength = currentMinimumStep ;
168 fGeometryLimitedStep = false ;
169 } else {
170 // Find whether the straight path intersects a volume
171 //
172 linearStepLength = fLinearNavigator->ComputeStep(startPosition,
173 startMomentumDir,
174 currentMinimumStep,
175 newSafety) ;
176 // Remember last safety origin & value.
177 //
178 fPreviousSftOrigin = startPosition ;
179 fPreviousSafety = newSafety ;
180
181
182 currentSafety = newSafety ;
183
184 fGeometryLimitedStep = (linearStepLength <= currentMinimumStep);
186 // The geometry limits the Step size (an intersection was found.)
187 geometryStepLength = linearStepLength ;
188 } else {
189 // The full Step is taken.
190 geometryStepLength = currentMinimumStep ;
191 }
192
193 fEndPointDistance = geometryStepLength ;
194
195 // Calculate final position
196 //
197 fTransportEndPosition = startPosition + geometryStepLength * startMomentumDir;
198
199 // Momentum direction, energy and polarisation are unchanged by transport
200 //
201 fTransportEndMomentumDir = startMomentumDir ;
202 fTransportEndKineticEnergy = primaryParticle->GetKineticEnergy() ;
203 fTransportEndSpin = track.GetPolarization();
204 fMomentumChanged = false ;
205 fEndGlobalTimeComputed = false ;
206 }
207
208
209 // If we are asked to go a step length of 0, and we are on a boundary
210 // then a boundary will also limit the step -> we must flag this.
211 //
212 if (currentMinimumStep == 0.0) {
213 if (currentSafety == 0.0) { fGeometryLimitedStep = true; }
214 }
215
216 // Update the safety starting from the end-point,
217 // if it will become negative at the end-point.
218 //
219 if (currentSafety < fEndPointDistance) {
220 if (particleCharge != 0.0) {
221 G4double endSafety =
223 currentSafety = endSafety ;
225 fPreviousSafety = currentSafety ;
226
227
228 // Because the Stepping Manager assumes it is from the start point,
229 // add the StepLength
230 //
231 currentSafety += fEndPointDistance ;
232
233#ifdef G4DEBUG_TRANSPORT
234 G4cout.precision(12) ;
235 G4cout << "***G4LongLivedNeutralTransportation::AlongStepGPIL ** " << G4endl ;
236 G4cout << " Called Navigator->ComputeSafety at " << fTransportEndPosition
237 << " and it returned safety= " << endSafety << G4endl ;
238 G4cout << " Adding endpoint distance " << fEndPointDistance
239 << " to obtain pseudo-safety= " << currentSafety << G4endl ;
240 } else {
241 G4cout << "***G4LongLivedNeutralTransportation::AlongStepGPIL ** " << G4endl ;
242 G4cout << " Avoiding call to ComputeSafety : " << G4endl;
243 G4cout << " charge = " << particleCharge << G4endl;
244 G4cout << " mag moment = " << magneticMoment << G4endl;
245#endif
246 }
247 }
248
249 fParticleChange.ProposeTrueStepLength(geometryStepLength) ;
250
251 return geometryStepLength ;
252}
G4bool fGeometryLimitedStep
Flag to determine whether a boundary was reached.

◆ AtRestDoIt()

G4VParticleChange * AtRestDoIt ( const G4Track &  ,
const G4Step &   
)
inline

No operation in AtRestDoIt.

Definition at line 177 of file G4LongLivedNeutralTransportation.h.

178 { return 0; }

◆ AtRestGetPhysicalInteractionLength()

G4double AtRestGetPhysicalInteractionLength ( const G4Track &  ,
G4ForceCondition *   
)
inline

No operation in AtRestGPIL.

Definition at line 170 of file G4LongLivedNeutralTransportation.h.

172 { return -1.0; }

◆ DoesGlobalFieldExist()

G4bool DoesGlobalFieldExist ( )
inlineprotected

Checks whether a field exists for the "global" field manager.

Definition at line 198 of file G4LongLivedNeutralTransportation.h.

199 {
200 G4TransportationManager* transportMgr;
201 transportMgr = G4TransportationManager::GetTransportationManager();
202
203 // fFieldExists= transportMgr->GetFieldManager()->DoesFieldExist();
204 // return fFieldExists;
205 return transportMgr->GetFieldManager()->DoesFieldExist();
206 }

◆ FieldExertedForce()

G4bool FieldExertedForce ( )
inline

References fFieldExertedForce.

Definition at line 109 of file G4LongLivedNeutralTransportation.h.

◆ GetMaxEnergyKilled()

G4double GetMaxEnergyKilled ( ) const
inline

Print values of looper thresholds.

Access fMaxEnergyKilled

◆ GetSilenceLooperWarnings()

G4bool GetSilenceLooperWarnings ( )
static

Do not throw exception about 'looping' particles.

Definition at line 515 of file G4LongLivedNeutralTransportation.cc.

516{
518}
static G4bool fSilenceLooperWarnings
Flag to Supress all 'looper' warnings.

◆ PostStepDoIt()

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

G4VProcess::PostStepDoIt() implementation,.

Parameters
track
stepData

Definition at line 346 of file G4LongLivedNeutralTransportation.cc.

348{
349 G4TouchableHandle retCurrentTouchable ; // The one to return
350 G4bool isLastStep = false;
351
352 // Initialize ParticleChange (by setting all its members equal
353 // to corresponding members in G4Track)
354 // fParticleChange.Initialize(track) ; // To initialise TouchableChange
355
356 fParticleChange.ProposeTrackStatus(track.GetTrackStatus()) ;
357
358 // If the Step was determined by the volume boundary,
359 // logically relocate the particle
360
362 // fCurrentTouchable will now become the previous touchable,
363 // and what was the previous will be freed.
364 // (Needed because the preStepPoint can point to the previous touchable)
365
366 fLinearNavigator->SetGeometricallyLimitedStep() ;
368 LocateGlobalPointAndUpdateTouchableHandle(track.GetPosition(),
369 track.GetMomentumDirection(),
371 true) ;
372 // Check whether the particle is out of the world volume
373 // If so it has exited and must be killed.
374 //
375 if (fCurrentTouchableHandle->GetVolume() == 0) {
376 fParticleChange.ProposeTrackStatus(fStopAndKill) ;
377 }
378 retCurrentTouchable = fCurrentTouchableHandle ;
379 fParticleChange.SetTouchableHandle(fCurrentTouchableHandle) ;
380
381 // Update the Step flag which identifies the Last Step in a volume
383 isLastStep = fLinearNavigator->ExitedMotherVolume()
384 or fLinearNavigator->EnteredDaughterVolume() ;
385 else
386 isLastStep = fFieldPropagator->IsLastStepInVolume();
387 } else { // fGeometryLimitedStep is false
388 // This serves only to move the Navigator's location
389 //
390 fLinearNavigator->LocateGlobalPointWithinVolume(track.GetPosition()) ;
391
392 // The value of the track's current Touchable is retained.
393 // (and it must be correct because we must use it below to
394 // overwrite the (unset) one in particle change)
395 // It must be fCurrentTouchable too ??
396 //
397 fParticleChange.SetTouchableHandle(track.GetTouchableHandle()) ;
398 retCurrentTouchable = track.GetTouchableHandle() ;
399
400 isLastStep = false;
401 } // endif ( fGeometryLimitedStep )
402 fLastStepInVolume = isLastStep;
403
404 fParticleChange.ProposeFirstStepInVolume(fFirstStepInVolume);
405 fParticleChange.ProposeLastStepInVolume(isLastStep);
406
407 const G4VPhysicalVolume* pNewVol = retCurrentTouchable->GetVolume() ;
408 G4Material* pNewMaterial = nullptr ;
409 G4VSensitiveDetector* pNewSensitiveDetector = nullptr ;
410
411 if (pNewVol != 0) {
412 pNewMaterial = pNewVol->GetLogicalVolume()->GetMaterial();
413 pNewSensitiveDetector = pNewVol->GetLogicalVolume()->GetSensitiveDetector();
414 }
415
416 fParticleChange.SetMaterialInTouchable(pNewMaterial) ;
417 fParticleChange.SetSensitiveDetectorInTouchable(pNewSensitiveDetector) ;
418
419 const G4MaterialCutsCouple* pNewMaterialCutsCouple = 0;
420 if (pNewVol != 0) {
421 pNewMaterialCutsCouple = pNewVol->GetLogicalVolume()->GetMaterialCutsCouple();
422 }
423
424 if (pNewVol != 0 && pNewMaterialCutsCouple != 0
425 && pNewMaterialCutsCouple->GetMaterial() != pNewMaterial) {
426 // for parametrized volume
427 //
428 pNewMaterialCutsCouple =
429 G4ProductionCutsTable::GetProductionCutsTable()
430 ->GetMaterialCutsCouple(pNewMaterial,
431 pNewMaterialCutsCouple->GetProductionCuts());
432 }
433 fParticleChange.SetMaterialCutsCoupleInTouchable(pNewMaterialCutsCouple);
434
435 // temporarily until Get/Set Material of ParticleChange,
436 // and StepPoint can be made const.
437 // Set the touchable in ParticleChange
438 // this must always be done because the particle change always
439 // uses this value to overwrite the current touchable pointer.
440 //
441 fParticleChange.SetTouchableHandle(retCurrentTouchable) ;
442
443 return &fParticleChange ;
444}

◆ PostStepGetPhysicalInteractionLength()

G4double PostStepGetPhysicalInteractionLength ( const G4Track &  track,
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.

337{
338 fFieldExertedForce = false; // Not known
339 *pForceCond = Forced ;
340 return DBL_MAX ; // was kInfinity ; but convention now is DBL_MAX
341}

◆ PrintStatistics()

void PrintStatistics ( std::ostream &  outStr) const

returns current logging info of the algorithm

Definition at line 94 of file G4LongLivedNeutralTransportation.cc.

95{
96 outStr << " G4LongLivedNeutralTransportation: Statistics for looping particles " << G4endl;
97 outStr << " No looping tracks found or killed. " << G4endl;
98
99}

◆ ProcessDescription()

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

G4LongLivedNeutralTransportation::ProcessDescription()

Parameters
outStreamOutput file with a description of process

Definition at line 558 of file G4LongLivedNeutralTransportation.cc.

562{
563 G4String indent = " "; // : "");
564 G4int oldPrec = outStream.precision(6);
565 // outStream << std::setprecision(6);
566 outStream << G4endl << indent << GetProcessName() << ": ";
567
568 outStream << " Parameters for looping particles: " << G4endl
569 << " warning-E = " << fThreshold_Warning_Energy / CLHEP::MeV << " MeV " << G4endl
570 << " important E = " << fThreshold_Important_Energy / CLHEP::MeV << " MeV " << G4endl
571 << " thresholdTrials " << fThresholdTrials << G4endl;
572 outStream.precision(oldPrec);
573}
G4int fThresholdTrials
Number of trials an important looper survives.
G4double fThreshold_Important_Energy
Give a few trial above this E for looping particle.
G4double fThreshold_Warning_Energy
Warn above this energy about looping particle.

◆ 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.

524{
525 // Restores the old high values -- potentially appropriate for energy-frontier
526 // HEP experiments.
527 // Caution: All tracks with E < 100 MeV that are found to loop are
528 SetThresholdWarningEnergy(100.0 * CLHEP::MeV); // Warn above this energy
529 SetThresholdImportantEnergy(250.0 * CLHEP::MeV); // Extra trial above this En
530
531 G4int maxTrials = 10;
532 SetThresholdTrials(maxTrials);
533
534
535
536}
void SetThresholdWarningEnergy(G4double newEnWarn)
Set fThreshold_Warning_Energy.
void SetThresholdImportantEnergy(G4double newEnImp)
Set fThreshold_Important_Energy.
void SetThresholdTrials(G4int newMaxTrials)
Set fThresholdTrials.

◆ SetLowLooperThresholds()

void SetLowLooperThresholds ( )

Shortcut method - old values (meant for HEP)

Definition at line 539 of file G4LongLivedNeutralTransportation.cc.

540{
541 // These values were the default in Geant4 10.5 - beta
542 SetThresholdWarningEnergy(1.0 * CLHEP::keV); // Warn above this En
543 SetThresholdImportantEnergy(1.0 * CLHEP::MeV); // Extra trials above it
544
545 G4int maxTrials = 30; // A new value - was 10
546 SetThresholdTrials(maxTrials);
547
548}

◆ SetSilenceLooperWarnings()

void SetSilenceLooperWarnings ( G4bool  val)
static

Do not warn about 'looping' particles.

Definition at line 507 of file G4LongLivedNeutralTransportation.cc.

508{
509 fSilenceLooperWarnings = val; // Flag to *Suppress* all 'looper' warnings
510 // G4CoupledTransportation::fSilenceLooperWarnings= val;
511}

◆ SetThresholdImportantEnergy()

void SetThresholdImportantEnergy ( G4double  newEnImp)
inline

Set fThreshold_Important_Energy.

Definition at line 127 of file G4LongLivedNeutralTransportation.h.

◆ SetThresholdTrials()

void SetThresholdTrials ( G4int  newMaxTrials)
inline

Set fThresholdTrials.

Definition at line 131 of file G4LongLivedNeutralTransportation.h.

◆ SetThresholdWarningEnergy()

void SetThresholdWarningEnergy ( G4double  newEnWarn)
inline

Set fThreshold_Warning_Energy.

Definition at line 123 of file G4LongLivedNeutralTransportation.h.

◆ StartTracking()

void StartTracking ( G4Track *  aTrack)

Reset state for new (potentially resumed) track.

Definition at line 453 of file G4LongLivedNeutralTransportation.cc.

454{
455 G4VProcess::StartTracking(aTrack);
456 fNewTrack = true;
457 fFirstStepInVolume = true;
458 fLastStepInVolume = false;
459
460 // The actions here are those that were taken in AlongStepGPIL
461 // when track.GetCurrentStepNumber()==1
462
463 // Whether field exists should be determined at run level -- TODO
465
466 // reset safety value and center
467 //
468 fPreviousSafety = 0.0 ;
469 fPreviousSftOrigin = G4ThreeVector(0., 0., 0.) ;
470
471 // reset looping counter -- for motion in field
472 fNoLooperTrials = 0;
473 // Must clear this state .. else it depends on last track's value
474 // --> a better solution would set this from state of suspended track TODO ?
475 // Was if( aTrack->GetCurrentStepNumber()==1 ) { .. }
476
477 // ChordFinder reset internal state
478 //
480 fFieldPropagator->ClearPropagatorState();
481 // Resets all state of field propagator class (ONLY) including safety
482 // values (in case of overlaps and to wipe for first track).
483 }
484
485 // Make sure to clear the chord finders of all fields (i.e. managers)
486 //
487 G4FieldManagerStore* fieldMgrStore = G4FieldManagerStore::GetInstance();
488 fieldMgrStore->ClearAllChordFindersState();
489
490 // Update the current touchable handle (from the track's)
491 //
492 fCurrentTouchableHandle = aTrack->GetTouchableHandle();
493
494 // Inform field propagator of new track
495 //
496 if (fFieldPropagator) fFieldPropagator->PrepareNewTrack();
497}
G4int fNoLooperTrials
Counter for steps in which particle reports 'looping', if it is above 'Important' Energy.

Friends And Related Function Documentation

◆ G4CoupledTransportation

friend class G4CoupledTransportation
friend

Definition at line 265 of file G4LongLivedNeutralTransportation.h.

Member Data Documentation

◆ fAbandonUnstableTrials

G4int fAbandonUnstableTrials = 0
private

Number of trials after which to abandon.

Definition at line 246 of file G4LongLivedNeutralTransportation.h.

◆ fAnyFieldExists

G4bool fAnyFieldExists = false
private

Flag for existing fields.

Definition at line 221 of file G4LongLivedNeutralTransportation.h.

◆ fCandidateEndGlobalTime

G4double fCandidateEndGlobalTime = 0.0
private

The particle's state after this Step, Store for DoIt.

Definition at line 219 of file G4LongLivedNeutralTransportation.h.

◆ fCurrentTouchableHandle

G4TouchableHandle fCurrentTouchableHandle
private

Current touchable handle.

Definition at line 230 of file G4LongLivedNeutralTransportation.h.

◆ fEndGlobalTimeComputed

G4bool fEndGlobalTimeComputed = false
private

The particle's state after this Step, Store for DoIt.

Definition at line 218 of file G4LongLivedNeutralTransportation.h.

◆ fEndPointDistance

G4double fEndPointDistance
private

Endpoint distance.

Definition at line 237 of file G4LongLivedNeutralTransportation.h.

◆ fFieldExertedForce

G4bool fFieldExertedForce = false
private

During current step.

Definition at line 228 of file G4LongLivedNeutralTransportation.h.

◆ fFieldPropagator

G4PropagatorInField* fFieldPropagator
private

Propagator used to transport the particle.

Definition at line 211 of file G4LongLivedNeutralTransportation.h.

◆ fFirstStepInVolume

G4bool fFirstStepInVolume = true
private

Flag first step in a geom.

volume

Definition at line 224 of file G4LongLivedNeutralTransportation.h.

◆ fGeometryLimitedStep

G4bool fGeometryLimitedStep = true
private

Flag to determine whether a boundary was reached.

Definition at line 226 of file G4LongLivedNeutralTransportation.h.

◆ fLastStepInVolume

G4bool fLastStepInVolume = false
private

Flag last step in a geom.

volume

Definition at line 225 of file G4LongLivedNeutralTransportation.h.

◆ fLinearNavigator

G4Navigator* fLinearNavigator
private

Propagator used to transport the particle.

Definition at line 210 of file G4LongLivedNeutralTransportation.h.

◆ fMomentumChanged

G4bool fMomentumChanged = true
private

The particle's state after this Step, Store for DoIt.

Definition at line 217 of file G4LongLivedNeutralTransportation.h.

◆ fNewTrack

G4bool fNewTrack = true
private

Flag from StartTracking.

Definition at line 223 of file G4LongLivedNeutralTransportation.h.

◆ fNoLooperTrials

G4int fNoLooperTrials = 0
private

Counter for steps in which particle reports 'looping', if it is above 'Important' Energy.

Definition at line 250 of file G4LongLivedNeutralTransportation.h.

◆ fParticleChange

G4ParticleChangeForTransport fParticleChange
private

New ParticleChange.

Definition at line 235 of file G4LongLivedNeutralTransportation.h.

◆ fPreviousSafety

G4double fPreviousSafety
private

Remember last safety value.

Definition at line 233 of file G4LongLivedNeutralTransportation.h.

◆ fPreviousSftOrigin

G4ThreeVector fPreviousSftOrigin
private

Remember last safety origin.

Definition at line 232 of file G4LongLivedNeutralTransportation.h.

◆ fpSafetyHelper

G4SafetyHelper* fpSafetyHelper
private

To pass it the safety value obtained.

Definition at line 261 of file G4LongLivedNeutralTransportation.h.

◆ fShortStepOptimisation

G4bool fShortStepOptimisation
private

Whether to avoid calling G4Navigator for short step ( < safety) If using it, the safety estimate for endpoint will likely be smaller.

Definition at line 259 of file G4LongLivedNeutralTransportation.h.

◆ fSilenceLooperWarnings

G4bool fSilenceLooperWarnings = false
staticprivate

Flag to Supress all 'looper' warnings.

Definition at line 268 of file G4LongLivedNeutralTransportation.h.

◆ fSumEnergyKilled

G4double fSumEnergyKilled = 0.0
private

Sum of abandoned looping tracks energies.

Definition at line 254 of file G4LongLivedNeutralTransportation.h.

◆ fThreshold_Important_Energy

G4double fThreshold_Important_Energy = 1.0 * CLHEP::MeV
private

Give a few trial above this E for looping particle.

Definition at line 242 of file G4LongLivedNeutralTransportation.h.

◆ fThreshold_Warning_Energy

G4double fThreshold_Warning_Energy = 1.0 * CLHEP::keV
private

Warn above this energy about looping particle.

Definition at line 241 of file G4LongLivedNeutralTransportation.h.

◆ fThresholdTrials

G4int fThresholdTrials = 10
private

Number of trials an important looper survives.

Definition at line 243 of file G4LongLivedNeutralTransportation.h.

◆ fTransportEndKineticEnergy

G4double fTransportEndKineticEnergy = 0.0
private

The particle's state after this Step, Store for DoIt.

Definition at line 215 of file G4LongLivedNeutralTransportation.h.

◆ fTransportEndMomentumDir

G4ThreeVector fTransportEndMomentumDir = G4ThreeVector(0.0, 0.0, 0.0)
private

The particle's state after this Step, Store for DoIt.

Definition at line 214 of file G4LongLivedNeutralTransportation.h.

◆ fTransportEndPosition

G4ThreeVector fTransportEndPosition = G4ThreeVector(0.0, 0.0, 0.0)
private

The particle's state after this Step, Store for DoIt.

Definition at line 213 of file G4LongLivedNeutralTransportation.h.

◆ fTransportEndSpin

G4ThreeVector fTransportEndSpin = G4ThreeVector(0.0, 0.0, 0.0)
private

The particle's state after this Step, Store for DoIt.

Definition at line 216 of file G4LongLivedNeutralTransportation.h.

◆ fUseGravity

G4bool fUseGravity
staticprivate

Flag take into account gravity.

Definition at line 267 of file G4LongLivedNeutralTransportation.h.

◆ fUseMagneticMoment

G4bool fUseMagneticMoment
staticprivate

Flag take into account magnetic moment.

Definition at line 266 of file G4LongLivedNeutralTransportation.h.


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