CMS 3D CMS Logo

MonopoleTransportation.cc
Go to the documentation of this file.
1 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
2 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
3 //
4 // This class is a process responsible for the transportation of
5 // magnetic monopoles, ie the geometrical propagation that encounters the
6 // geometrical sub-volumes of the detectors.
7 //
8 // =======================================================================
9 // Created: 3 May 2010, J. Apostolakis, B. Bozsogi
10 // G4MonopoleTransportation class for
11 // Geant4 extended example "monopole"
12 //
13 // Adopted for CMSSW by V.Ivanchenko 30 April 2018
14 //
15 // =======================================================================
16 
21 
22 #include "G4ProductionCutsTable.hh"
23 #include "G4ParticleTable.hh"
24 #include "G4ChordFinder.hh"
25 #include "G4SafetyHelper.hh"
26 #include "G4FieldManagerStore.hh"
27 #include "G4TransportationProcessType.hh"
28 #include "G4SystemOfUnits.hh"
29 
30 class G4VSensitiveDetector;
31 
32 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
33 
35  sim::ChordFinderSetter* chordFinderSetter,
36  G4int verb)
37  : G4VProcess( G4String("MonopoleTransportation"), fTransportation ),
38  fParticleDef(mpl),
39  fChordFinderSetter(chordFinderSetter),
40  fieldMgrCMS(nullptr),
41  fLinearNavigator(nullptr),
42  fFieldPropagator(nullptr),
43  fParticleIsLooping( false ),
44  fPreviousSftOrigin (0.,0.,0.),
45  fPreviousSafety ( 0.0 ),
46  fThreshold_Warning_Energy( 100 * MeV ),
47  fThreshold_Important_Energy( 250 * MeV ),
48  fThresholdTrials( 10 ),
49  fNoLooperTrials(0),
50  fSumEnergyKilled( 0.0 ), fMaxEnergyKilled( 0.0 ),
51  fShortStepOptimisation(false), // Old default: true (=fast short steps)
52  fpSafetyHelper(nullptr)
53 {
54  verboseLevel = verb;
55 
56  // set Process Sub Type
57  SetProcessSubType(TRANSPORTATION);
58 
59 #ifdef G4MULTITHREADED
60  // Do not finalize the MonopoleTransportation class
61  if (G4Threading::IsMasterThread())
62  {
63  return;
64  }
65 #endif
66 
67  G4TransportationManager* transportMgr =
68  G4TransportationManager::GetTransportationManager();
69 
70  fLinearNavigator = transportMgr->GetNavigatorForTracking() ;
71 
72  fFieldPropagator = transportMgr->GetPropagatorInField() ;
73  fpSafetyHelper = transportMgr->GetSafetyHelper();
74 
75  // Cannot determine whether a field exists here,
76  // because it would only work if the field manager has informed
77  // about the detector's field before this transportation process
78  // is constructed.
79  // Instead later the method DoesGlobalFieldExist() is called
80 
81  static G4TouchableHandle nullTouchableHandle; // Points to (G4VTouchable*) 0
82  fCurrentTouchableHandle = nullTouchableHandle;
83 
84  fEndGlobalTimeComputed = false;
86 }
87 
88 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
89 
91 {
92  if( (verboseLevel > 0) && (fSumEnergyKilled > 0.0 ) ){
93  /*
94  G4cout << " MonopoleTransportation: Statistics for looping particles "
95  << G4endl;
96  G4cout << " Sum of energy of loopers killed: " << fSumEnergyKilled << G4endl;
97  G4cout << " Max energy of loopers killed: " << fMaxEnergyKilled << G4endl;
98  */
99  }
100 }
101 
102 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
103 //
104 // Responsibilities:
105 // Find whether the geometry limits the Step, and to what length
106 // Calculate the new value of the safety and return it.
107 // Store the final time, position and momentum.
108 
111  G4double, // previousStepSize
112  G4double currentMinimumStep,
113  G4double& currentSafety,
114  G4GPILSelection* selection )
115 {
116  // change to monopole equation
118 
119  G4double geometryStepLength, newSafety ;
121 
122  // Initial actions moved to StartTrack()
123  // --------------------------------------
124  // Note: in case another process changes touchable handle
125  // it will be necessary to add here (for all steps)
126  // fCurrentTouchableHandle = aTrack->GetTouchableHandle();
127 
128  // GPILSelection is set to defaule value of CandidateForSelection
129  // It is a return value
130  //
131  *selection = CandidateForSelection ;
132 
133  // Get initial Energy/Momentum of the track
134  //
135  const G4DynamicParticle* pParticle = track.GetDynamicParticle() ;
136  const G4ThreeVector& startMomentumDir = pParticle->GetMomentumDirection() ;
137  G4ThreeVector startPosition = track.GetPosition() ;
138 
139  // The Step Point safety can be limited by other geometries and/or the
140  // assumptions of any process - it's not always the geometrical safety.
141  // We calculate the starting point's isotropic safety here.
142  //
143  G4ThreeVector OriginShift = startPosition - fPreviousSftOrigin ;
144  G4double MagSqShift = OriginShift.mag2() ;
145  if( MagSqShift >= sqr(fPreviousSafety) )
146  {
147  currentSafety = 0.0 ;
148  }
149  else
150  {
151  currentSafety = fPreviousSafety - std::sqrt(MagSqShift) ;
152  }
153 
154  // Is the monopole charged ?
155  //
156  G4double particleMagneticCharge = fParticleDef->MagneticCharge() ;
157  G4double particleElectricCharge = pParticle->GetCharge();
158 
160 
161  // There is no need to locate the current volume. It is Done elsewhere:
162  // On track construction
163  // By the tracking, after all AlongStepDoIts, in "Relocation"
164 
165  // Check whether the particle have an (EM) field force exerting upon it
166  G4bool fieldExertsForce = false ;
167 
168  if(particleMagneticCharge != 0.0 && fieldMgrCMS) {
169  // Message the field Manager, to configure it for this track
170  fieldMgrCMS->ConfigureForTrack( &track );
171  // Moved here, in order to allow a transition
172  // from a zero-field status (with fieldMgr->(field)0
173  // to a finite field status
174 
175  // If the field manager has no field, there is no field !
176  fieldExertsForce = (fieldMgrCMS->GetDetectorField() != nullptr);
177  }
178 
179  // G4cout << " G4Transport: field exerts force= " << fieldExertsForce
180  // << " fieldMgr= " << fieldMgr << G4endl;
181 
182  // Choose the calculation of the transportation: Field or not
183  //
184  if( !fieldExertsForce )
185  {
186  G4double linearStepLength ;
187  if( fShortStepOptimisation && (currentMinimumStep <= currentSafety) )
188  {
189  // The Step is guaranteed to be taken
190  //
191  geometryStepLength = currentMinimumStep ;
193  }
194  else
195  {
196  // Find whether the straight path intersects a volume
197  //
198  linearStepLength = fLinearNavigator->ComputeStep( startPosition,
199  startMomentumDir,
200  currentMinimumStep,
201  newSafety) ;
202  // Remember last safety origin & value.
203  //
204  fPreviousSftOrigin = startPosition ;
205  fPreviousSafety = newSafety ;
206  // fpSafetyHelper->SetCurrentSafety( newSafety, startPosition);
207 
208  // The safety at the initial point has been re-calculated:
209  //
210  currentSafety = newSafety ;
211 
212  fGeometryLimitedStep= (linearStepLength <= currentMinimumStep);
214  {
215  // The geometry limits the Step size (an intersection was found.)
216  geometryStepLength = linearStepLength ;
217  }
218  else
219  {
220  // The full Step is taken.
221  geometryStepLength = currentMinimumStep ;
222  }
223  }
224  endpointDistance = geometryStepLength ;
225 
226  // Calculate final position
227  //
228  fTransportEndPosition = startPosition+geometryStepLength*startMomentumDir ;
229 
230  // Momentum direction, energy and polarisation are unchanged by transport
231  //
232  fTransportEndMomentumDir = startMomentumDir ;
233  fTransportEndKineticEnergy = track.GetKineticEnergy() ;
234  fTransportEndSpin = track.GetPolarization();
238  }
239  else // A field exerts force
240  {
241  G4double momentumMagnitude = pParticle->GetTotalMomentum() ;
242  G4ThreeVector EndUnitMomentum ;
243  G4double lengthAlongCurve ;
244  G4double restMass = fParticleDef->GetPDGMass() ;
245 
246  G4ChargeState chargeState(particleElectricCharge, // The charge can change (dynamic)
247  fParticleDef->GetPDGSpin(),
248  0, // Magnetic moment: pParticleDef->GetMagneticMoment(),
249  0, // Electric Dipole moment - not in Particle Definition
250  particleMagneticCharge ); // in Mev/c
251 
252  G4EquationOfMotion* equationOfMotion =
253  fFieldPropagator->GetChordFinder()->GetIntegrationDriver()->GetEquationOfMotion();
254 
255  equationOfMotion
256  ->SetChargeMomentumMass( chargeState, // Was particleMagneticCharge - in Mev/c
257  momentumMagnitude, // Was particleElectricCharge
258  restMass ) ;
259  // SetChargeMomentumMass now passes both the electric and magnetic charge - in chargeState
260 
261  G4ThreeVector spin = track.GetPolarization() ;
262  G4FieldTrack aFieldTrack = G4FieldTrack( startPosition,
263  track.GetMomentumDirection(),
264  0.0,
265  track.GetKineticEnergy(),
266  restMass,
267  track.GetVelocity(),
268  track.GetGlobalTime(), // Lab.
269  track.GetProperTime(), // Part.
270  &spin ) ;
271  if( currentMinimumStep > 0 )
272  {
273  // Do the Transport in the field (non recti-linear)
274  //
275  lengthAlongCurve = fFieldPropagator->ComputeStep( aFieldTrack,
276  currentMinimumStep,
277  currentSafety,
278  track.GetVolume() ) ;
279  fGeometryLimitedStep= lengthAlongCurve < currentMinimumStep;
280  if( fGeometryLimitedStep ) {
281  geometryStepLength = lengthAlongCurve ;
282  } else {
283  geometryStepLength = currentMinimumStep ;
284  }
285  }
286  else
287  {
288  geometryStepLength = lengthAlongCurve= 0.0 ;
290  }
291 
292  // Remember last safety origin & value.
293  //
294  fPreviousSftOrigin = startPosition ;
295  fPreviousSafety = currentSafety ;
296  // fpSafetyHelper->SetCurrentSafety( newSafety, startPosition);
297 
298  // Get the End-Position and End-Momentum (Dir-ection)
299  //
300  fTransportEndPosition = aFieldTrack.GetPosition() ;
301 
302  // Momentum: Magnitude and direction can be changed too now ...
303  //
305  fTransportEndMomentumDir = aFieldTrack.GetMomentumDir() ;
306 
307  fTransportEndKineticEnergy = aFieldTrack.GetKineticEnergy() ;
308 
309  fCandidateEndGlobalTime = aFieldTrack.GetLabTimeOfFlight();
310  fEndGlobalTimeComputed = true;
311 
312  fTransportEndSpin = aFieldTrack.GetSpin();
313  fParticleIsLooping = fFieldPropagator->IsParticleLooping() ;
314  endpointDistance = (fTransportEndPosition - startPosition).mag() ;
315  }
316 
317  // If we are asked to go a step length of 0, and we are on a boundary
318  // then a boundary will also limit the step -> we must flag this.
319  //
320  if( currentMinimumStep == 0.0 )
321  {
322  if( currentSafety == 0.0 ) fGeometryLimitedStep = true ;
323  }
324 
325  // Update the safety starting from the end-point,
326  // if it will become negative at the end-point.
327  //
328  if( currentSafety < endpointDistance )
329  {
330  // if( particleMagneticCharge == 0.0 )
331  // G4cout << " Avoiding call to ComputeSafety : charge = 0.0 " << G4endl;
332 
333  if( particleMagneticCharge != 0.0 ) {
334 
335  G4double endSafety =
336  fLinearNavigator->ComputeSafety( fTransportEndPosition) ;
337  currentSafety = endSafety ;
338  fPreviousSftOrigin = fTransportEndPosition ;
339  fPreviousSafety = currentSafety ;
340  fpSafetyHelper->SetCurrentSafety( currentSafety, fTransportEndPosition);
341 
342  // Because the Stepping Manager assumes it is from the start point,
343  // add the StepLength
344  //
345  currentSafety += endpointDistance ;
346 
347 #ifdef G4DEBUG_TRANSPORT
348  G4cout.precision(12) ;
349  G4cout << "***MonopoleTransportation::AlongStepGPIL ** " << G4endl ;
350  G4cout << " Called Navigator->ComputeSafety at " << fTransportEndPosition
351  << " and it returned safety= " << endSafety << G4endl ;
352  G4cout << " Adding endpoint distance " << endpointDistance
353  << " to obtain pseudo-safety= " << currentSafety << G4endl ;
354 #endif
355  }
356  }
357 
358  fParticleChange.ProposeTrueStepLength(geometryStepLength) ;
359 
360  return geometryStepLength ;
361 }
362 
363 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
364 //
365 // Initialize ParticleChange (by setting all its members equal
366 // to corresponding members in G4Track)
367 
368 G4VParticleChange* MonopoleTransportation::AlongStepDoIt( const G4Track& track,
369  const G4Step& stepData )
370 {
371  fParticleChange.Initialize(track) ;
372 
373  // Code for specific process
374  //
375  fParticleChange.ProposePosition(fTransportEndPosition) ;
376  fParticleChange.ProposeMomentumDirection(fTransportEndMomentumDir) ;
378  fParticleChange.SetMomentumChanged(fMomentumChanged) ;
379 
380  fParticleChange.ProposePolarization(fTransportEndSpin);
381 
382  G4double deltaTime = 0.0 ;
383 
384  // Calculate Lab Time of Flight (ONLY if field Equations used it!)
385  // G4double endTime = fCandidateEndGlobalTime;
386  // G4double delta_time = endTime - startTime;
387 
388  G4double startTime = track.GetGlobalTime() ;
389 
391  {
392  // The time was not integrated .. make the best estimate possible
393  //
394  G4double finalVelocity = track.GetVelocity() ;
395  G4double initialVelocity = stepData.GetPreStepPoint()->GetVelocity() ;
396  G4double stepLength = track.GetStepLength() ;
397 
398  deltaTime= 0.0; // in case initialVelocity = 0
399  //const G4DynamicParticle* fpDynamicParticle = track.GetDynamicParticle();
400  if (finalVelocity > 0.0)
401  {
402  G4double meanInverseVelocity ;
403  meanInverseVelocity = 0.5
404  * ( 1.0 / initialVelocity + 1.0 / finalVelocity ) ;
405  deltaTime = stepLength * meanInverseVelocity ;
406  }
407  else if( initialVelocity > 0.0 )
408  {
409  deltaTime = stepLength/initialVelocity ;
410  }
411  fCandidateEndGlobalTime = startTime + deltaTime ;
412  }
413  else
414  {
415  deltaTime = fCandidateEndGlobalTime - startTime ;
416  }
417 
418  fParticleChange.ProposeGlobalTime( fCandidateEndGlobalTime ) ;
419 
420  // Now Correct by Lorentz factor to get "proper" deltaTime
421 
422  G4double restMass = track.GetDynamicParticle()->GetMass() ;
423  G4double deltaProperTime = deltaTime*( restMass/track.GetTotalEnergy() ) ;
424 
425  fParticleChange.ProposeProperTime(track.GetProperTime() + deltaProperTime) ;
426  //fParticleChange. ProposeTrueStepLength( track.GetStepLength() ) ;
427 
428  // If the particle is caught looping or is stuck (in very difficult
429  // boundaries) in a magnetic field (doing many steps)
430  // THEN this kills it ...
431  //
432  if ( fParticleIsLooping )
433  {
434  G4double endEnergy= fTransportEndKineticEnergy;
435 
436  if( (endEnergy < fThreshold_Important_Energy)
438  // Kill the looping particle
439  //
440  fParticleChange.ProposeTrackStatus( fStopAndKill ) ;
441 
442  // 'Bare' statistics
443  fSumEnergyKilled += endEnergy;
444  if( endEnergy > fMaxEnergyKilled) { fMaxEnergyKilled= endEnergy; }
445 
446 #ifdef G4VERBOSE
447  if( (verboseLevel > 1) ||
448  ( endEnergy > fThreshold_Warning_Energy ) ) {
449  G4cout << " MonopoleTransportation is killing track that is looping or stuck "
450  << G4endl
451  << " This track has " << track.GetKineticEnergy() / MeV
452  << " MeV energy." << G4endl;
453  G4cout << " Number of trials = " << fNoLooperTrials
454  << " No of calls to AlongStepDoIt = " << noCalls
455  << G4endl;
456  }
457 #endif
458  fNoLooperTrials=0;
459  }
460  else{
461  fNoLooperTrials ++;
462 #ifdef G4VERBOSE
463  if( (verboseLevel > 2) ){
464  G4cout << " MonopoleTransportation::AlongStepDoIt(): Particle looping - "
465  << " Number of trials = " << fNoLooperTrials
466  << " No of calls to = " << noCalls
467  << G4endl;
468  }
469 #endif
470  }
471  }else{
472  fNoLooperTrials=0;
473  }
474 
475  // Another (sometimes better way) is to use a user-limit maximum Step size
476  // to alleviate this problem ..
477 
478  // Introduce smooth curved trajectories to particle-change
479  //
480  fParticleChange.SetPointerToVectorOfAuxiliaryPoints
481  (fFieldPropagator->GimmeTrajectoryVectorAndForgetIt() );
482 
483  return &fParticleChange ;
484 }
485 
486 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
487 //
488 // This ensures that the PostStep action is always called,
489 // so that it can do the relocation if it is needed.
490 //
491 
494  G4double, // previousStepSize
495  G4ForceCondition* pForceCond )
496 {
497  *pForceCond = Forced ;
498  return DBL_MAX ; // was kInfinity ; but convention now is DBL_MAX
499 }
500 
501 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
502 
503 G4VParticleChange* MonopoleTransportation::PostStepDoIt( const G4Track& track,
504  const G4Step& )
505 {
506  G4TouchableHandle retCurrentTouchable ; // The one to return
507 
508  fParticleChange.ProposeTrackStatus(track.GetTrackStatus()) ;
509 
510  // If the Step was determined by the volume boundary,
511  // logically relocate the particle
512 
514  {
515  // fCurrentTouchable will now become the previous touchable,
516  // and what was the previous will be freed.
517  // (Needed because the preStepPoint can point to the previous touchable)
518 
519  fLinearNavigator->SetGeometricallyLimitedStep() ;
521  LocateGlobalPointAndUpdateTouchableHandle( track.GetPosition(),
522  track.GetMomentumDirection(),
524  true ) ;
525  // Check whether the particle is out of the world volume
526  // If so it has exited and must be killed.
527  //
528  if( fCurrentTouchableHandle->GetVolume() == nullptr )
529  {
530  fParticleChange.ProposeTrackStatus( fStopAndKill ) ;
531  }
532  retCurrentTouchable = fCurrentTouchableHandle ;
533  fParticleChange.SetTouchableHandle( fCurrentTouchableHandle ) ;
534  }
535  else // fGeometryLimitedStep is false
536  {
537  // This serves only to move the Navigator's location
538  //
539  fLinearNavigator->LocateGlobalPointWithinVolume( track.GetPosition() ) ;
540 
541  // The value of the track's current Touchable is retained.
542  // (and it must be correct because we must use it below to
543  // overwrite the (unset) one in particle change)
544  // It must be fCurrentTouchable too ??
545  //
546  fParticleChange.SetTouchableHandle( track.GetTouchableHandle() ) ;
547  retCurrentTouchable = track.GetTouchableHandle() ;
548  } // endif ( fGeometryLimitedStep )
549 
550  const G4VPhysicalVolume* pNewVol = retCurrentTouchable->GetVolume() ;
551  const G4Material* pNewMaterial = nullptr ;
552  const G4VSensitiveDetector* pNewSensitiveDetector = nullptr ;
553 
554  if( pNewVol != nullptr )
555  {
556  pNewMaterial= pNewVol->GetLogicalVolume()->GetMaterial();
557  pNewSensitiveDetector= pNewVol->GetLogicalVolume()->GetSensitiveDetector();
558  }
559 
560  // ( <const_cast> pNewMaterial ) ;
561  // ( <const_cast> pNewSensitiveDetector) ;
562 
563  fParticleChange.SetMaterialInTouchable(
564  (G4Material *) pNewMaterial ) ;
565  fParticleChange.SetSensitiveDetectorInTouchable(
566  (G4VSensitiveDetector *) pNewSensitiveDetector ) ;
567 
568  const G4MaterialCutsCouple* pNewMaterialCutsCouple = nullptr;
569  if( pNewVol != nullptr )
570  {
571  pNewMaterialCutsCouple=pNewVol->GetLogicalVolume()->GetMaterialCutsCouple();
572  }
573 
574  if( pNewVol!=nullptr && pNewMaterialCutsCouple!=nullptr &&
575  pNewMaterialCutsCouple->GetMaterial()!=pNewMaterial )
576  {
577  // for parametrized volume
578  //
579  pNewMaterialCutsCouple =
580  G4ProductionCutsTable::GetProductionCutsTable()
581  ->GetMaterialCutsCouple(pNewMaterial,
582  pNewMaterialCutsCouple->GetProductionCuts());
583  }
584  fParticleChange.SetMaterialCutsCoupleInTouchable( pNewMaterialCutsCouple );
585 
586  // temporarily until Get/Set Material of ParticleChange,
587  // and StepPoint can be made const.
588  // Set the touchable in ParticleChange
589  // this must always be done because the particle change always
590  // uses this value to overwrite the current touchable pointer.
591  //
592  fParticleChange.SetTouchableHandle(retCurrentTouchable) ;
593 
594  // change to normal equation
596 
597  return &fParticleChange ;
598 }
599 
600 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
601 
602 // New method takes over the responsibility to reset the state
603 // of MonopoleTransportation object at the start of a new track
604 // or the resumption of a suspended track.
605 
606 void
608 {
609  // initialise pointer
610  if(!fieldMgrCMS) {
611  G4FieldManager* fieldMgr =
612  fFieldPropagator->FindAndSetFieldManager(aTrack->GetVolume());
613  fieldMgrCMS = static_cast<CMSFieldManager*>(fieldMgr);
614  }
615 
616  G4VProcess::StartTracking(aTrack);
617 
618  // The actions here are those that were taken in AlongStepGPIL
619  // when track.GetCurrentStepNumber()==1
620 
621  // reset safety value and center
622  //
623  fPreviousSafety = 0.0 ;
624  fPreviousSftOrigin = G4ThreeVector(0.,0.,0.) ;
625 
626  // reset looping counter -- for motion in field
627  fNoLooperTrials= 0;
628  // Must clear this state .. else it depends on last track's value
629  // --> a better solution would set this from state of suspended track TODO ?
630  // Was if( aTrack->GetCurrentStepNumber()==1 ) { .. }
631 
632  // ChordFinder reset internal state
633  //
634  if( DoesGlobalFieldExist() ) {
635  fFieldPropagator->ClearPropagatorState();
636  // Resets all state of field propagator class (ONLY)
637  // including safety values (in case of overlaps and to wipe for first track).
638 
639  G4ChordFinder* chordF= fFieldPropagator->GetChordFinder();
640  if( chordF ) chordF->ResetStepEstimate();
641  }
642 
643  // Make sure to clear the chord finders of all fields (ie managers)
644  G4FieldManagerStore* fieldMgrStore= G4FieldManagerStore::GetInstance();
645  fieldMgrStore->ClearAllChordFindersState();
646 
647  // Update the current touchable handle (from the track's)
648  //
649  fCurrentTouchableHandle = aTrack->GetTouchableHandle();
650 }
651 
652 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
653 
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
void StartTracking(G4Track *aTrack) override
selection
main part
Definition: corrVsCorr.py:98
MonopoleTransportation(const Monopole *p, sim::ChordFinderSetter *cfs, G4int verbosityLevel=1)
#define nullptr
const double MeV
T sqrt(T t)
Definition: SSEVec.h:18
void ConfigureForTrack(const G4Track *) override
G4ParticleChangeForTransport fParticleChange
G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection) override
G4double PostStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4ForceCondition *pForceCond) override
G4PropagatorInField * fFieldPropagator
G4TouchableHandle fCurrentTouchableHandle
G4double MagneticCharge() const
Definition: Monopole.h:19
Square< F >::type sqr(const F &f)
Definition: Square.h:13
G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData) override
G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData) override
void SetMonopoleTracking(G4bool)
float spin(float ph)