11#include <framework/logging/Logger.h>
13#include "G4TransportationProcessType.hh"
15#include "G4PhysicalConstants.hh"
16#include "G4SystemOfUnits.hh"
17#include "G4ProductionCutsTable.hh"
18#include "G4PrimaryParticle.hh"
20#include "G4FieldManagerStore.hh"
21#include <simulation/longlivedneutral/G4LongLivedNeutralTransportation.h>
23class G4VSensitiveDetector;
34 : G4VProcess(G4String(
"LongLivedNeutralTransportation"), fTransportation),
35 fFieldExertedForce(false),
36 fPreviousSftOrigin(0., 0., 0.),
38 fEndPointDistance(-1.0),
39 fShortStepOptimisation(false)
41 SetProcessSubType(
static_cast<G4int
>(TRANSPORTATION));
43 SetVerboseLevel(verbosity);
45 G4TransportationManager* transportMgr ;
47 transportMgr = G4TransportationManager::GetTransportationManager() ;
66 static G4ThreadLocal G4TouchableHandle* pNullTouchableHandle = 0;
67 pNullTouchableHandle =
new G4TouchableHandle;
73 if (verboseLevel > 0) {
74 G4cout <<
" G4LongLivedNeutralTransportation constructor> set fShortStepOptimisation to ";
76 else { G4cout <<
"false" << G4endl; }
96 outStr <<
" G4LongLivedNeutralTransportation: Statistics for looping particles " << G4endl;
97 outStr <<
" No looping tracks found or killed. " << G4endl;
111 G4double currentMinimumStep,
112 G4double& currentSafety,
113 G4GPILSelection* selection)
115 G4double geometryStepLength = -1.0, newSafety = -1.0;
126 *selection = CandidateForSelection ;
136 const G4DynamicParticle* pParticle = track.GetDynamicParticle() ;
137 const G4PrimaryParticle* primaryParticle = pParticle->GetPrimaryParticle();
138 G4ThreeVector startMomentumDir = primaryParticle->GetMomentumDirection() ;
139 G4ThreeVector startPosition = track.GetPosition() ;
146 G4double MagSqShift = OriginShift.mag2() ;
148 currentSafety = 0.0 ;
155 G4double particleCharge = 0 ;
163 G4double linearStepLength ;
167 geometryStepLength = currentMinimumStep ;
182 currentSafety = newSafety ;
187 geometryStepLength = linearStepLength ;
190 geometryStepLength = currentMinimumStep ;
212 if (currentMinimumStep == 0.0) {
220 if (particleCharge != 0.0) {
223 currentSafety = endSafety ;
233#ifdef G4DEBUG_TRANSPORT
234 G4cout.precision(12) ;
235 G4cout <<
"***G4LongLivedNeutralTransportation::AlongStepGPIL ** " << G4endl ;
237 <<
" and it returned safety= " << endSafety << G4endl ;
239 <<
" to obtain pseudo-safety= " << currentSafety << G4endl ;
241 G4cout <<
"***G4LongLivedNeutralTransportation::AlongStepGPIL ** " << G4endl ;
242 G4cout <<
" Avoiding call to ComputeSafety : " << G4endl;
243 G4cout <<
" charge = " << particleCharge << G4endl;
244 G4cout <<
" mag moment = " << magneticMoment << G4endl;
251 return geometryStepLength ;
260 const G4Step& stepData)
276 G4double deltaTime = 0.0 ;
282 G4double startTime = track.GetGlobalTime() ;
284 const G4DynamicParticle* pParticle = track.GetDynamicParticle() ;
285 const G4PrimaryParticle* primaryParticle = pParticle->GetPrimaryParticle();
294 B2DEBUG(0,
"Velocity calculated from Step data: " << stepData.GetPreStepPoint()->GetVelocity());
296 G4double initialVelocity = c_light * primaryParticle->GetTotalMomentum() / primaryParticle->GetTotalEnergy();
297 G4double stepLength = track.GetStepLength();
300 if (initialVelocity > 0.0) { deltaTime = stepLength / initialVelocity; }
310 G4double restMass = primaryParticle->GetMass() ;
311 G4double deltaProperTime = deltaTime * (restMass / primaryParticle->GetTotalEnergy()) ;
313 fParticleChange.ProposeProperTime(track.GetProperTime() + deltaProperTime) ;
336 G4ForceCondition* pForceCond)
339 *pForceCond = Forced ;
349 G4TouchableHandle retCurrentTouchable ;
350 G4bool isLastStep =
false;
368 LocateGlobalPointAndUpdateTouchableHandle(track.GetPosition(),
369 track.GetMomentumDirection(),
398 retCurrentTouchable = track.GetTouchableHandle() ;
407 const G4VPhysicalVolume* pNewVol = retCurrentTouchable->GetVolume() ;
408 G4Material* pNewMaterial = nullptr ;
409 G4VSensitiveDetector* pNewSensitiveDetector = nullptr ;
412 pNewMaterial = pNewVol->GetLogicalVolume()->GetMaterial();
413 pNewSensitiveDetector = pNewVol->GetLogicalVolume()->GetSensitiveDetector();
417 fParticleChange.SetSensitiveDetectorInTouchable(pNewSensitiveDetector) ;
419 const G4MaterialCutsCouple* pNewMaterialCutsCouple = 0;
421 pNewMaterialCutsCouple = pNewVol->GetLogicalVolume()->GetMaterialCutsCouple();
424 if (pNewVol != 0 && pNewMaterialCutsCouple != 0
425 && pNewMaterialCutsCouple->GetMaterial() != pNewMaterial) {
428 pNewMaterialCutsCouple =
429 G4ProductionCutsTable::GetProductionCutsTable()
430 ->GetMaterialCutsCouple(pNewMaterial,
431 pNewMaterialCutsCouple->GetProductionCuts());
433 fParticleChange.SetMaterialCutsCoupleInTouchable(pNewMaterialCutsCouple);
455 G4VProcess::StartTracking(aTrack);
487 G4FieldManagerStore* fieldMgrStore = G4FieldManagerStore::GetInstance();
488 fieldMgrStore->ClearAllChordFindersState();
531 G4int maxTrials = 10;
545 G4int maxTrials = 30;
563 G4String indent =
" ";
564 G4int oldPrec = outStream.precision(6);
566 outStream << G4endl << indent << GetProcessName() <<
": ";
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);
G4bool fFirstStepInVolume
Flag first step in a geom.
G4TouchableHandle fCurrentTouchableHandle
Current touchable handle.
void PrintStatistics(std::ostream &outStr) const
returns current logging info of the algorithm
G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData)
G4VProcess::PostStepDoIt() implementation,.
G4bool fAnyFieldExists
Flag for existing fields.
G4double fCandidateEndGlobalTime
The particle's state after this Step, Store for DoIt.
G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double ¤tSafety, G4GPILSelection *selection)
G4VProcess::AlongStepGetPhysicalInteractionLength() implementation,.
void SetThresholdWarningEnergy(G4double newEnWarn)
Set fThreshold_Warning_Energy.
void SetThresholdImportantEnergy(G4double newEnImp)
Set fThreshold_Important_Energy.
virtual void ProcessDescription(std::ostream &outStream) const
G4LongLivedNeutralTransportation::ProcessDescription()
G4bool fGeometryLimitedStep
Flag to determine whether a boundary was reached.
G4bool fFieldExertedForce
During current step.
G4ThreeVector fTransportEndMomentumDir
The particle's state after this Step, Store for DoIt.
void SetThresholdTrials(G4int newMaxTrials)
Set fThresholdTrials.
G4double fPreviousSafety
Remember last safety value.
G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData)
G4VProcess::AlongStepDoIt() implementation,.
G4double fEndPointDistance
Endpoint distance.
G4bool fLastStepInVolume
Flag last step in a geom.
void SetHighLooperThresholds()
Get/Set parameters for killing loopers: Above 'important' energy a 'looping' particle in field will N...
G4double fSumEnergyKilled
Sum of abandoned looping tracks energies.
G4ThreeVector fPreviousSftOrigin
Remember last safety origin.
static void SetSilenceLooperWarnings(G4bool val)
Do not warn about 'looping' particles.
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *pForceCond)
G4VProcess::PostStepGetPhysicalInteractionLength() implementation.
G4ThreeVector fTransportEndSpin
The particle's state after this Step, Store for DoIt.
G4Navigator * fLinearNavigator
Propagator used to transport the particle.
static G4bool fSilenceLooperWarnings
Flag to Supress all 'looper' warnings.
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 fShortStepOptimisation
Whether to avoid calling G4Navigator for short step ( < safety) If using it, the safety estimate for ...
G4int fNoLooperTrials
Counter for steps in which particle reports 'looping', if it is above 'Important' Energy.
G4bool fMomentumChanged
The particle's state after this Step, Store for DoIt.
G4bool fNewTrack
Flag from StartTracking.
G4LongLivedNeutralTransportation(G4int verbosityLevel=1)
Constructor.
void StartTracking(G4Track *aTrack)
Reset state for new (potentially resumed) track.
G4ParticleChangeForTransport fParticleChange
New ParticleChange.
void SetLowLooperThresholds()
Shortcut method - old values (meant for HEP)
G4bool DoesGlobalFieldExist()
Checks whether a field exists for the "global" field manager.
~G4LongLivedNeutralTransportation()
Destructor.
static G4bool GetSilenceLooperWarnings()
Do not throw exception about 'looping' particles.
G4PropagatorInField * fFieldPropagator
Propagator used to transport the particle.
G4double fTransportEndKineticEnergy
The particle's state after this Step, Store for DoIt.
Abstract base class for different kinds of events.