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