11 #include <simulation/monopoles/G4MonopoleTransportation.h>
12 #include <simulation/monopoles/G4Monopole.h>
14 #include <G4ProductionCutsTable.hh>
15 #include <G4ParticleTable.hh>
16 #include <G4ChordFinder.hh>
17 #include <G4SafetyHelper.hh>
18 #include <G4FieldManagerStore.hh>
19 #include <G4TransportationProcessType.hh>
20 #include <CLHEP/Units/SystemOfUnits.h>
21 #include <framework/logging/Logger.h>
23 class G4VSensitiveDetector;
27 using namespace Belle2::Monopoles;
28 using namespace CLHEP;
30 G4MonopoleTransportation::G4MonopoleTransportation(
const G4Monopole* mpl,
32 : G4VProcess(G4String(
"MonopoleTransportation"), fTransportation),
37 fParticleIsLooping(false),
38 fPreviousSftOrigin(0., 0., 0.),
40 fThreshold_Trap_Energy(10 * MeV),
41 fThreshold_Warning_Energy(100 * MeV),
42 fThreshold_Important_Energy(250 * MeV),
45 fSumEnergyKilled(0.0), fMaxEnergyKilled(0.0),
46 fShortStepOptimisation(false),
52 SetProcessSubType(TRANSPORTATION);
56 G4TransportationManager* transportMgr ;
58 transportMgr = G4TransportationManager::GetTransportationManager() ;
74 static G4TouchableHandle nullTouchableHandle;
84 B2DEBUG(25,
" G4MonopoleTransportation: Statistics for looping particles ");
98 G4double currentMinimumStep,
99 G4double& currentSafety,
100 G4GPILSelection* selection)
105 G4double geometryStepLength, newSafety ;
117 *selection = CandidateForSelection ;
121 const G4DynamicParticle* pParticle = track.GetDynamicParticle() ;
122 G4ThreeVector startMomentumDir = pParticle->GetMomentumDirection() ;
123 G4ThreeVector startPosition = track.GetPosition() ;
132 G4double MagSqShift = OriginShift.mag2() ;
134 currentSafety = 0.0 ;
142 G4double particleElectricCharge = pParticle->GetCharge();
153 G4FieldManager* fieldMgr =
nullptr;
154 G4bool fieldExertsForce = false ;
156 if ((particleMagneticCharge != 0.0)) {
160 fieldMgr->ConfigureForTrack(&track);
166 fieldExertsForce = (fieldMgr->GetDetectorField() != 0);
172 if (!fieldExertsForce) {
173 G4double linearStepLength ;
177 geometryStepLength = currentMinimumStep ;
194 currentSafety = newSafety ;
199 geometryStepLength = linearStepLength ;
202 geometryStepLength = currentMinimumStep ;
220 G4double momentumMagnitude = pParticle->GetTotalMomentum() ;
221 G4ThreeVector EndUnitMomentum ;
222 G4double lengthAlongCurve ;
225 G4ChargeState chargeState(particleElectricCharge,
229 particleMagneticCharge);
231 G4EquationOfMotion* equationOfMotion =
233 ->GetEquationOfMotion();
236 ->SetChargeMomentumMass(chargeState,
241 G4ThreeVector spin = track.GetPolarization() ;
242 G4FieldTrack aFieldTrack = G4FieldTrack(startPosition,
243 track.GetMomentumDirection(),
245 track.GetKineticEnergy(),
248 track.GetGlobalTime(),
249 track.GetProperTime(),
251 if (currentMinimumStep > 0) {
260 geometryStepLength = lengthAlongCurve ;
262 geometryStepLength = currentMinimumStep ;
265 geometryStepLength = lengthAlongCurve = 0.0 ;
297 if (currentMinimumStep == 0.0) {
308 if (particleMagneticCharge != 0.0) {
312 currentSafety = endSafety ;
322 #ifdef G4DEBUG_TRANSPORT
323 G4cout.precision(12) ;
324 G4cout <<
"***G4MonopoleTransportation::AlongStepGPIL ** " << G4endl ;
326 <<
" and it returned safety= " << endSafety << G4endl ;
328 <<
" to obtain pseudo-safety= " << currentSafety << G4endl ;
338 return geometryStepLength ;
346 const G4Step& stepData)
348 static G4int noCalls = 0;
349 static const G4ParticleDefinition* fOpticalPhoton =
350 G4ParticleTable::GetParticleTable()->FindParticle(
"opticalphoton");
365 G4double deltaTime = 0.0 ;
371 G4double startTime = track.GetGlobalTime() ;
376 G4double finalVelocity = track.GetVelocity() ;
377 G4double initialVelocity = stepData.GetPreStepPoint()->GetVelocity() ;
378 G4double stepLength = track.GetStepLength() ;
381 const G4DynamicParticle* fpDynamicParticle = track.GetDynamicParticle();
382 if (fpDynamicParticle->GetDefinition() == fOpticalPhoton) {
385 deltaTime = stepLength / finalVelocity ;
386 }
else if (finalVelocity > 0.0) {
387 G4double meanInverseVelocity ;
389 meanInverseVelocity = 0.5
390 * (1.0 / initialVelocity + 1.0 / finalVelocity) ;
391 deltaTime = stepLength * meanInverseVelocity ;
392 }
else if (initialVelocity > 0.0) {
393 deltaTime = stepLength / initialVelocity ;
404 G4double restMass = track.GetDynamicParticle()->GetMass() ;
405 G4double deltaProperTime = deltaTime * (restMass / track.GetTotalEnergy()) ;
407 fParticleChange.ProposeProperTime(track.GetProperTime() + deltaProperTime) ;
428 if ((verboseLevel > 1) ||
430 B2DEBUG(25,
" G4MonopoleTransportation is killing track that is looping or stuck ");
431 B2DEBUG(25,
" This track has " << track.GetKineticEnergy() / MeV <<
" MeV energy.");
433 <<
" No of calls to AlongStepDoIt = " << noCalls);
440 if ((verboseLevel > 2)) {
441 B2DEBUG(25,
" G4MonopoleTransportation::AlongStepDoIt(): Particle looping - "
443 <<
" No of calls to = " << noCalls);
476 G4ForceCondition* pForceCond)
478 *pForceCond = Forced ;
485 G4TouchableHandle retCurrentTouchable ;
503 LocateGlobalPointAndUpdateTouchableHandle(track.GetPosition(),
504 track.GetMomentumDirection(),
526 retCurrentTouchable = track.GetTouchableHandle() ;
529 const G4VPhysicalVolume* pNewVol = retCurrentTouchable->GetVolume() ;
530 G4Material* pNewMaterial = nullptr ;
531 G4VSensitiveDetector* pNewSensitiveDetector = nullptr ;
534 pNewMaterial = pNewVol->GetLogicalVolume()->GetMaterial();
535 pNewSensitiveDetector = pNewVol->GetLogicalVolume()->GetSensitiveDetector();
541 pNewSensitiveDetector) ;
543 const G4MaterialCutsCouple* pNewMaterialCutsCouple = 0;
545 pNewMaterialCutsCouple = pNewVol->GetLogicalVolume()->GetMaterialCutsCouple();
548 if (pNewVol != 0 && pNewMaterialCutsCouple != 0 &&
549 pNewMaterialCutsCouple->GetMaterial() != pNewMaterial) {
552 pNewMaterialCutsCouple =
553 G4ProductionCutsTable::GetProductionCutsTable()
554 ->GetMaterialCutsCouple(pNewMaterial,
555 pNewMaterialCutsCouple->GetProductionCuts());
557 fParticleChange.SetMaterialCutsCoupleInTouchable(pNewMaterialCutsCouple);
577 G4VProcess::StartTracking(aTrack);
605 static G4FieldManagerStore* fieldMgrStore = G4FieldManagerStore::GetInstance();
606 fieldMgrStore->ClearAllChordFindersState();
void SwitchChordFinder(G4int val)
Switches chord finder between 1 - basf2 FullSim chord finder 2 - monopole chord finder Since monopole...
static G4MonopoleFieldSetup * GetMonopoleFieldSetup()
Returns G4MonopoleFieldSetup instance.
G4TouchableHandle fCurrentTouchableHandle
Current touchable handle.
G4double endpointDistance
Endpint distance.
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData)
G4VProcess::PostStepDoIt() implementation.
G4double fCandidateEndGlobalTime
The particle's state after this Step, Store for DoIt.
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double ¤tSafety, G4GPILSelection *selection)
G4VProcess::AlongStepGetPhysicalInteractionLength() implementation.
G4int fThresholdTrials
Nubmer of trials for looping particles.
G4double fMaxEnergyKilled
Max of abandoned looping tracks energies.
G4bool fGeometryLimitedStep
Flag to determine whether a boundary was reached.
~G4MonopoleTransportation()
Destructor.
G4ThreeVector fTransportEndMomentumDir
The particle's state after this Step, Store for DoIt.
G4double fPreviousSafety
Remember last safety value.
virtual G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData)
G4VProcess::AlongStepDoIt() implementation, Proposes changes during step to fParticleChange of this c...
G4SafetyHelper * fpSafetyHelper
To pass it the safety value obtained.
G4double fSumEnergyKilled
Sum of abandoned looping tracks energies.
G4ThreeVector fPreviousSftOrigin
Remember last safety origin.
virtual 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.
G4double fThreshold_Important_Energy
Hesitate above this about looping particle for a certain no of trials.
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.
virtual void StartTracking(G4Track *aTrack)
Reset state for new (potentially resumed) track.
G4ParticleChangeForTransport fParticleChange
New ParticleChange.
G4double fThreshold_Warning_Energy
Warn above this energy about looping particle.
G4bool DoesGlobalFieldExist()
Checks whether a field exists for the "global" field manager.
const G4Monopole * fParticleDef
Monopole definition for charge and mass reference.
G4bool fParticleIsLooping
Is the monopole stuck in looping.
G4double fThreshold_Trap_Energy
Assume monopoles below this can bound to material.
G4PropagatorInField * fFieldPropagator
Propagator used to transport the particle.
G4double fTransportEndKineticEnergy
The particle's state after this Step, Store for DoIt.
G4MonopoleFieldSetup * fMagSetup
Monpole field setup.
A class to hold monopole description as a particle.
G4double MagneticCharge() const
Returns magnetic charge of the monopole.
Abstract base class for different kinds of events.