43 #include "G4MonopoleTransportation.hh" 44 #include "G4ProductionCutsTable.hh" 45 #include "G4ParticleTable.hh" 46 #include "G4ChordFinder.hh" 47 #include "G4SafetyHelper.hh" 48 #include "G4FieldManagerStore.hh" 49 #include "SimG4Core/Physics/interface/G4Monopole.hh" 53 class G4VSensitiveDetector;
60 static const G4TouchableHandle nullTouchableHandle;
63 G4MonopoleTransportation::G4MonopoleTransportation(
const G4Monopole* mpl,
66 : G4VProcess( G4String(
"MonopoleTransportation"), fTransportation ),
67 fChordFinderSetter(chordFinderSetter),
68 fParticleIsLooping(
false ),
69 fPreviousSftOrigin (0.,0.,0.),
70 fPreviousSafety ( 0.0 ),
71 fThreshold_Warning_Energy( 100 *
MeV ),
72 fThreshold_Important_Energy( 250 *
MeV ),
73 fThresholdTrials( 10 ),
74 fUnimportant_Energy( 1 *
MeV ),
76 fSumEnergyKilled( 0.0 ), fMaxEnergyKilled( 0.0 ),
77 fShortStepOptimisation(
false),
84 G4TransportationManager* transportMgr ;
86 transportMgr = G4TransportationManager::GetTransportationManager() ;
88 fLinearNavigator = transportMgr->GetNavigatorForTracking() ;
92 fFieldPropagator = transportMgr->GetPropagatorInField() ;
94 fpSafetyHelper = transportMgr->GetSafetyHelper();
102 fCurrentTouchableHandle = nullTouchableHandle;
104 fEndGlobalTimeComputed =
false;
105 fCandidateEndGlobalTime = 0;
110 G4MonopoleTransportation::~G4MonopoleTransportation()
112 if( (fVerboseLevel > 0) && (fSumEnergyKilled > 0.0 ) ){
113 G4cout <<
" G4MonopoleTransportation: Statistics for looping particles " << G4endl;
114 G4cout <<
" Sum of energy of loopers killed: " << fSumEnergyKilled << G4endl;
115 G4cout <<
" Max energy of loopers killed: " << fMaxEnergyKilled << G4endl;
126 G4double G4MonopoleTransportation::
127 AlongStepGetPhysicalInteractionLength(
const G4Track&
track,
129 G4double currentMinimumStep,
130 G4double& currentSafety,
134 G4FieldManager* fieldMgr = fFieldPropagator->FindAndSetFieldManager(track.GetVolume());
138 G4double geometryStepLength, newSafety ;
139 fParticleIsLooping =
false ;
150 *selection = CandidateForSelection ;
154 const G4DynamicParticle* pParticle = track.GetDynamicParticle() ;
155 G4ThreeVector startMomentumDir = pParticle->GetMomentumDirection() ;
156 G4ThreeVector startPosition = track.GetPosition() ;
164 G4ThreeVector OriginShift = startPosition - fPreviousSftOrigin ;
165 G4double MagSqShift = OriginShift.mag2() ;
166 if( MagSqShift >=
sqr(fPreviousSafety) )
168 currentSafety = 0.0 ;
172 currentSafety = fPreviousSafety -
std::sqrt(MagSqShift) ;
177 G4double particleMagneticCharge = pParticleDef->MagneticCharge() ;
178 G4double particleElectricCharge = pParticle->GetCharge();
180 fGeometryLimitedStep =
false ;
189 G4bool fieldExertsForce =
false ;
191 if( (particleMagneticCharge != 0.0) )
195 fieldMgr->ConfigureForTrack( &track );
201 fieldExertsForce = (fieldMgr->GetDetectorField() != 0);
210 if( !fieldExertsForce )
212 G4double linearStepLength ;
213 if( fShortStepOptimisation && (currentMinimumStep <= currentSafety) )
217 geometryStepLength = currentMinimumStep ;
218 fGeometryLimitedStep =
false ;
224 linearStepLength = fLinearNavigator->ComputeStep( startPosition,
230 fPreviousSftOrigin = startPosition ;
231 fPreviousSafety = newSafety ;
236 currentSafety = newSafety ;
238 fGeometryLimitedStep= (linearStepLength <= currentMinimumStep);
239 if( fGeometryLimitedStep )
242 geometryStepLength = linearStepLength ;
247 geometryStepLength = currentMinimumStep ;
250 endpointDistance = geometryStepLength ;
254 fTransportEndPosition = startPosition+geometryStepLength*startMomentumDir ;
258 fTransportEndMomentumDir = startMomentumDir ;
259 fTransportEndKineticEnergy = track.GetKineticEnergy() ;
260 fTransportEndSpin = track.GetPolarization();
261 fParticleIsLooping =
false ;
262 fMomentumChanged =
false ;
263 fEndGlobalTimeComputed =
false ;
268 G4ThreeVector EndUnitMomentum ;
269 G4double lengthAlongCurve ;
270 G4double restMass = pParticleDef->GetPDGMass();
271 G4double momentumMagnitude = pParticle->GetTotalMomentum();
276 G4ChargeState chargeState(particleElectricCharge,
277 pParticleDef->GetPDGSpin(),
280 particleMagneticCharge );
282 G4EquationOfMotion* equationOfMotion =
283 (fFieldPropagator->GetChordFinder()->GetIntegrationDriver()->GetStepper())
284 ->GetEquationOfMotion();
286 equationOfMotion->SetChargeMomentumMass( chargeState,
293 G4ThreeVector
spin = track.GetPolarization() ;
294 G4FieldTrack aFieldTrack = G4FieldTrack( startPosition,
295 track.GetMomentumDirection(),
297 track.GetKineticEnergy(),
300 track.GetGlobalTime(),
301 track.GetProperTime(),
303 if( currentMinimumStep > 0 )
307 lengthAlongCurve = fFieldPropagator->ComputeStep( aFieldTrack,
310 track.GetVolume() ) ;
311 fGeometryLimitedStep= lengthAlongCurve < currentMinimumStep;
312 if( fGeometryLimitedStep ) {
313 geometryStepLength = lengthAlongCurve ;
315 geometryStepLength = currentMinimumStep ;
320 geometryStepLength = lengthAlongCurve= 0.0 ;
321 fGeometryLimitedStep =
false ;
326 fPreviousSftOrigin = startPosition ;
327 fPreviousSafety = currentSafety ;
332 fTransportEndPosition = aFieldTrack.GetPosition() ;
336 fMomentumChanged =
true ;
337 fTransportEndMomentumDir = aFieldTrack.GetMomentumDir() ;
339 fTransportEndKineticEnergy = aFieldTrack.GetKineticEnergy() ;
346 fCandidateEndGlobalTime = aFieldTrack.GetLabTimeOfFlight();
347 fEndGlobalTimeComputed =
true;
414 fTransportEndSpin = aFieldTrack.GetSpin();
415 fParticleIsLooping = fFieldPropagator->IsParticleLooping() ;
416 endpointDistance = (fTransportEndPosition - startPosition).
mag() ;
422 if( currentMinimumStep == 0.0 )
424 if( currentSafety == 0.0 ) fGeometryLimitedStep =
true ;
430 if( currentSafety < endpointDistance )
435 if( particleMagneticCharge != 0.0 ) {
438 fLinearNavigator->ComputeSafety( fTransportEndPosition) ;
439 currentSafety = endSafety ;
440 fPreviousSftOrigin = fTransportEndPosition ;
441 fPreviousSafety = currentSafety ;
442 fpSafetyHelper->SetCurrentSafety( currentSafety, fTransportEndPosition);
447 currentSafety += endpointDistance ;
449 #ifdef G4DEBUG_TRANSPORT 451 G4cout <<
"***G4MonopoleTransportation::AlongStepGPIL ** " << G4endl ;
452 G4cout <<
" Called Navigator->ComputeSafety at " << fTransportEndPosition
453 <<
" and it returned safety= " << endSafety << G4endl ;
454 G4cout <<
" Adding endpoint distance " << endpointDistance
455 <<
" to obtain pseudo-safety= " << currentSafety << G4endl ;
460 fParticleChange.ProposeTrueStepLength(geometryStepLength) ;
462 return geometryStepLength ;
470 G4VParticleChange* G4MonopoleTransportation::AlongStepDoIt(
const G4Track& track,
471 const G4Step& stepData )
473 static const G4ParticleDefinition* fOpticalPhoton =
474 G4ParticleTable::GetParticleTable()->FindParticle(
"opticalphoton");
477 static thread_local G4int noCalls=0;
481 fParticleChange.Initialize(track) ;
485 fParticleChange.ProposePosition(fTransportEndPosition) ;
486 fParticleChange.ProposeMomentumDirection(fTransportEndMomentumDir) ;
487 fParticleChange.ProposeEnergy(fTransportEndKineticEnergy) ;
488 fParticleChange.SetMomentumChanged(fMomentumChanged) ;
490 fParticleChange.ProposePolarization(fTransportEndSpin);
492 G4double deltaTime = 0.0 ;
498 G4double startTime = track.GetGlobalTime() ;
500 if (!fEndGlobalTimeComputed)
504 G4double finalVelocity = track.GetVelocity() ;
505 G4double initialVelocity = stepData.GetPreStepPoint()->GetVelocity() ;
506 G4double stepLength = track.GetStepLength() ;
509 const G4DynamicParticle* fpDynamicParticle = track.GetDynamicParticle();
510 if (fpDynamicParticle->GetDefinition()== fOpticalPhoton)
514 deltaTime = stepLength/finalVelocity ;
516 else if (finalVelocity > 0.0)
518 G4double meanInverseVelocity ;
520 meanInverseVelocity = 0.5
521 * ( 1.0 / initialVelocity + 1.0 / finalVelocity ) ;
522 deltaTime = stepLength * meanInverseVelocity ;
524 else if( initialVelocity > 0.0 )
526 deltaTime = stepLength/initialVelocity ;
528 fCandidateEndGlobalTime = startTime + deltaTime ;
532 deltaTime = fCandidateEndGlobalTime - startTime ;
535 fParticleChange.ProposeGlobalTime( fCandidateEndGlobalTime ) ;
539 G4double restMass = track.GetDynamicParticle()->GetMass() ;
540 G4double deltaProperTime = deltaTime*( restMass/track.GetTotalEnergy() ) ;
542 fParticleChange.ProposeProperTime(track.GetProperTime() + deltaProperTime) ;
549 if ( fParticleIsLooping )
551 G4double endEnergy= fTransportEndKineticEnergy;
553 if( (endEnergy < fThreshold_Important_Energy)
554 || (fNoLooperTrials >= fThresholdTrials ) ){
557 fParticleChange.ProposeTrackStatus( fStopAndKill ) ;
560 fSumEnergyKilled += endEnergy;
561 if( endEnergy > fMaxEnergyKilled) { fMaxEnergyKilled= endEnergy; }
564 if( (fVerboseLevel > 1) ||
565 ( endEnergy > fThreshold_Warning_Energy ) ) {
566 G4cout <<
" G4MonopoleTransportation is killing track that is looping or stuck " 568 <<
" This track has " << track.GetKineticEnergy() /
MeV 569 <<
" MeV energy." << G4endl;
570 G4cout <<
" Number of trials = " << fNoLooperTrials
571 <<
" No of calls to AlongStepDoIt = " << noCalls
580 if( (fVerboseLevel > 2) ){
581 G4cout <<
" G4MonopoleTransportation::AlongStepDoIt(): Particle looping - " 582 <<
" Number of trials = " << fNoLooperTrials
583 <<
" No of calls to = " << noCalls
597 fParticleChange.SetPointerToVectorOfAuxiliaryPoints
598 (fFieldPropagator->GimmeTrajectoryVectorAndForgetIt() );
600 return &fParticleChange ;
609 G4double G4MonopoleTransportation::
610 PostStepGetPhysicalInteractionLength(
const G4Track& track,
612 G4ForceCondition* pForceCond )
615 G4FieldManager* fieldMgr=fFieldPropagator->FindAndSetFieldManager( track.GetVolume() );
619 *pForceCond = Forced ;
626 G4VParticleChange* G4MonopoleTransportation::PostStepDoIt(
const G4Track& track,
629 G4TouchableHandle retCurrentTouchable ;
635 fParticleChange.ProposeTrackStatus(track.GetTrackStatus()) ;
640 if(fGeometryLimitedStep)
646 fLinearNavigator->SetGeometricallyLimitedStep() ;
648 LocateGlobalPointAndUpdateTouchableHandle( track.GetPosition(),
649 track.GetMomentumDirection(),
650 fCurrentTouchableHandle,
655 if( fCurrentTouchableHandle->GetVolume() == 0 )
657 fParticleChange.ProposeTrackStatus( fStopAndKill ) ;
659 retCurrentTouchable = fCurrentTouchableHandle ;
660 fParticleChange.SetTouchableHandle( fCurrentTouchableHandle ) ;
666 fLinearNavigator->LocateGlobalPointWithinVolume( track.GetPosition() ) ;
673 fParticleChange.SetTouchableHandle( track.GetTouchableHandle() ) ;
674 retCurrentTouchable = track.GetTouchableHandle() ;
677 const G4VPhysicalVolume* pNewVol = retCurrentTouchable->GetVolume() ;
678 const G4Material* pNewMaterial = 0 ;
679 const G4VSensitiveDetector* pNewSensitiveDetector = 0 ;
683 pNewMaterial= pNewVol->GetLogicalVolume()->GetMaterial();
684 pNewSensitiveDetector= pNewVol->GetLogicalVolume()->GetSensitiveDetector();
690 fParticleChange.SetMaterialInTouchable( (G4Material *) pNewMaterial ) ;
691 fParticleChange.SetSensitiveDetectorInTouchable( (G4VSensitiveDetector *) pNewSensitiveDetector ) ;
693 const G4MaterialCutsCouple* pNewMaterialCutsCouple = 0;
696 pNewMaterialCutsCouple=pNewVol->GetLogicalVolume()->GetMaterialCutsCouple();
699 if( pNewVol!=0 && pNewMaterialCutsCouple!=0 && pNewMaterialCutsCouple->GetMaterial()!=pNewMaterial )
703 pNewMaterialCutsCouple =
704 G4ProductionCutsTable::GetProductionCutsTable()
705 ->GetMaterialCutsCouple(pNewMaterial,
706 pNewMaterialCutsCouple->GetProductionCuts());
708 fParticleChange.SetMaterialCutsCoupleInTouchable( pNewMaterialCutsCouple );
716 fParticleChange.SetTouchableHandle(retCurrentTouchable) ;
718 return &fParticleChange ;
725 G4MonopoleTransportation::StartTracking(G4Track* aTrack)
727 G4VProcess::StartTracking(aTrack);
734 fPreviousSafety = 0.0 ;
735 fPreviousSftOrigin = G4ThreeVector(0.,0.,0.) ;
745 if( DoesGlobalFieldExist() ) {
746 fFieldPropagator->ClearPropagatorState();
755 G4FieldManagerStore* fieldMgrStore= G4FieldManagerStore::GetInstance();
756 fieldMgrStore->ClearAllChordFindersState();
760 fCurrentTouchableHandle = aTrack->GetTouchableHandle();
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
Square< F >::type sqr(const F &f)
void SetMonopoleTracking(G4bool)