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