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"
51 class G4VSensitiveDetector;
58 static const G4TouchableHandle nullTouchableHandle;
61 G4MonopoleTransportation::G4MonopoleTransportation(
const G4Monopole* mpl,
64 : G4VProcess( G4String(
"MonopoleTransportation"), fTransportation ),
65 fFieldBuilder(fieldBuilder),
66 fParticleIsLooping(
false ),
67 fPreviousSftOrigin (0.,0.,0.),
68 fPreviousSafety ( 0.0 ),
69 fThreshold_Warning_Energy( 100 * MeV ),
70 fThreshold_Important_Energy( 250 * MeV ),
71 fThresholdTrials( 10 ),
72 fUnimportant_Energy( 1 * MeV ),
74 fSumEnergyKilled( 0.0 ), fMaxEnergyKilled( 0.0 ),
75 fShortStepOptimisation(
false),
82 G4TransportationManager* transportMgr ;
84 transportMgr = G4TransportationManager::GetTransportationManager() ;
86 fLinearNavigator = transportMgr->GetNavigatorForTracking() ;
90 fFieldPropagator = transportMgr->GetPropagatorInField() ;
92 fpSafetyHelper = transportMgr->GetSafetyHelper();
100 fCurrentTouchableHandle = nullTouchableHandle;
102 fEndGlobalTimeComputed =
false;
103 fCandidateEndGlobalTime = 0;
108 G4MonopoleTransportation::~G4MonopoleTransportation()
110 if( (fVerboseLevel > 0) && (fSumEnergyKilled > 0.0 ) ){
111 G4cout <<
" G4MonopoleTransportation: Statistics for looping particles " << G4endl;
112 G4cout <<
" Sum of energy of loopers killed: " << fSumEnergyKilled << G4endl;
113 G4cout <<
" Max energy of loopers killed: " << fMaxEnergyKilled << G4endl;
124 G4double G4MonopoleTransportation::
125 AlongStepGetPhysicalInteractionLength(
const G4Track& track,
127 G4double currentMinimumStep,
128 G4double& currentSafety,
134 G4FieldManager* fieldMgrX=fFieldPropagator->FindAndSetFieldManager( track.GetVolume() );
135 fFieldBuilder->setStepperAndChordFinder (fieldMgrX, 1);
137 G4double geometryStepLength, newSafety ;
138 fParticleIsLooping =
false ;
149 *selection = CandidateForSelection ;
153 const G4DynamicParticle* pParticle = track.GetDynamicParticle() ;
154 G4ThreeVector startMomentumDir = pParticle->GetMomentumDirection() ;
155 G4ThreeVector startPosition = track.GetPosition() ;
163 G4ThreeVector OriginShift = startPosition - fPreviousSftOrigin ;
164 G4double MagSqShift = OriginShift.mag2() ;
165 if( MagSqShift >=
sqr(fPreviousSafety) )
167 currentSafety = 0.0 ;
171 currentSafety = fPreviousSafety -
std::sqrt(MagSqShift) ;
176 G4double particleMagneticCharge = pParticleDef->MagneticCharge() ;
177 G4double particleElectricCharge = pParticle->GetCharge();
179 fGeometryLimitedStep =
false ;
188 G4FieldManager* fieldMgr=0;
189 G4bool fieldExertsForce =
false ;
191 if( (particleMagneticCharge != 0.0) )
193 fieldMgr= fFieldPropagator->FindAndSetFieldManager( track.GetVolume() );
196 fieldMgr->ConfigureForTrack( &track );
202 fieldExertsForce = (fieldMgr->GetDetectorField() != 0);
211 if( !fieldExertsForce )
213 G4double linearStepLength ;
214 if( fShortStepOptimisation && (currentMinimumStep <= currentSafety) )
218 geometryStepLength = currentMinimumStep ;
219 fGeometryLimitedStep =
false ;
225 linearStepLength = fLinearNavigator->ComputeStep( startPosition,
231 fPreviousSftOrigin = startPosition ;
232 fPreviousSafety = newSafety ;
237 currentSafety = newSafety ;
239 fGeometryLimitedStep= (linearStepLength <= currentMinimumStep);
240 if( fGeometryLimitedStep )
243 geometryStepLength = linearStepLength ;
248 geometryStepLength = currentMinimumStep ;
251 endpointDistance = geometryStepLength ;
255 fTransportEndPosition = startPosition+geometryStepLength*startMomentumDir ;
259 fTransportEndMomentumDir = startMomentumDir ;
260 fTransportEndKineticEnergy = track.GetKineticEnergy() ;
261 fTransportEndSpin = track.GetPolarization();
262 fParticleIsLooping =
false ;
263 fMomentumChanged =
false ;
264 fEndGlobalTimeComputed =
false ;
269 G4ThreeVector EndUnitMomentum ;
270 G4double lengthAlongCurve ;
271 G4double restMass = pParticleDef->GetPDGMass() ;
273 fFieldPropagator->SetChargeMomentumMass( particleMagneticCharge,
274 particleElectricCharge,
281 G4ThreeVector
spin = track.GetPolarization() ;
282 G4FieldTrack aFieldTrack = G4FieldTrack( startPosition,
283 track.GetMomentumDirection(),
285 track.GetKineticEnergy(),
288 track.GetGlobalTime(),
289 track.GetProperTime(),
291 if( currentMinimumStep > 0 )
295 lengthAlongCurve = fFieldPropagator->ComputeStep( aFieldTrack,
298 track.GetVolume() ) ;
299 fGeometryLimitedStep= lengthAlongCurve < currentMinimumStep;
300 if( fGeometryLimitedStep ) {
301 geometryStepLength = lengthAlongCurve ;
303 geometryStepLength = currentMinimumStep ;
308 geometryStepLength = lengthAlongCurve= 0.0 ;
309 fGeometryLimitedStep =
false ;
314 fPreviousSftOrigin = startPosition ;
315 fPreviousSafety = currentSafety ;
320 fTransportEndPosition = aFieldTrack.GetPosition() ;
324 fMomentumChanged =
true ;
325 fTransportEndMomentumDir = aFieldTrack.GetMomentumDir() ;
327 fTransportEndKineticEnergy = aFieldTrack.GetKineticEnergy() ;
334 fCandidateEndGlobalTime = aFieldTrack.GetLabTimeOfFlight();
335 fEndGlobalTimeComputed =
true;
402 fTransportEndSpin = aFieldTrack.GetSpin();
403 fParticleIsLooping = fFieldPropagator->IsParticleLooping() ;
404 endpointDistance = (fTransportEndPosition - startPosition).
mag() ;
410 if( currentMinimumStep == 0.0 )
412 if( currentSafety == 0.0 ) fGeometryLimitedStep =
true ;
418 if( currentSafety < endpointDistance )
423 if( particleMagneticCharge != 0.0 ) {
426 fLinearNavigator->ComputeSafety( fTransportEndPosition) ;
427 currentSafety = endSafety ;
428 fPreviousSftOrigin = fTransportEndPosition ;
429 fPreviousSafety = currentSafety ;
430 fpSafetyHelper->SetCurrentSafety( currentSafety, fTransportEndPosition);
435 currentSafety += endpointDistance ;
437 #ifdef G4DEBUG_TRANSPORT
438 G4cout.precision(12) ;
439 G4cout <<
"***G4MonopoleTransportation::AlongStepGPIL ** " << G4endl ;
440 G4cout <<
" Called Navigator->ComputeSafety at " << fTransportEndPosition
441 <<
" and it returned safety= " << endSafety << G4endl ;
442 G4cout <<
" Adding endpoint distance " << endpointDistance
443 <<
" to obtain pseudo-safety= " << currentSafety << G4endl ;
448 fParticleChange.ProposeTrueStepLength(geometryStepLength) ;
450 return geometryStepLength ;
458 G4VParticleChange* G4MonopoleTransportation::AlongStepDoIt(
const G4Track& track,
459 const G4Step& stepData )
461 static const G4ParticleDefinition* fOpticalPhoton =
462 G4ParticleTable::GetParticleTable()->FindParticle(
"opticalphoton");
465 static G4int noCalls=0;
469 fParticleChange.Initialize(track) ;
473 fParticleChange.ProposePosition(fTransportEndPosition) ;
474 fParticleChange.ProposeMomentumDirection(fTransportEndMomentumDir) ;
475 fParticleChange.ProposeEnergy(fTransportEndKineticEnergy) ;
476 fParticleChange.SetMomentumChanged(fMomentumChanged) ;
478 fParticleChange.ProposePolarization(fTransportEndSpin);
480 G4double deltaTime = 0.0 ;
486 G4double startTime = track.GetGlobalTime() ;
488 if (!fEndGlobalTimeComputed)
492 G4double finalVelocity = track.GetVelocity() ;
493 G4double initialVelocity = stepData.GetPreStepPoint()->GetVelocity() ;
494 G4double stepLength = track.GetStepLength() ;
497 const G4DynamicParticle* fpDynamicParticle = track.GetDynamicParticle();
498 if (fpDynamicParticle->GetDefinition()== fOpticalPhoton)
502 deltaTime = stepLength/finalVelocity ;
504 else if (finalVelocity > 0.0)
506 G4double meanInverseVelocity ;
508 meanInverseVelocity = 0.5
509 * ( 1.0 / initialVelocity + 1.0 / finalVelocity ) ;
510 deltaTime = stepLength * meanInverseVelocity ;
512 else if( initialVelocity > 0.0 )
514 deltaTime = stepLength/initialVelocity ;
516 fCandidateEndGlobalTime = startTime + deltaTime ;
520 deltaTime = fCandidateEndGlobalTime - startTime ;
523 fParticleChange.ProposeGlobalTime( fCandidateEndGlobalTime ) ;
527 G4double restMass = track.GetDynamicParticle()->GetMass() ;
528 G4double deltaProperTime = deltaTime*( restMass/track.GetTotalEnergy() ) ;
530 fParticleChange.ProposeProperTime(track.GetProperTime() + deltaProperTime) ;
537 if ( fParticleIsLooping )
539 G4double endEnergy= fTransportEndKineticEnergy;
541 if( (endEnergy < fThreshold_Important_Energy)
542 || (fNoLooperTrials >= fThresholdTrials ) ){
545 fParticleChange.ProposeTrackStatus( fStopAndKill ) ;
548 fSumEnergyKilled += endEnergy;
549 if( endEnergy > fMaxEnergyKilled) { fMaxEnergyKilled= endEnergy; }
552 if( (fVerboseLevel > 1) ||
553 ( endEnergy > fThreshold_Warning_Energy ) ) {
554 G4cout <<
" G4MonopoleTransportation is killing track that is looping or stuck "
556 <<
" This track has " << track.GetKineticEnergy() / MeV
557 <<
" MeV energy." << G4endl;
558 G4cout <<
" Number of trials = " << fNoLooperTrials
559 <<
" No of calls to AlongStepDoIt = " << noCalls
568 if( (fVerboseLevel > 2) ){
569 G4cout <<
" G4MonopoleTransportation::AlongStepDoIt(): Particle looping - "
570 <<
" Number of trials = " << fNoLooperTrials
571 <<
" No of calls to = " << noCalls
585 fParticleChange.SetPointerToVectorOfAuxiliaryPoints
586 (fFieldPropagator->GimmeTrajectoryVectorAndForgetIt() );
588 return &fParticleChange ;
597 G4double G4MonopoleTransportation::
598 PostStepGetPhysicalInteractionLength(
const G4Track& track,
600 G4ForceCondition* pForceCond )
604 G4FieldManager* fieldMgr=fFieldPropagator->FindAndSetFieldManager( track.GetVolume() );
605 fFieldBuilder->setStepperAndChordFinder (fieldMgr, 0);
607 *pForceCond = Forced ;
614 G4VParticleChange* G4MonopoleTransportation::PostStepDoIt(
const G4Track& track,
617 G4TouchableHandle retCurrentTouchable ;
623 fParticleChange.ProposeTrackStatus(track.GetTrackStatus()) ;
628 if(fGeometryLimitedStep)
634 fLinearNavigator->SetGeometricallyLimitedStep() ;
636 LocateGlobalPointAndUpdateTouchableHandle( track.GetPosition(),
637 track.GetMomentumDirection(),
638 fCurrentTouchableHandle,
643 if( fCurrentTouchableHandle->GetVolume() == 0 )
645 fParticleChange.ProposeTrackStatus( fStopAndKill ) ;
647 retCurrentTouchable = fCurrentTouchableHandle ;
648 fParticleChange.SetTouchableHandle( fCurrentTouchableHandle ) ;
654 fLinearNavigator->LocateGlobalPointWithinVolume( track.GetPosition() ) ;
661 fParticleChange.SetTouchableHandle( track.GetTouchableHandle() ) ;
662 retCurrentTouchable = track.GetTouchableHandle() ;
665 const G4VPhysicalVolume* pNewVol = retCurrentTouchable->GetVolume() ;
666 const G4Material* pNewMaterial = 0 ;
667 const G4VSensitiveDetector* pNewSensitiveDetector = 0 ;
671 pNewMaterial= pNewVol->GetLogicalVolume()->GetMaterial();
672 pNewSensitiveDetector= pNewVol->GetLogicalVolume()->GetSensitiveDetector();
678 fParticleChange.SetMaterialInTouchable( (G4Material *) pNewMaterial ) ;
679 fParticleChange.SetSensitiveDetectorInTouchable( (G4VSensitiveDetector *) pNewSensitiveDetector ) ;
681 const G4MaterialCutsCouple* pNewMaterialCutsCouple = 0;
684 pNewMaterialCutsCouple=pNewVol->GetLogicalVolume()->GetMaterialCutsCouple();
687 if( pNewVol!=0 && pNewMaterialCutsCouple!=0 && pNewMaterialCutsCouple->GetMaterial()!=pNewMaterial )
691 pNewMaterialCutsCouple =
692 G4ProductionCutsTable::GetProductionCutsTable()
693 ->GetMaterialCutsCouple(pNewMaterial,
694 pNewMaterialCutsCouple->GetProductionCuts());
696 fParticleChange.SetMaterialCutsCoupleInTouchable( pNewMaterialCutsCouple );
704 fParticleChange.SetTouchableHandle(retCurrentTouchable) ;
706 return &fParticleChange ;
713 G4MonopoleTransportation::StartTracking(G4Track* aTrack)
715 G4VProcess::StartTracking(aTrack);
722 fPreviousSafety = 0.0 ;
723 fPreviousSftOrigin = G4ThreeVector(0.,0.,0.) ;
733 if( DoesGlobalFieldExist() ) {
734 fFieldPropagator->ClearPropagatorState();
743 static G4FieldManagerStore* fieldMgrStore= G4FieldManagerStore::GetInstance();
744 fieldMgrStore->ClearAllChordFindersState();
748 fCurrentTouchableHandle = aTrack->GetTouchableHandle();
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
Square< F >::type sqr(const F &f)
volatile std::atomic< bool > shutdown_flag false