CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros 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 //
26 // GEANT4 tag $Name: $
27 //
28 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
29 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
30 //
31 // This class is a process responsible for the transportation of
32 // magnetic monopoles, ie the geometrical propagation that encounters the
33 // geometrical sub-volumes of the detectors.
34 //
35 // For monopoles, uses a different equation of motion and ignores energy
36 // conservation.
37 //
38 
39 // =======================================================================
40 // Created: 3 May 2010, J. Apostolakis, B. Bozsogi
41 // =======================================================================
42 
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"
50 
51 class G4VSensitiveDetector;
52 
54 //
55 // Constructor
56 
57 namespace {
58  static const G4TouchableHandle nullTouchableHandle; // Points to (G4VTouchable*) 0
59 }
60 
61 G4MonopoleTransportation::G4MonopoleTransportation(const G4Monopole* mpl,
62  sim::FieldBuilder* fieldBuilder,
63  G4int verb)
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 ),
73  fNoLooperTrials(0),
74  fSumEnergyKilled( 0.0 ), fMaxEnergyKilled( 0.0 ),
75  fShortStepOptimisation(false), // Old default: true (=fast short steps)
76  fVerboseLevel( verb )
77 {
78  pParticleDef = mpl;
79 
80  //magSetup = G4MonopoleFieldSetup::GetMonopoleFieldSetup();
81 
82  G4TransportationManager* transportMgr ;
83 
84  transportMgr = G4TransportationManager::GetTransportationManager() ;
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  fCurrentTouchableHandle = nullTouchableHandle;
101 
102  fEndGlobalTimeComputed = false;
103  fCandidateEndGlobalTime = 0;
104 }
105 
107 
108 G4MonopoleTransportation::~G4MonopoleTransportation()
109 {
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;
114  }
115 }
116 
118 //
119 // Responsibilities:
120 // Find whether the geometry limits the Step, and to what length
121 // Calculate the new value of the safety and return it.
122 // Store the final time, position and momentum.
123 
124 G4double G4MonopoleTransportation::
125 AlongStepGetPhysicalInteractionLength( const G4Track& track,
126  G4double, // previousStepSize
127  G4double currentMinimumStep,
128  G4double& currentSafety,
129  G4GPILSelection* selection )
130 {
131 
132  //magSetup->SetStepperAndChordFinder(1);
133  // change to monopole equation
134  G4FieldManager* fieldMgrX=fFieldPropagator->FindAndSetFieldManager( track.GetVolume() );
135  fFieldBuilder->setStepperAndChordFinder (fieldMgrX, 1);
136 
137  G4double geometryStepLength, newSafety ;
138  fParticleIsLooping = false ;
139 
140  // Initial actions moved to StartTrack()
141  // --------------------------------------
142  // Note: in case another process changes touchable handle
143  // it will be necessary to add here (for all steps)
144  // fCurrentTouchableHandle = aTrack->GetTouchableHandle();
145 
146  // GPILSelection is set to defaule value of CandidateForSelection
147  // It is a return value
148  //
149  *selection = CandidateForSelection ;
150 
151  // Get initial Energy/Momentum of the track
152  //
153  const G4DynamicParticle* pParticle = track.GetDynamicParticle() ;
154  G4ThreeVector startMomentumDir = pParticle->GetMomentumDirection() ;
155  G4ThreeVector startPosition = track.GetPosition() ;
156 
157  // G4double theTime = track.GetGlobalTime() ;
158 
159  // The Step Point safety can be limited by other geometries and/or the
160  // assumptions of any process - it's not always the geometrical safety.
161  // We calculate the starting point's isotropic safety here.
162  //
163  G4ThreeVector OriginShift = startPosition - fPreviousSftOrigin ;
164  G4double MagSqShift = OriginShift.mag2() ;
165  if( MagSqShift >= sqr(fPreviousSafety) )
166  {
167  currentSafety = 0.0 ;
168  }
169  else
170  {
171  currentSafety = fPreviousSafety - std::sqrt(MagSqShift) ;
172  }
173 
174  // Is the monopole charged ?
175  //
176  G4double particleMagneticCharge = pParticleDef->MagneticCharge() ;
177  G4double particleElectricCharge = pParticle->GetCharge();
178 
179  fGeometryLimitedStep = false ;
180  // fEndGlobalTimeComputed = false ;
181 
182  // There is no need to locate the current volume. It is Done elsewhere:
183  // On track construction
184  // By the tracking, after all AlongStepDoIts, in "Relocation"
185 
186  // Check whether the particle have an (EM) field force exerting upon it
187  //
188  G4FieldManager* fieldMgr=0;
189  G4bool fieldExertsForce = false ;
190 
191  if( (particleMagneticCharge != 0.0) )
192  {
193  fieldMgr= fFieldPropagator->FindAndSetFieldManager( track.GetVolume() );
194  if (fieldMgr != 0) {
195  // Message the field Manager, to configure it for this track
196  fieldMgr->ConfigureForTrack( &track );
197  // Moved here, in order to allow a transition
198  // from a zero-field status (with fieldMgr->(field)0
199  // to a finite field status
200 
201  // If the field manager has no field, there is no field !
202  fieldExertsForce = (fieldMgr->GetDetectorField() != 0);
203  }
204  }
205 
206  // G4cout << " G4Transport: field exerts force= " << fieldExertsForce
207  // << " fieldMgr= " << fieldMgr << G4endl;
208 
209  // Choose the calculation of the transportation: Field or not
210  //
211  if( !fieldExertsForce )
212  {
213  G4double linearStepLength ;
214  if( fShortStepOptimisation && (currentMinimumStep <= currentSafety) )
215  {
216  // The Step is guaranteed to be taken
217  //
218  geometryStepLength = currentMinimumStep ;
219  fGeometryLimitedStep = false ;
220  }
221  else
222  {
223  // Find whether the straight path intersects a volume
224  //
225  linearStepLength = fLinearNavigator->ComputeStep( startPosition,
226  startMomentumDir,
227  currentMinimumStep,
228  newSafety) ;
229  // Remember last safety origin & value.
230  //
231  fPreviousSftOrigin = startPosition ;
232  fPreviousSafety = newSafety ;
233  // fpSafetyHelper->SetCurrentSafety( newSafety, startPosition);
234 
235  // The safety at the initial point has been re-calculated:
236  //
237  currentSafety = newSafety ;
238 
239  fGeometryLimitedStep= (linearStepLength <= currentMinimumStep);
240  if( fGeometryLimitedStep )
241  {
242  // The geometry limits the Step size (an intersection was found.)
243  geometryStepLength = linearStepLength ;
244  }
245  else
246  {
247  // The full Step is taken.
248  geometryStepLength = currentMinimumStep ;
249  }
250  }
251  endpointDistance = geometryStepLength ;
252 
253  // Calculate final position
254  //
255  fTransportEndPosition = startPosition+geometryStepLength*startMomentumDir ;
256 
257  // Momentum direction, energy and polarisation are unchanged by transport
258  //
259  fTransportEndMomentumDir = startMomentumDir ;
260  fTransportEndKineticEnergy = track.GetKineticEnergy() ;
261  fTransportEndSpin = track.GetPolarization();
262  fParticleIsLooping = false ;
263  fMomentumChanged = false ;
264  fEndGlobalTimeComputed = false ;
265  }
266  else // A field exerts force
267  {
268  // G4double momentumMagnitude = pParticle->GetTotalMomentum() ;
269  G4ThreeVector EndUnitMomentum ;
270  G4double lengthAlongCurve ;
271  G4double restMass = pParticleDef->GetPDGMass() ;
272 
273  fFieldPropagator->SetChargeMomentumMass( particleMagneticCharge, // in Mev/c
274  particleElectricCharge, // in e+ units
275  restMass ) ;
276 
277  // SetChargeMomentumMass is _not_ used here as it would in everywhere else,
278  // it's just a workaround to pass the electric charge as well.
279 
280 
281  G4ThreeVector spin = track.GetPolarization() ;
282  G4FieldTrack aFieldTrack = G4FieldTrack( startPosition,
283  track.GetMomentumDirection(),
284  0.0,
285  track.GetKineticEnergy(),
286  restMass,
287  track.GetVelocity(),
288  track.GetGlobalTime(), // Lab.
289  track.GetProperTime(), // Part.
290  &spin ) ;
291  if( currentMinimumStep > 0 )
292  {
293  // Do the Transport in the field (non recti-linear)
294  //
295  lengthAlongCurve = fFieldPropagator->ComputeStep( aFieldTrack,
296  currentMinimumStep,
297  currentSafety,
298  track.GetVolume() ) ;
299  fGeometryLimitedStep= lengthAlongCurve < currentMinimumStep;
300  if( fGeometryLimitedStep ) {
301  geometryStepLength = lengthAlongCurve ;
302  } else {
303  geometryStepLength = currentMinimumStep ;
304  }
305  }
306  else
307  {
308  geometryStepLength = lengthAlongCurve= 0.0 ;
309  fGeometryLimitedStep = false ;
310  }
311 
312  // Remember last safety origin & value.
313  //
314  fPreviousSftOrigin = startPosition ;
315  fPreviousSafety = currentSafety ;
316  // fpSafetyHelper->SetCurrentSafety( newSafety, startPosition);
317 
318  // Get the End-Position and End-Momentum (Dir-ection)
319  //
320  fTransportEndPosition = aFieldTrack.GetPosition() ;
321 
322  // Momentum: Magnitude and direction can be changed too now ...
323  //
324  fMomentumChanged = true ;
325  fTransportEndMomentumDir = aFieldTrack.GetMomentumDir() ;
326 
327  fTransportEndKineticEnergy = aFieldTrack.GetKineticEnergy() ;
328 
329 // if( fFieldPropagator->GetCurrentFieldManager()->DoesFieldChangeEnergy() )
330 // {
331 // // If the field can change energy, then the time must be integrated
332 // // - so this should have been updated
333 
334  fCandidateEndGlobalTime = aFieldTrack.GetLabTimeOfFlight();
335  fEndGlobalTimeComputed = true;
336 
337 // // was ( fCandidateEndGlobalTime != track.GetGlobalTime() );
338 // // a cleaner way is to have FieldTrack knowing whether time is updated.
339 // }
340 // else
341 // {
342 // // The energy should be unchanged by field transport,
343 // // - so the time changed will be calculated elsewhere
344 // //
345 // fEndGlobalTimeComputed = false;
346 //
347 // // Check that the integration preserved the energy
348 // // - and if not correct this!
349 // G4double startEnergy= track.GetKineticEnergy();
350 // G4double endEnergy= fTransportEndKineticEnergy;
351 //
352 // static G4int no_inexact_steps=0, no_large_ediff;
353 // G4double absEdiff = std::fabs(startEnergy- endEnergy);
354 // if( absEdiff > perMillion * endEnergy )
355 // {
356 // no_inexact_steps++;
357 // // Possible statistics keeping here ...
358 // }
359 // if( fVerboseLevel > 1 )
360 // {
361 // if( std::fabs(startEnergy- endEnergy) > perThousand * endEnergy )
362 // {
363 // static G4int no_warnings= 0, warnModulo=1, moduloFactor= 10;
364 // no_large_ediff ++;
365 // if( (no_large_ediff% warnModulo) == 0 )
366 // {
367 // no_warnings++;
368 // G4cout << "WARNING - G4MonopoleTransportation::AlongStepGetPIL() "
369 // << " Energy change in Step is above 1^-3 relative value. " << G4endl
370 // << " Relative change in 'tracking' step = "
371 // << std::setw(15) << (endEnergy-startEnergy)/startEnergy << G4endl
372 // << " Starting E= " << std::setw(12) << startEnergy / MeV << " MeV " << G4endl
373 // << " Ending E= " << std::setw(12) << endEnergy / MeV << " MeV " << G4endl;
374 // G4cout << " Energy has been corrected -- however, review"
375 // << " field propagation parameters for accuracy." << G4endl;
376 // if( (fVerboseLevel > 2 ) || (no_warnings<4) || (no_large_ediff == warnModulo * moduloFactor) ){
377 // G4cout << " These include EpsilonStepMax(/Min) in G4FieldManager "
378 // << " which determine fractional error per step for integrated quantities. " << G4endl
379 // << " Note also the influence of the permitted number of integration steps."
380 // << G4endl;
381 // }
382 // G4cerr << "ERROR - G4MonopoleTransportation::AlongStepGetPIL()" << G4endl
383 // << " Bad 'endpoint'. Energy change detected"
384 // << " and corrected. "
385 // << " Has occurred already "
386 // << no_large_ediff << " times." << G4endl;
387 // if( no_large_ediff == warnModulo * moduloFactor )
388 // {
389 // warnModulo *= moduloFactor;
390 // }
391 // }
392 // }
393 // } // end of if (fVerboseLevel)
394 //
395 // // Correct the energy for fields that conserve it
396 // // This - hides the integration error
397 // // - but gives a better physical answer
398 // fTransportEndKineticEnergy= track.GetKineticEnergy();
399 // }
400 
401 
402  fTransportEndSpin = aFieldTrack.GetSpin();
403  fParticleIsLooping = fFieldPropagator->IsParticleLooping() ;
404  endpointDistance = (fTransportEndPosition - startPosition).mag() ;
405  }
406 
407  // If we are asked to go a step length of 0, and we are on a boundary
408  // then a boundary will also limit the step -> we must flag this.
409  //
410  if( currentMinimumStep == 0.0 )
411  {
412  if( currentSafety == 0.0 ) fGeometryLimitedStep = true ;
413  }
414 
415  // Update the safety starting from the end-point,
416  // if it will become negative at the end-point.
417  //
418  if( currentSafety < endpointDistance )
419  {
420  // if( particleMagneticCharge == 0.0 )
421  // G4cout << " Avoiding call to ComputeSafety : charge = 0.0 " << G4endl;
422 
423  if( particleMagneticCharge != 0.0 ) {
424 
425  G4double endSafety =
426  fLinearNavigator->ComputeSafety( fTransportEndPosition) ;
427  currentSafety = endSafety ;
428  fPreviousSftOrigin = fTransportEndPosition ;
429  fPreviousSafety = currentSafety ;
430  fpSafetyHelper->SetCurrentSafety( currentSafety, fTransportEndPosition);
431 
432  // Because the Stepping Manager assumes it is from the start point,
433  // add the StepLength
434  //
435  currentSafety += endpointDistance ;
436 
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 ;
444 #endif
445  }
446  }
447 
448  fParticleChange.ProposeTrueStepLength(geometryStepLength) ;
449 
450  return geometryStepLength ;
451 }
452 
454 //
455 // Initialize ParticleChange (by setting all its members equal
456 // to corresponding members in G4Track)
457 
458 G4VParticleChange* G4MonopoleTransportation::AlongStepDoIt( const G4Track& track,
459  const G4Step& stepData )
460 {
461  static const G4ParticleDefinition* fOpticalPhoton =
462  G4ParticleTable::GetParticleTable()->FindParticle("opticalphoton");
463 
464 #ifdef G4VERBOSE
465  static G4int noCalls=0;
466  noCalls++;
467 #endif
468 
469  fParticleChange.Initialize(track) ;
470 
471  // Code for specific process
472  //
473  fParticleChange.ProposePosition(fTransportEndPosition) ;
474  fParticleChange.ProposeMomentumDirection(fTransportEndMomentumDir) ;
475  fParticleChange.ProposeEnergy(fTransportEndKineticEnergy) ;
476  fParticleChange.SetMomentumChanged(fMomentumChanged) ;
477 
478  fParticleChange.ProposePolarization(fTransportEndSpin);
479 
480  G4double deltaTime = 0.0 ;
481 
482  // Calculate Lab Time of Flight (ONLY if field Equations used it!)
483  // G4double endTime = fCandidateEndGlobalTime;
484  // G4double delta_time = endTime - startTime;
485 
486  G4double startTime = track.GetGlobalTime() ;
487 
488  if (!fEndGlobalTimeComputed)
489  {
490  // The time was not integrated .. make the best estimate possible
491  //
492  G4double finalVelocity = track.GetVelocity() ;
493  G4double initialVelocity = stepData.GetPreStepPoint()->GetVelocity() ;
494  G4double stepLength = track.GetStepLength() ;
495 
496  deltaTime= 0.0; // in case initialVelocity = 0
497  const G4DynamicParticle* fpDynamicParticle = track.GetDynamicParticle();
498  if (fpDynamicParticle->GetDefinition()== fOpticalPhoton)
499  {
500  // A photon is in the medium of the final point
501  // during the step, so it has the final velocity.
502  deltaTime = stepLength/finalVelocity ;
503  }
504  else if (finalVelocity > 0.0)
505  {
506  G4double meanInverseVelocity ;
507  // deltaTime = stepLength/finalVelocity ;
508  meanInverseVelocity = 0.5
509  * ( 1.0 / initialVelocity + 1.0 / finalVelocity ) ;
510  deltaTime = stepLength * meanInverseVelocity ;
511  }
512  else if( initialVelocity > 0.0 )
513  {
514  deltaTime = stepLength/initialVelocity ;
515  }
516  fCandidateEndGlobalTime = startTime + deltaTime ;
517  }
518  else
519  {
520  deltaTime = fCandidateEndGlobalTime - startTime ;
521  }
522 
523  fParticleChange.ProposeGlobalTime( fCandidateEndGlobalTime ) ;
524 
525  // Now Correct by Lorentz factor to get "proper" deltaTime
526 
527  G4double restMass = track.GetDynamicParticle()->GetMass() ;
528  G4double deltaProperTime = deltaTime*( restMass/track.GetTotalEnergy() ) ;
529 
530  fParticleChange.ProposeProperTime(track.GetProperTime() + deltaProperTime) ;
531  //fParticleChange. ProposeTrueStepLength( track.GetStepLength() ) ;
532 
533  // If the particle is caught looping or is stuck (in very difficult
534  // boundaries) in a magnetic field (doing many steps)
535  // THEN this kills it ...
536  //
537  if ( fParticleIsLooping )
538  {
539  G4double endEnergy= fTransportEndKineticEnergy;
540 
541  if( (endEnergy < fThreshold_Important_Energy)
542  || (fNoLooperTrials >= fThresholdTrials ) ){
543  // Kill the looping particle
544  //
545  fParticleChange.ProposeTrackStatus( fStopAndKill ) ;
546 
547  // 'Bare' statistics
548  fSumEnergyKilled += endEnergy;
549  if( endEnergy > fMaxEnergyKilled) { fMaxEnergyKilled= endEnergy; }
550 
551 #ifdef G4VERBOSE
552  if( (fVerboseLevel > 1) ||
553  ( endEnergy > fThreshold_Warning_Energy ) ) {
554  G4cout << " G4MonopoleTransportation is killing track that is looping or stuck "
555  << G4endl
556  << " This track has " << track.GetKineticEnergy() / MeV
557  << " MeV energy." << G4endl;
558  G4cout << " Number of trials = " << fNoLooperTrials
559  << " No of calls to AlongStepDoIt = " << noCalls
560  << G4endl;
561  }
562 #endif
563  fNoLooperTrials=0;
564  }
565  else{
566  fNoLooperTrials ++;
567 #ifdef G4VERBOSE
568  if( (fVerboseLevel > 2) ){
569  G4cout << " G4MonopoleTransportation::AlongStepDoIt(): Particle looping - "
570  << " Number of trials = " << fNoLooperTrials
571  << " No of calls to = " << noCalls
572  << G4endl;
573  }
574 #endif
575  }
576  }else{
577  fNoLooperTrials=0;
578  }
579 
580  // Another (sometimes better way) is to use a user-limit maximum Step size
581  // to alleviate this problem ..
582 
583  // Introduce smooth curved trajectories to particle-change
584  //
585  fParticleChange.SetPointerToVectorOfAuxiliaryPoints
586  (fFieldPropagator->GimmeTrajectoryVectorAndForgetIt() );
587 
588  return &fParticleChange ;
589 }
590 
592 //
593 // This ensures that the PostStep action is always called,
594 // so that it can do the relocation if it is needed.
595 //
596 
597 G4double G4MonopoleTransportation::
598 PostStepGetPhysicalInteractionLength( const G4Track& track,
599  G4double, // previousStepSize
600  G4ForceCondition* pForceCond )
601 {
602  //magSetup->SetStepperAndChordFinder(0);
603  // change back to usual equation
604  G4FieldManager* fieldMgr=fFieldPropagator->FindAndSetFieldManager( track.GetVolume() );
605  fFieldBuilder->setStepperAndChordFinder (fieldMgr, 0);
606 
607  *pForceCond = Forced ;
608  return DBL_MAX ; // was kInfinity ; but convention now is DBL_MAX
609 }
610 
612 //
613 
614 G4VParticleChange* G4MonopoleTransportation::PostStepDoIt( const G4Track& track,
615  const G4Step& )
616 {
617  G4TouchableHandle retCurrentTouchable ; // The one to return
618 
619  // Initialize ParticleChange (by setting all its members equal
620  // to corresponding members in G4Track)
621  // fParticleChange.Initialize(track) ; // To initialise TouchableChange
622 
623  fParticleChange.ProposeTrackStatus(track.GetTrackStatus()) ;
624 
625  // If the Step was determined by the volume boundary,
626  // logically relocate the particle
627 
628  if(fGeometryLimitedStep)
629  {
630  // fCurrentTouchable will now become the previous touchable,
631  // and what was the previous will be freed.
632  // (Needed because the preStepPoint can point to the previous touchable)
633 
634  fLinearNavigator->SetGeometricallyLimitedStep() ;
635  fLinearNavigator->
636  LocateGlobalPointAndUpdateTouchableHandle( track.GetPosition(),
637  track.GetMomentumDirection(),
638  fCurrentTouchableHandle,
639  true ) ;
640  // Check whether the particle is out of the world volume
641  // If so it has exited and must be killed.
642  //
643  if( fCurrentTouchableHandle->GetVolume() == 0 )
644  {
645  fParticleChange.ProposeTrackStatus( fStopAndKill ) ;
646  }
647  retCurrentTouchable = fCurrentTouchableHandle ;
648  fParticleChange.SetTouchableHandle( fCurrentTouchableHandle ) ;
649  }
650  else // fGeometryLimitedStep is false
651  {
652  // This serves only to move the Navigator's location
653  //
654  fLinearNavigator->LocateGlobalPointWithinVolume( track.GetPosition() ) ;
655 
656  // The value of the track's current Touchable is retained.
657  // (and it must be correct because we must use it below to
658  // overwrite the (unset) one in particle change)
659  // It must be fCurrentTouchable too ??
660  //
661  fParticleChange.SetTouchableHandle( track.GetTouchableHandle() ) ;
662  retCurrentTouchable = track.GetTouchableHandle() ;
663  } // endif ( fGeometryLimitedStep )
664 
665  const G4VPhysicalVolume* pNewVol = retCurrentTouchable->GetVolume() ;
666  const G4Material* pNewMaterial = 0 ;
667  const G4VSensitiveDetector* pNewSensitiveDetector = 0 ;
668 
669  if( pNewVol != 0 )
670  {
671  pNewMaterial= pNewVol->GetLogicalVolume()->GetMaterial();
672  pNewSensitiveDetector= pNewVol->GetLogicalVolume()->GetSensitiveDetector();
673  }
674 
675  // ( <const_cast> pNewMaterial ) ;
676  // ( <const_cast> pNewSensitiveDetector) ;
677 
678  fParticleChange.SetMaterialInTouchable( (G4Material *) pNewMaterial ) ;
679  fParticleChange.SetSensitiveDetectorInTouchable( (G4VSensitiveDetector *) pNewSensitiveDetector ) ;
680 
681  const G4MaterialCutsCouple* pNewMaterialCutsCouple = 0;
682  if( pNewVol != 0 )
683  {
684  pNewMaterialCutsCouple=pNewVol->GetLogicalVolume()->GetMaterialCutsCouple();
685  }
686 
687  if( pNewVol!=0 && pNewMaterialCutsCouple!=0 && pNewMaterialCutsCouple->GetMaterial()!=pNewMaterial )
688  {
689  // for parametrized volume
690  //
691  pNewMaterialCutsCouple =
692  G4ProductionCutsTable::GetProductionCutsTable()
693  ->GetMaterialCutsCouple(pNewMaterial,
694  pNewMaterialCutsCouple->GetProductionCuts());
695  }
696  fParticleChange.SetMaterialCutsCoupleInTouchable( pNewMaterialCutsCouple );
697 
698  // temporarily until Get/Set Material of ParticleChange,
699  // and StepPoint can be made const.
700  // Set the touchable in ParticleChange
701  // this must always be done because the particle change always
702  // uses this value to overwrite the current touchable pointer.
703  //
704  fParticleChange.SetTouchableHandle(retCurrentTouchable) ;
705 
706  return &fParticleChange ;
707 }
708 
709 // New method takes over the responsibility to reset the state of G4MonopoleTransportation
710 // object at the start of a new track or the resumption of a suspended track.
711 
712 void
713 G4MonopoleTransportation::StartTracking(G4Track* aTrack)
714 {
715  G4VProcess::StartTracking(aTrack);
716 
717 // The actions here are those that were taken in AlongStepGPIL
718 // when track.GetCurrentStepNumber()==1
719 
720  // reset safety value and center
721  //
722  fPreviousSafety = 0.0 ;
723  fPreviousSftOrigin = G4ThreeVector(0.,0.,0.) ;
724 
725  // reset looping counter -- for motion in field
726  fNoLooperTrials= 0;
727  // Must clear this state .. else it depends on last track's value
728  // --> a better solution would set this from state of suspended track TODO ?
729  // Was if( aTrack->GetCurrentStepNumber()==1 ) { .. }
730 
731  // ChordFinder reset internal state
732  //
733  if( DoesGlobalFieldExist() ) {
734  fFieldPropagator->ClearPropagatorState();
735  // Resets all state of field propagator class (ONLY)
736  // including safety values (in case of overlaps and to wipe for first track).
737 
738  // G4ChordFinder* chordF= fFieldPropagator->GetChordFinder();
739  // if( chordF ) chordF->ResetStepEstimate();
740  }
741 
742  // Make sure to clear the chord finders of all fields (ie managers)
743  static G4FieldManagerStore* fieldMgrStore= G4FieldManagerStore::GetInstance();
744  fieldMgrStore->ClearAllChordFindersState();
745 
746  // Update the current touchable handle (from the track's)
747  //
748  fCurrentTouchableHandle = aTrack->GetTouchableHandle();
749 }
750 
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
selection
main part
Definition: corrVsCorr.py:98
T sqrt(T t)
Definition: SSEVec.h:48
Square< F >::type sqr(const F &f)
Definition: Square.h:13
volatile std::atomic< bool > shutdown_flag false
float spin(float ph)