Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4Transportation.cc
이 파일의 문서화 페이지로 가기
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 //
27 // $Id: G4Transportation.cc 2011/06/10 16:19:46 japost Exp japost $
28 //
29 // ------------------------------------------------------------
30 // GEANT 4 include file implementation
31 //
32 // ------------------------------------------------------------
33 //
34 // This class is a process responsible for the transportation of
35 // a particle, ie the geometrical propagation that encounters the
36 // geometrical sub-volumes of the detectors.
37 //
38 // It is also tasked with the key role of proposing the "isotropic safety",
39 // which will be used to update the post-step point's safety.
40 //
41 // =======================================================================
42 // Modified:
43 // 10 Jan 2015, M.Kelsey: Use G4DynamicParticle mass, NOT PDGMass
44 // 28 Oct 2011, P.Gumpl./J.Ap: Detect gravity field, use magnetic moment
45 // 20 Nov 2008, J.Apostolakis: Push safety to helper - after ComputeSafety
46 // 9 Nov 2007, J.Apostolakis: Flag for short steps, push safety to helper
47 // 19 Jan 2006, P.MoraDeFreitas: Fix for suspended tracks (StartTracking)
48 // 11 Aug 2004, M.Asai: Add G4VSensitiveDetector* for updating stepPoint.
49 // 21 June 2003, J.Apostolakis: Calling field manager with
50 // track, to enable it to configure its accuracy
51 // 13 May 2003, J.Apostolakis: Zero field areas now taken into
52 // account correclty in all cases (thanks to W Pokorski).
53 // 29 June 2001, J.Apostolakis, D.Cote-Ahern, P.Gumplinger:
54 // correction for spin tracking
55 // 20 Febr 2001, J.Apostolakis: update for new FieldTrack
56 // 22 Sept 2000, V.Grichine: update of Kinetic Energy
57 // Created: 19 March 1997, J. Apostolakis
58 // =======================================================================
59 
60 #include "G4Transportation.hh"
63 
64 #include "G4PhysicalConstants.hh"
65 #include "G4SystemOfUnits.hh"
66 #include "G4ProductionCutsTable.hh"
67 #include "G4ParticleTable.hh"
68 
69 #include "G4ChargeState.hh"
70 #include "G4EquationOfMotion.hh"
71 
72 #include "G4FieldManagerStore.hh"
74 
76 
78 
80 //
81 // Constructor
82 
84  : G4VProcess( G4String("Transportation"), fTransportation ),
85  fTransportEndPosition( 0.0, 0.0, 0.0 ),
86  fTransportEndMomentumDir( 0.0, 0.0, 0.0 ),
87  fTransportEndKineticEnergy( 0.0 ),
88  fTransportEndSpin( 0.0, 0.0, 0.0 ),
89  fMomentumChanged(true),
90  fEndGlobalTimeComputed(false),
91  fCandidateEndGlobalTime(0.0),
92  fParticleIsLooping( false ),
93  fNewTrack( true ),
94  fFirstStepInVolume( true ),
95  fLastStepInVolume( false ),
96  fGeometryLimitedStep(true),
97  fFieldExertedForce( false ),
98  fPreviousSftOrigin( 0.,0.,0. ),
99  fPreviousSafety( 0.0 ),
100  fEndPointDistance( -1.0 ),
101  fThreshold_Warning_Energy( 1.0 * CLHEP::kiloelectronvolt ),
102  fThreshold_Important_Energy( 1.0 * MeV ),
103  fThresholdTrials( 10 ),
104  fNoLooperTrials( 0 ),
105  fSumEnergyKilled( 0.0 ), fMaxEnergyKilled( 0.0 ),
106  fShortStepOptimisation( false ) // Old default: true (=fast short steps)
107 {
108  SetProcessSubType(static_cast<G4int>(TRANSPORTATION));
109  pParticleChange= &fParticleChange; // Required to conform to G4VProcess
110  SetVerboseLevel(verbosity);
111 
112  G4TransportationManager* transportMgr ;
113 
115 
116  fLinearNavigator = transportMgr->GetNavigatorForTracking() ;
117 
118  fFieldPropagator = transportMgr->GetPropagatorInField() ;
119 
120  fpSafetyHelper = transportMgr->GetSafetyHelper(); // New
121 
122  fpLogger = new G4TransportationLogger("G4Transportation", verbosity);
124  GetThresholdTrials() );
125 
126  // Cannot determine whether a field exists here, as it would
127  // depend on the relative order of creating the detector's
128  // field and this process. That order is not guaranted.
129  // Instead later the method DoesGlobalFieldExist() is called
130 
131  static G4ThreadLocal G4TouchableHandle* pNullTouchableHandle = 0;
132  if ( !pNullTouchableHandle) { pNullTouchableHandle = new G4TouchableHandle; }
133  fCurrentTouchableHandle = *pNullTouchableHandle;
134  // Points to (G4VTouchable*) 0
135 
136 #ifdef G4VERBOSE
137  if( verboseLevel > 0)
138  {
139  G4cout << " G4Transportation constructor> set fShortStepOptimisation to ";
140  if ( fShortStepOptimisation ) G4cout << "true" << G4endl;
141  else G4cout << "false" << G4endl;
142  }
143 #endif
144 }
145 
147 
149 {
150  delete fpLogger;
151  if( (verboseLevel > 0) && (fSumEnergyKilled > 0.0 ) )
152  {
153  G4cout << " G4Transportation: Statistics for looping particles " << G4endl;
154  G4cout << " Sum of energy of loopers killed: "
155  << fSumEnergyKilled << G4endl;
156  G4cout << " Max energy of loopers killed: "
157  << fMaxEnergyKilled << G4endl;
158  }
159 }
160 
162 //
163 // Responsibilities:
164 // Find whether the geometry limits the Step, and to what length
165 // Calculate the new value of the safety and return it.
166 // Store the final time, position and momentum.
167 
170  G4double, // previousStepSize
171  G4double currentMinimumStep,
172  G4double& currentSafety,
173  G4GPILSelection* selection )
174 {
175  G4double geometryStepLength= -1.0, newSafety= -1.0;
176  fParticleIsLooping = false ;
177 
178  // Initial actions moved to StartTrack()
179  // --------------------------------------
180  // Note: in case another process changes touchable handle
181  // it will be necessary to add here (for all steps)
182  // fCurrentTouchableHandle = aTrack->GetTouchableHandle();
183 
184  // GPILSelection is set to defaule value of CandidateForSelection
185  // It is a return value
186  //
187  *selection = CandidateForSelection ;
188 
190  fLastStepInVolume= false;
191  fNewTrack = false;
192 
194 
195  // Get initial Energy/Momentum of the track
196  //
197  const G4DynamicParticle* pParticle = track.GetDynamicParticle() ;
198  const G4ParticleDefinition* pParticleDef = pParticle->GetDefinition() ;
199  G4ThreeVector startMomentumDir = pParticle->GetMomentumDirection() ;
200  G4ThreeVector startPosition = track.GetPosition() ;
201 
202  // G4double theTime = track.GetGlobalTime() ;
203 
204  // The Step Point safety can be limited by other geometries and/or the
205  // assumptions of any process - it's not always the geometrical safety.
206  // We calculate the starting point's isotropic safety here.
207  //
208  G4ThreeVector OriginShift = startPosition - fPreviousSftOrigin ;
209  G4double MagSqShift = OriginShift.mag2() ;
210  if( MagSqShift >= sqr(fPreviousSafety) )
211  {
212  currentSafety = 0.0 ;
213  }
214  else
215  {
216  currentSafety = fPreviousSafety - std::sqrt(MagSqShift) ;
217  }
218 
219  // Is the particle charged or has it a magnetic moment?
220  //
221  G4double particleCharge = pParticle->GetCharge() ;
222  G4double magneticMoment = pParticle->GetMagneticMoment() ;
223  G4double restMass = pParticle->GetMass() ;
224 
225  fGeometryLimitedStep = false ;
226  // fEndGlobalTimeComputed = false ;
227 
228  // There is no need to locate the current volume. It is Done elsewhere:
229  // On track construction
230  // By the tracking, after all AlongStepDoIts, in "Relocation"
231 
232  // Check if the particle has a force, EM or gravitational, exerted on it
233  //
234  G4FieldManager* fieldMgr=0;
235  G4bool fieldExertsForce = false ;
236 
237  G4bool gravityOn = false;
238  G4bool fieldExists= false; // Field is not 0 (null pointer)
239 
240  fieldMgr = fFieldPropagator->FindAndSetFieldManager( track.GetVolume() );
241  if( fieldMgr != 0 )
242  {
243  // Message the field Manager, to configure it for this track
244  fieldMgr->ConfigureForTrack( &track );
245  // Is here to allow a transition from no-field pointer
246  // to finite field (non-zero pointer).
247 
248  // If the field manager has no field ptr, the field is zero
249  // by definition ( = there is no field ! )
250  const G4Field* ptrField= fieldMgr->GetDetectorField();
251  fieldExists = (ptrField!=0) ;
252  if( fieldExists )
253  {
254  gravityOn= ptrField->IsGravityActive();
255 
256  if( (particleCharge != 0.0)
257  || (fUseMagneticMoment && (magneticMoment != 0.0) )
258  || (gravityOn && (restMass != 0.0) )
259  )
260  {
261  fieldExertsForce = fieldExists;
262  }
263  }
264  }
265  fFieldExertedForce = fieldExertsForce;
266 
267  if( !fieldExertsForce )
268  {
269  G4double linearStepLength ;
270  if( fShortStepOptimisation && (currentMinimumStep <= currentSafety) )
271  {
272  // The Step is guaranteed to be taken
273  //
274  geometryStepLength = currentMinimumStep ;
275  fGeometryLimitedStep = false ;
276  }
277  else
278  {
279  // Find whether the straight path intersects a volume
280  //
281  linearStepLength = fLinearNavigator->ComputeStep( startPosition,
282  startMomentumDir,
283  currentMinimumStep,
284  newSafety) ;
285  // Remember last safety origin & value.
286  //
287  fPreviousSftOrigin = startPosition ;
288  fPreviousSafety = newSafety ;
289  fpSafetyHelper->SetCurrentSafety( newSafety, startPosition);
290 
291  currentSafety = newSafety ;
292 
293  fGeometryLimitedStep= (linearStepLength <= currentMinimumStep);
295  {
296  // The geometry limits the Step size (an intersection was found.)
297  geometryStepLength = linearStepLength ;
298  }
299  else
300  {
301  // The full Step is taken.
302  geometryStepLength = currentMinimumStep ;
303  }
304  }
305  fEndPointDistance = geometryStepLength ;
306 
307  // Calculate final position
308  //
309  fTransportEndPosition = startPosition+geometryStepLength*startMomentumDir ;
310 
311  // Momentum direction, energy and polarisation are unchanged by transport
312  //
313  fTransportEndMomentumDir = startMomentumDir ;
316  fParticleIsLooping = false ;
317  fMomentumChanged = false ;
318  fEndGlobalTimeComputed = false ;
319  }
320  else // A field exerts force
321  {
322  G4double momentumMagnitude = pParticle->GetTotalMomentum() ;
323  G4ThreeVector EndUnitMomentum ;
324  G4double lengthAlongCurve ;
325 
326  G4ChargeState chargeState(particleCharge, // The charge can change (dynamic)
327  magneticMoment,
328  pParticleDef->GetPDGSpin() );
329  // For insurance, could set it again
330  // chargeState.SetPDGSpin(pParticleDef->GetPDGSpin() ); // Provisionally in same object
331 
332  auto equationOfMotion = fFieldPropagator->GetCurrentEquationOfMotion();
333 
334  equationOfMotion->SetChargeMomentumMass( chargeState,
335  momentumMagnitude,
336  restMass);
337 
338  G4FieldTrack aFieldTrack = G4FieldTrack( startPosition,
339  track.GetGlobalTime(), // Lab.
340  // track.GetProperTime(), // Particle rest frame
341  track.GetMomentumDirection(),
342  track.GetKineticEnergy(),
343  restMass,
344  particleCharge,
345  track.GetPolarization(),
346  pParticleDef->GetPDGMagneticMoment(),
347  0.0, // Length along track
348  pParticleDef->GetPDGSpin()
349  ) ;
350 
351  if( currentMinimumStep > 0 )
352  {
353  // Do the Transport in the field (non recti-linear)
354  //
355  lengthAlongCurve = fFieldPropagator->ComputeStep( aFieldTrack,
356  currentMinimumStep,
357  currentSafety,
358  track.GetVolume() ) ;
359 
361  // It is possible that step was reduced in PropagatorInField due to previous zero steps
362  // To cope with case that reduced step is taken in full, we must rely on PiF to obtain this
363  // value.
364 
365  geometryStepLength = std::min( lengthAlongCurve, currentMinimumStep );
366 
367  // Remember last safety origin & value.
368  //
369  fPreviousSftOrigin = startPosition ;
370  fPreviousSafety = currentSafety ;
371  fpSafetyHelper->SetCurrentSafety( currentSafety, startPosition);
372  }
373  else
374  {
375  geometryStepLength = lengthAlongCurve= 0.0 ;
376  fGeometryLimitedStep = false ;
377  }
378 
379  // Get the End-Position and End-Momentum (Dir-ection)
380  //
381  fTransportEndPosition = aFieldTrack.GetPosition() ;
382 
383  // Momentum: Magnitude and direction can be changed too now ...
384  //
385  fMomentumChanged = true ;
386  fTransportEndMomentumDir = aFieldTrack.GetMomentumDir() ;
387 
389 
391  {
392  // If the field can change energy, then the time must be integrated
393  // - so this should have been updated
394  //
396  fEndGlobalTimeComputed = true;
397 
398  // was ( fCandidateEndGlobalTime != track.GetGlobalTime() );
399  // a cleaner way is to have FieldTrack knowing whether time is updated.
400  }
401  else
402  {
403  // The energy should be unchanged by field transport,
404  // - so the time changed will be calculated elsewhere
405  //
406  fEndGlobalTimeComputed = false;
407 
408  // Check that the integration preserved the energy
409  // - and if not correct this!
410  G4double startEnergy= track.GetKineticEnergy();
412 
413  static G4ThreadLocal G4int no_inexact_steps=0, no_large_ediff;
414  G4double absEdiff = std::fabs(startEnergy- endEnergy);
415  if( absEdiff > perMillion * endEnergy )
416  {
417  no_inexact_steps++;
418  // Possible statistics keeping here ...
419  }
420  if( verboseLevel > 1 )
421  {
422  if( std::fabs(startEnergy- endEnergy) > perThousand * endEnergy )
423  {
424  static G4ThreadLocal G4int no_warnings= 0, warnModulo=1,
425  moduloFactor= 10;
426  no_large_ediff ++;
427  if( (no_large_ediff% warnModulo) == 0 )
428  {
429  no_warnings++;
430  G4cout << "WARNING - G4Transportation::AlongStepGetPIL() "
431  << " Energy change in Step is above 1^-3 relative value. " << G4endl
432  << " Relative change in 'tracking' step = "
433  << std::setw(15) << (endEnergy-startEnergy)/startEnergy << G4endl
434  << " Starting E= " << std::setw(12) << startEnergy / MeV << " MeV " << G4endl
435  << " Ending E= " << std::setw(12) << endEnergy / MeV << " MeV " << G4endl;
436  G4cout << " Energy has been corrected -- however, review"
437  << " field propagation parameters for accuracy." << G4endl;
438  if( (verboseLevel > 2 ) || (no_warnings<4) || (no_large_ediff == warnModulo * moduloFactor) )
439  {
440  G4cout << " These include EpsilonStepMax(/Min) in G4FieldManager "
441  << " which determine fractional error per step for integrated quantities. " << G4endl
442  << " Note also the influence of the permitted number of integration steps."
443  << G4endl;
444  }
445  G4cerr << "ERROR - G4Transportation::AlongStepGetPIL()" << G4endl
446  << " Bad 'endpoint'. Energy change detected"
447  << " and corrected. "
448  << " Has occurred already "
449  << no_large_ediff << " times." << G4endl;
450  if( no_large_ediff == warnModulo * moduloFactor )
451  {
452  warnModulo *= moduloFactor;
453  }
454  }
455  }
456  } // end of if (verboseLevel)
457 
458  // Correct the energy for fields that conserve it
459  // This - hides the integration error
460  // - but gives a better physical answer
462  }
463 
464  fTransportEndSpin = aFieldTrack.GetSpin();
466  fEndPointDistance = (fTransportEndPosition - startPosition).mag() ;
467  }
468 
469  // If we are asked to go a step length of 0, and we are on a boundary
470  // then a boundary will also limit the step -> we must flag this.
471  //
472  if( currentMinimumStep == 0.0 )
473  {
474  if( currentSafety == 0.0 ) { fGeometryLimitedStep = true; }
475  }
476 
477  // Update the safety starting from the end-point,
478  // if it will become negative at the end-point.
479  //
480  if( currentSafety < fEndPointDistance )
481  {
482  if( particleCharge != 0.0 )
483  {
484  G4double endSafety =
486  currentSafety = endSafety ;
487  fPreviousSftOrigin = fTransportEndPosition ;
488  fPreviousSafety = currentSafety ;
490 
491  // Because the Stepping Manager assumes it is from the start point,
492  // add the StepLength
493  //
494  currentSafety += fEndPointDistance ;
495 
496 #ifdef G4DEBUG_TRANSPORT
497  G4cout.precision(12) ;
498  G4cout << "***G4Transportation::AlongStepGPIL ** " << G4endl ;
499  G4cout << " Called Navigator->ComputeSafety at " << fTransportEndPosition
500  << " and it returned safety= " << endSafety << G4endl ;
501  G4cout << " Adding endpoint distance " << fEndPointDistance
502  << " to obtain pseudo-safety= " << currentSafety << G4endl ;
503  }
504  else
505  {
506  G4cout << "***G4Transportation::AlongStepGPIL ** " << G4endl ;
507  G4cout << " Avoiding call to ComputeSafety : " << G4endl;
508  G4cout << " charge = " << particleCharge << G4endl;
509  G4cout << " mag moment = " << magneticMoment << G4endl;
510 #endif
511  }
512  }
513 
514  fParticleChange.ProposeTrueStepLength(geometryStepLength) ;
515 
516  return geometryStepLength ;
517 }
518 
520 //
521 // Initialize ParticleChange (by setting all its members equal
522 // to corresponding members in G4Track)
523 
525  const G4Step& stepData )
526 {
527  static G4ThreadLocal G4long noCallsASDI=0;
528  const char *methodName= "AlongStepDoIt";
529  noCallsASDI++;
530 
531  fParticleChange.Initialize(track) ;
532 
533  // Code for specific process
534  //
539 
541 
542  G4double deltaTime = 0.0 ;
543 
544  // Calculate Lab Time of Flight (ONLY if field Equations used it!)
545  // G4double endTime = fCandidateEndGlobalTime;
546  // G4double delta_time = endTime - startTime;
547 
548  G4double startTime = track.GetGlobalTime() ;
549 
551  {
552  // The time was not integrated .. make the best estimate possible
553  //
554  G4double initialVelocity = stepData.GetPreStepPoint()->GetVelocity();
556 
557  deltaTime= 0.0; // in case initialVelocity = 0
558  if ( initialVelocity > 0.0 ) { deltaTime = stepLength/initialVelocity; }
559 
560  fCandidateEndGlobalTime = startTime + deltaTime ;
561  fParticleChange.ProposeLocalTime( track.GetLocalTime() + deltaTime) ;
562  }
563  else
564  {
565  deltaTime = fCandidateEndGlobalTime - startTime ;
567  }
568 
569 
570  // Now Correct by Lorentz factor to get delta "proper" Time
571 
572  G4double restMass = track.GetDynamicParticle()->GetMass() ;
573  G4double deltaProperTime = deltaTime*( restMass/track.GetTotalEnergy() ) ;
574 
575  fParticleChange.ProposeProperTime(track.GetProperTime() + deltaProperTime) ;
576  //fParticleChange. ProposeTrueStepLength( track.GetStepLength() ) ;
577 
578  // If the particle is caught looping or is stuck (in very difficult
579  // boundaries) in a magnetic field (doing many steps) THEN this kills it ...
580  //
581  if ( fParticleIsLooping )
582  {
584  fNoLooperTrials ++;
585 
586  if( (endEnergy < fThreshold_Important_Energy)
588  {
589  // Kill the looping particle
590  //
592 
593  // 'Bare' statistics
594  //
595  fSumEnergyKilled += endEnergy;
596  if( endEnergy > fMaxEnergyKilled) { fMaxEnergyKilled= endEnergy; }
597 
598  if( endEnergy > fThreshold_Warning_Energy )
599  {
600  fpLogger->ReportLoopingTrack( track, stepData, fNoLooperTrials, noCallsASDI, methodName );
601  // const char* fullMethodName= "G4Transportation::AlongStepDotIt()";
602  // ReportLoopingTrack( track, stepData, noCallsASDI, fullMethodName );
603  }
604  fNoLooperTrials=0;
605  }
606  else
607  {
608 #ifdef G4VERBOSE
609  if( verboseLevel > 2 )
610  {
611  G4cout << " G4Transportation::AlongStepDoIt(): Particle looping - "
612  << G4endl
613  << " Number of trials = " << fNoLooperTrials
614  << G4endl
615  << " No of calls to = " << noCallsASDI << G4endl;
616  }
617 #endif
618  }
619  }
620  else
621  {
622  fNoLooperTrials=0;
623  }
624 
625  // Another (sometimes better way) is to use a user-limit maximum Step size
626  // to alleviate this problem ..
627 
628  // Introduce smooth curved trajectories to particle-change
629  //
632 
633  return &fParticleChange ;
634 }
635 
637 //
638 // This ensures that the PostStep action is always called,
639 // so that it can do the relocation if it is needed.
640 //
641 
644  G4double, // previousStepSize
645  G4ForceCondition* pForceCond )
646 {
647  fFieldExertedForce = false; // Not known
648  *pForceCond = Forced ;
649  return DBL_MAX ; // was kInfinity ; but convention now is DBL_MAX
650 }
651 
653 //
654 
656  const G4Step& )
657 {
658  G4TouchableHandle retCurrentTouchable ; // The one to return
659  G4bool isLastStep= false;
660 
661  // Initialize ParticleChange (by setting all its members equal
662  // to corresponding members in G4Track)
663  // fParticleChange.Initialize(track) ; // To initialise TouchableChange
664 
666 
667  // If the Step was determined by the volume boundary,
668  // logically relocate the particle
669 
671  {
672  // fCurrentTouchable will now become the previous touchable,
673  // and what was the previous will be freed.
674  // (Needed because the preStepPoint can point to the previous touchable)
675 
678  LocateGlobalPointAndUpdateTouchableHandle( track.GetPosition(),
679  track.GetMomentumDirection(),
681  true ) ;
682  // Check whether the particle is out of the world volume
683  // If so it has exited and must be killed.
684  //
685  if( fCurrentTouchableHandle->GetVolume() == 0 )
686  {
688  }
689  retCurrentTouchable = fCurrentTouchableHandle ;
691 
692  // Update the Step flag which identifies the Last Step in a volume
693  if( !fFieldExertedForce )
694  isLastStep = fLinearNavigator->ExitedMotherVolume()
696  else
697  isLastStep = fFieldPropagator->IsLastStepInVolume();
698  }
699  else // fGeometryLimitedStep is false
700  {
701  // This serves only to move the Navigator's location
702  //
704 
705  // The value of the track's current Touchable is retained.
706  // (and it must be correct because we must use it below to
707  // overwrite the (unset) one in particle change)
708  // It must be fCurrentTouchable too ??
709  //
711  retCurrentTouchable = track.GetTouchableHandle() ;
712 
713  isLastStep= false;
714  } // endif ( fGeometryLimitedStep )
715  fLastStepInVolume= isLastStep;
716 
719 
720  const G4VPhysicalVolume* pNewVol = retCurrentTouchable->GetVolume() ;
721  const G4Material* pNewMaterial = 0 ;
722  const G4VSensitiveDetector* pNewSensitiveDetector = 0 ;
723 
724  if( pNewVol != 0 )
725  {
726  pNewMaterial= pNewVol->GetLogicalVolume()->GetMaterial();
727  pNewSensitiveDetector= pNewVol->GetLogicalVolume()->GetSensitiveDetector();
728  }
729 
730  // ( <const_cast> pNewMaterial ) ;
731  // ( <const_cast> pNewSensitiveDetector) ;
732 
735 
736  const G4MaterialCutsCouple* pNewMaterialCutsCouple = 0;
737  if( pNewVol != 0 )
738  {
739  pNewMaterialCutsCouple=pNewVol->GetLogicalVolume()->GetMaterialCutsCouple();
740  }
741 
742  if( pNewVol!=0 && pNewMaterialCutsCouple!=0 && pNewMaterialCutsCouple->GetMaterial()!=pNewMaterial )
743  {
744  // for parametrized volume
745  //
746  pNewMaterialCutsCouple =
748  ->GetMaterialCutsCouple(pNewMaterial,
749  pNewMaterialCutsCouple->GetProductionCuts());
750  }
751  fParticleChange.SetMaterialCutsCoupleInTouchable( pNewMaterialCutsCouple );
752 
753  // temporarily until Get/Set Material of ParticleChange,
754  // and StepPoint can be made const.
755  // Set the touchable in ParticleChange
756  // this must always be done because the particle change always
757  // uses this value to overwrite the current touchable pointer.
758  //
759  fParticleChange.SetTouchableHandle(retCurrentTouchable) ;
760 
761  return &fParticleChange ;
762 }
763 
764 // New method takes over the responsibility to reset the state of G4Transportation
765 // object at the start of a new track or the resumption of a suspended track.
766 
767 void
769 {
771  fNewTrack= true;
772  fFirstStepInVolume= true;
773  fLastStepInVolume= false;
774 
775  // The actions here are those that were taken in AlongStepGPIL
776  // when track.GetCurrentStepNumber()==1
777 
778  // reset safety value and center
779  //
780  fPreviousSafety = 0.0 ;
781  fPreviousSftOrigin = G4ThreeVector(0.,0.,0.) ;
782 
783  // reset looping counter -- for motion in field
784  fNoLooperTrials= 0;
785  // Must clear this state .. else it depends on last track's value
786  // --> a better solution would set this from state of suspended track TODO ?
787  // Was if( aTrack->GetCurrentStepNumber()==1 ) { .. }
788 
789  // ChordFinder reset internal state
790  //
791  if( DoesGlobalFieldExist() )
792  {
794  // Resets all state of field propagator class (ONLY)
795  // including safety values (in case of overlaps and to wipe for first track).
796 
797  // G4ChordFinder* chordF= fFieldPropagator->GetChordFinder();
798  // if( chordF ) chordF->ResetStepEstimate();
799  }
800 
801  // Make sure to clear the chord finders of all fields (ie managers)
802  //
804  fieldMgrStore->ClearAllChordFindersState();
805 
806  // Update the current touchable handle (from the track's)
807  //
809 
810  // Inform field propagator of new track
812 }
813 
814 // ------------------------===========================---------------------------
815 
817 {
818  G4bool lastValue= fUseMagneticMoment;
819  fUseMagneticMoment= useMoment;
821  return lastValue;
822 }
G4double GetLocalTime() const
G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData)
G4double GetKineticEnergy() const
virtual void ConfigureForTrack(const G4Track *)
virtual G4double ComputeStep(const G4ThreeVector &pGlobalPoint, const G4ThreeVector &pDirection, const G4double pCurrentProposedStepLength, G4double &pNewSafety)
Definition: G4Navigator.cc:753
CLHEP::Hep3Vector G4ThreeVector
G4SafetyHelper * GetSafetyHelper() const
void ProposePolarization(G4double Px, G4double Py, G4double Pz)
void ProposeGlobalTime(G4double t)
G4LogicalVolume * GetLogicalVolume() const
G4bool DoesFieldChangeEnergy() const
static constexpr double MeV
Definition: G4SIunits.hh:214
G4VSensitiveDetector * GetSensitiveDetector() const
void SetCurrentSafety(G4double val, const G4ThreeVector &pos)
static G4FieldManagerStore * GetInstance()
void SetThresholds(G4double newEnWarn, G4double importantEnergy, G4int newMaxTrials)
static G4bool EnableUseMagneticMoment(G4bool useMoment=true)
G4StepPoint * GetPreStepPoint() const
G4bool EnteredDaughterVolume() const
virtual void LocateGlobalPointWithinVolume(const G4ThreeVector &position)
Definition: G4Navigator.cc:588
G4double GetStepLength() const
#define G4endl
Definition: G4ios.hh:61
const G4ThreeVector & GetMomentumDirection() const
G4ReferenceCountedHandle< G4VTouchable > G4TouchableHandle
G4EquationOfMotion * GetCurrentEquationOfMotion()
#define fFieldExertedForce
void SetSensitiveDetectorInTouchable(G4VSensitiveDetector *fSensitiveDetector)
G4FieldManager * GetCurrentFieldManager()
#define fPreviousSafety
G4ThreeVector fPreviousSftOrigin
const G4TouchableHandle & GetTouchableHandle() const
const G4ThreeVector & GetPolarization() const
static constexpr double perMillion
Definition: G4SIunits.hh:334
G4ParticleChangeForTransport fParticleChange
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:44
G4Navigator * GetNavigatorForTracking() const
void SetPointerToVectorOfAuxiliaryPoints(std::vector< G4ThreeVector > *theNewVectorPointer)
G4bool DoesGlobalFieldExist()
G4double ComputeStep(G4FieldTrack &pFieldTrack, G4double pCurrentProposedStepLength, G4double &pNewSafety, G4VPhysicalVolume *pPhysVol=0)
G4double fTransportEndKineticEnergy
#define G4ThreadLocal
Definition: tls.hh:69
G4Material * GetMaterial() const
G4double GetVelocity() const
G4double GetProperTime() const
virtual void Initialize(const G4Track &)
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:440
G4ThreeVector GetSpin() const
G4ThreeVector GetPosition() const
void ProposeLocalTime(G4double t)
virtual void StartTracking(G4Track *)
Definition: G4VProcess.cc:102
G4double GetPDGSpin() const
G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData)
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
G4ThreeVector fTransportEndSpin
G4double GetThresholdWarningEnergy() const
G4FieldManager * FindAndSetFieldManager(G4VPhysicalVolume *pCurrentPhysVol)
G4ParticleDefinition * GetDefinition() const
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
void ProposeLastStepInVolume(G4bool flag)
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
G4double PostStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4ForceCondition *pForceCond)
#define fPreviousSftOrigin
G4double GetKineticEnergy() const
long G4long
Definition: G4Types.hh:80
G4ThreeVector fTransportEndPosition
void ProposePosition(G4double x, G4double y, G4double z)
G4double GetMagneticMoment() const
G4GPILSelection
static G4bool fUseMagneticMoment
G4double GetGlobalTime() const
static constexpr double kiloelectronvolt
Definition: G4SIunits.hh:207
G4double GetCharge() const
G4ThreeVector fTransportEndMomentumDir
#define fNewTrack
const G4ThreeVector & GetPosition() const
const G4Field * GetDetectorField() const
double mag2() const
Definition: G4Step.hh:76
void ProposeEnergy(G4double finalEnergy)
std::vector< G4ThreeVector > * GimmeTrajectoryVectorAndForgetIt() const
G4Transportation(G4int verbosityLevel=1)
void ProposeTrueStepLength(G4double truePathLength)
G4GLOB_DLL std::ostream G4cerr
G4double GetPDGMagneticMoment() const
G4TouchableHandle fCurrentTouchableHandle
static G4TransportationManager * GetTransportationManager()
G4double fCandidateEndGlobalTime
virtual G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4double pProposedMaxLength=DBL_MAX, const G4bool keepState=true)
static G4ProductionCutsTable * GetProductionCutsTable()
int G4int
Definition: G4Types.hh:78
void ProposeProperTime(G4double finalProperTime)
G4PropagatorInField * fFieldPropagator
G4int verboseLevel
Definition: G4VProcess.hh:371
void ProposeFirstStepInVolume(G4bool flag)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4bool ExitedMotherVolume() const
const G4ThreeVector & GetMomentumDir() const
void SetMaterialInTouchable(G4Material *fMaterial)
G4bool IsParticleLooping() const
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:435
G4ForceCondition
G4GLOB_DLL std::ostream G4cout
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
G4ProductionCuts * GetProductionCuts() const
G4double GetThresholdImportantEnergy() const
G4SafetyHelper * fpSafetyHelper
G4PropagatorInField * GetPropagatorInField() const
void ReportLoopingTrack(const G4Track &track, const G4Step &stepInfo, G4int numTrials, G4long noCalls, const char *methodName) const
G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
const G4Material * GetMaterial() const
G4double GetMass() const
T sqr(const T &x)
Definition: templates.hh:145
G4double GetTotalMomentum() const
G4Navigator * fLinearNavigator
void SetTouchableHandle(const G4TouchableHandle &fTouchable)
void StartTracking(G4Track *aTrack)
void SetMaterialCutsCoupleInTouchable(const G4MaterialCutsCouple *fMaterialCutsCouple)
G4TrackStatus GetTrackStatus() const
G4double fThreshold_Warning_Energy
G4int GetThresholdTrials() const
const G4ThreeVector & GetMomentumDirection() const
#define DBL_MAX
Definition: templates.hh:83
G4TransportationLogger * fpLogger
void ProposeTrackStatus(G4TrackStatus status)
G4double fThreshold_Important_Energy
G4double GetLabTimeOfFlight() const
const G4DynamicParticle * GetDynamicParticle() const
virtual void SetChargeMomentumMass(G4ChargeState particleCharge, G4double MomentumXc, G4double MassXc2)=0
static constexpr double perThousand
Definition: G4SIunits.hh:333
void SetMomentumChanged(G4bool b)
G4bool IsGravityActive() const
Definition: G4Field.hh:98
void SetGeometricallyLimitedStep()
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4VPhysicalVolume * GetVolume() const
G4double GetTotalEnergy() const