CMS 3D CMS Logo

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