Geant4  10.01.p01
G4Transportation.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 //
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"
62 
63 #include "G4PhysicalConstants.hh"
64 #include "G4SystemOfUnits.hh"
65 #include "G4ProductionCutsTable.hh"
66 #include "G4ParticleTable.hh"
67 
68 #include "G4ChargeState.hh"
69 #include "G4EquationOfMotion.hh"
70 
71 #include "G4FieldManagerStore.hh"
72 
74 
77 //
78 // Constructor
79 
81  : G4VProcess( G4String("Transportation"), fTransportation ),
82  fTransportEndPosition( 0.0, 0.0, 0.0 ),
83  fTransportEndMomentumDir( 0.0, 0.0, 0.0 ),
84  fTransportEndKineticEnergy( 0.0 ),
85  fTransportEndSpin( 0.0, 0.0, 0.0 ),
86  fMomentumChanged(true),
87  fEndGlobalTimeComputed(false),
88  fCandidateEndGlobalTime(0.0),
89  fParticleIsLooping( false ),
90  fGeometryLimitedStep(true),
91  fPreviousSftOrigin( 0.,0.,0. ),
92  fPreviousSafety( 0.0 ),
93  // fParticleChange(),
94  fEndPointDistance( -1.0 ),
95  fThreshold_Warning_Energy( 100 * MeV ),
96  fThreshold_Important_Energy( 250 * MeV ),
97  fThresholdTrials( 10 ),
98  fNoLooperTrials( 0 ),
99  fSumEnergyKilled( 0.0 ), fMaxEnergyKilled( 0.0 ),
100  fShortStepOptimisation( false ), // Old default: true (=fast short steps)
101  fVerboseLevel( verbosity )
102 {
103  // set Process Sub Type
104  SetProcessSubType(static_cast<G4int>(TRANSPORTATION));
105  pParticleChange= &fParticleChange; // Required to conform to G4VProcess
106 
107  G4TransportationManager* transportMgr ;
108 
110 
111  fLinearNavigator = transportMgr->GetNavigatorForTracking() ;
112 
113  fFieldPropagator = transportMgr->GetPropagatorInField() ;
114 
115  fpSafetyHelper = transportMgr->GetSafetyHelper(); // New
116 
117  // Cannot determine whether a field exists here, as it would
118  // depend on the relative order of creating the detector's
119  // field and this process. That order is not guaranted.
120  // Instead later the method DoesGlobalFieldExist() is called
121 
122  static G4ThreadLocal G4TouchableHandle* pNullTouchableHandle = 0;
123  if ( !pNullTouchableHandle) { pNullTouchableHandle = new G4TouchableHandle; }
124  fCurrentTouchableHandle = *pNullTouchableHandle;
125  // Points to (G4VTouchable*) 0
126 
127 #ifdef G4VERBOSE
128  if( fVerboseLevel > 0)
129  {
130  G4cout << " G4Transportation constructor> set fShortStepOptimisation to ";
131  if ( fShortStepOptimisation ) G4cout << "true" << G4endl;
132  else G4cout << "false" << G4endl;
133  }
134 #endif
135 }
136 
138 
140 {
141  if( (fVerboseLevel > 0) && (fSumEnergyKilled > 0.0 ) )
142  {
143  G4cout << " G4Transportation: Statistics for looping particles " << G4endl;
144  G4cout << " Sum of energy of loopers killed: " << fSumEnergyKilled << G4endl;
145  G4cout << " Max energy of loopers killed: " << fMaxEnergyKilled << G4endl;
146  }
147 }
148 
150 //
151 // Responsibilities:
152 // Find whether the geometry limits the Step, and to what length
153 // Calculate the new value of the safety and return it.
154 // Store the final time, position and momentum.
155 
158  G4double, // previousStepSize
159  G4double currentMinimumStep,
160  G4double& currentSafety,
161  G4GPILSelection* selection )
162 {
163  G4double geometryStepLength= -1.0, newSafety= -1.0;
164  fParticleIsLooping = false ;
165 
166  // Initial actions moved to StartTrack()
167  // --------------------------------------
168  // Note: in case another process changes touchable handle
169  // it will be necessary to add here (for all steps)
170  // fCurrentTouchableHandle = aTrack->GetTouchableHandle();
171 
172  // GPILSelection is set to defaule value of CandidateForSelection
173  // It is a return value
174  //
175  *selection = CandidateForSelection ;
176 
177  // Get initial Energy/Momentum of the track
178  //
179  const G4DynamicParticle* pParticle = track.GetDynamicParticle() ;
180  const G4ParticleDefinition* pParticleDef = pParticle->GetDefinition() ;
181  G4ThreeVector startMomentumDir = pParticle->GetMomentumDirection() ;
182  G4ThreeVector startPosition = track.GetPosition() ;
183 
184  // G4double theTime = track.GetGlobalTime() ;
185 
186  // The Step Point safety can be limited by other geometries and/or the
187  // assumptions of any process - it's not always the geometrical safety.
188  // We calculate the starting point's isotropic safety here.
189  //
190  G4ThreeVector OriginShift = startPosition - fPreviousSftOrigin ;
191  G4double MagSqShift = OriginShift.mag2() ;
192  if( MagSqShift >= sqr(fPreviousSafety) )
193  {
194  currentSafety = 0.0 ;
195  }
196  else
197  {
198  currentSafety = fPreviousSafety - std::sqrt(MagSqShift) ;
199  }
200 
201  // Is the particle charged or has it a magnetic moment?
202  //
203  G4double particleCharge = pParticle->GetCharge() ;
204  G4double magneticMoment = pParticle->GetMagneticMoment() ;
205  G4double restMass = pParticle->GetMass() ;
206 
207  fGeometryLimitedStep = false ;
208  // fEndGlobalTimeComputed = false ;
209 
210  // There is no need to locate the current volume. It is Done elsewhere:
211  // On track construction
212  // By the tracking, after all AlongStepDoIts, in "Relocation"
213 
214  // Check if the particle has a force, EM or gravitational, exerted on it
215  //
216  G4FieldManager* fieldMgr=0;
217  G4bool fieldExertsForce = false ;
218 
219  G4bool gravityOn = false;
220  G4bool fieldExists= false; // Field is not 0 (null pointer)
221 
222  fieldMgr = fFieldPropagator->FindAndSetFieldManager( track.GetVolume() );
223  if( fieldMgr != 0 )
224  {
225  // Message the field Manager, to configure it for this track
226  fieldMgr->ConfigureForTrack( &track );
227  // Is here to allow a transition from no-field pointer
228  // to finite field (non-zero pointer).
229 
230  // If the field manager has no field ptr, the field is zero
231  // by definition ( = there is no field ! )
232  const G4Field* ptrField= fieldMgr->GetDetectorField();
233  fieldExists = (ptrField!=0) ;
234  if( fieldExists )
235  {
236  gravityOn= ptrField->IsGravityActive();
237 
238  if( (particleCharge != 0.0)
239  || (fUseMagneticMoment && (magneticMoment != 0.0) )
240  || (gravityOn && (restMass != 0.0) )
241  )
242  {
243  fieldExertsForce = fieldExists;
244  }
245  }
246  }
247  // G4cout << " G4Transport: field exerts force= " << fieldExertsForce
248  // << " fieldMgr= " << fieldMgr << G4endl;
249 
250  if( !fieldExertsForce )
251  {
252  G4double linearStepLength ;
253  if( fShortStepOptimisation && (currentMinimumStep <= currentSafety) )
254  {
255  // The Step is guaranteed to be taken
256  //
257  geometryStepLength = currentMinimumStep ;
258  fGeometryLimitedStep = false ;
259  }
260  else
261  {
262  // Find whether the straight path intersects a volume
263  //
264  linearStepLength = fLinearNavigator->ComputeStep( startPosition,
265  startMomentumDir,
266  currentMinimumStep,
267  newSafety) ;
268  // Remember last safety origin & value.
269  //
270  fPreviousSftOrigin = startPosition ;
271  fPreviousSafety = newSafety ;
272  fpSafetyHelper->SetCurrentSafety( newSafety, startPosition);
273 
274  currentSafety = newSafety ;
275 
276  fGeometryLimitedStep= (linearStepLength <= currentMinimumStep);
278  {
279  // The geometry limits the Step size (an intersection was found.)
280  geometryStepLength = linearStepLength ;
281  }
282  else
283  {
284  // The full Step is taken.
285  geometryStepLength = currentMinimumStep ;
286  }
287  }
288  fEndPointDistance = geometryStepLength ;
289 
290  // Calculate final position
291  //
292  fTransportEndPosition = startPosition+geometryStepLength*startMomentumDir ;
293 
294  // Momentum direction, energy and polarisation are unchanged by transport
295  //
296  fTransportEndMomentumDir = startMomentumDir ;
299  fParticleIsLooping = false ;
300  fMomentumChanged = false ;
301  fEndGlobalTimeComputed = false ;
302  }
303  else // A field exerts force
304  {
305  G4double momentumMagnitude = pParticle->GetTotalMomentum() ;
306  G4ThreeVector EndUnitMomentum ;
307  G4double lengthAlongCurve ;
308 
309  G4ChargeState chargeState(particleCharge, // The charge can change (dynamic)
310  magneticMoment,
311  pParticleDef->GetPDGSpin() );
312  // For insurance, could set it again
313  // chargeState.SetPDGSpin(pParticleDef->GetPDGSpin() ); // Provisionally in same object
314 
315  G4EquationOfMotion* equationOfMotion =
317  ->GetEquationOfMotion();
318 
319 // equationOfMotion->SetChargeMomentumMass( particleCharge,
320  equationOfMotion->SetChargeMomentumMass( chargeState,
321  momentumMagnitude,
322  restMass);
323 
324  G4FieldTrack aFieldTrack = G4FieldTrack( startPosition,
325  track.GetGlobalTime(), // Lab.
326  // track.GetProperTime(), // Particle rest frame
327  track.GetMomentumDirection(),
328  track.GetKineticEnergy(),
329  restMass,
330  particleCharge,
331  track.GetPolarization(),
332  pParticleDef->GetPDGMagneticMoment(),
333  0.0, // Length along track
334  pParticleDef->GetPDGSpin()
335  ) ;
336 
337  if( currentMinimumStep > 0 )
338  {
339  // Do the Transport in the field (non recti-linear)
340  //
341  lengthAlongCurve = fFieldPropagator->ComputeStep( aFieldTrack,
342  currentMinimumStep,
343  currentSafety,
344  track.GetVolume() ) ;
345  fGeometryLimitedStep= lengthAlongCurve < currentMinimumStep;
347  {
348  geometryStepLength = lengthAlongCurve ;
349  }
350  else
351  {
352  geometryStepLength = currentMinimumStep ;
353  }
354 
355  // Remember last safety origin & value.
356  //
357  fPreviousSftOrigin = startPosition ;
358  fPreviousSafety = currentSafety ;
359  fpSafetyHelper->SetCurrentSafety( currentSafety, startPosition);
360  }
361  else
362  {
363  geometryStepLength = lengthAlongCurve= 0.0 ;
364  fGeometryLimitedStep = false ;
365  }
366 
367  // Get the End-Position and End-Momentum (Dir-ection)
368  //
369  fTransportEndPosition = aFieldTrack.GetPosition() ;
370 
371  // Momentum: Magnitude and direction can be changed too now ...
372  //
373  fMomentumChanged = true ;
374  fTransportEndMomentumDir = aFieldTrack.GetMomentumDir() ;
375 
377 
379  {
380  // If the field can change energy, then the time must be integrated
381  // - so this should have been updated
382  //
384  fEndGlobalTimeComputed = true;
385 
386  // was ( fCandidateEndGlobalTime != track.GetGlobalTime() );
387  // a cleaner way is to have FieldTrack knowing whether time is updated.
388  }
389  else
390  {
391  // The energy should be unchanged by field transport,
392  // - so the time changed will be calculated elsewhere
393  //
394  fEndGlobalTimeComputed = false;
395 
396  // Check that the integration preserved the energy
397  // - and if not correct this!
398  G4double startEnergy= track.GetKineticEnergy();
400 
401  static G4ThreadLocal G4int no_inexact_steps=0, no_large_ediff;
402  G4double absEdiff = std::fabs(startEnergy- endEnergy);
403  if( absEdiff > perMillion * endEnergy )
404  {
405  no_inexact_steps++;
406  // Possible statistics keeping here ...
407  }
408  if( fVerboseLevel > 1 )
409  {
410  if( std::fabs(startEnergy- endEnergy) > perThousand * endEnergy )
411  {
412  static G4ThreadLocal G4int no_warnings= 0, warnModulo=1,
413  moduloFactor= 10;
414  no_large_ediff ++;
415  if( (no_large_ediff% warnModulo) == 0 )
416  {
417  no_warnings++;
418  G4cout << "WARNING - G4Transportation::AlongStepGetPIL() "
419  << " Energy change in Step is above 1^-3 relative value. " << G4endl
420  << " Relative change in 'tracking' step = "
421  << std::setw(15) << (endEnergy-startEnergy)/startEnergy << G4endl
422  << " Starting E= " << std::setw(12) << startEnergy / MeV << " MeV " << G4endl
423  << " Ending E= " << std::setw(12) << endEnergy / MeV << " MeV " << G4endl;
424  G4cout << " Energy has been corrected -- however, review"
425  << " field propagation parameters for accuracy." << G4endl;
426  if( (fVerboseLevel > 2 ) || (no_warnings<4) || (no_large_ediff == warnModulo * moduloFactor) )
427  {
428  G4cout << " These include EpsilonStepMax(/Min) in G4FieldManager "
429  << " which determine fractional error per step for integrated quantities. " << G4endl
430  << " Note also the influence of the permitted number of integration steps."
431  << G4endl;
432  }
433  G4cerr << "ERROR - G4Transportation::AlongStepGetPIL()" << G4endl
434  << " Bad 'endpoint'. Energy change detected"
435  << " and corrected. "
436  << " Has occurred already "
437  << no_large_ediff << " times." << G4endl;
438  if( no_large_ediff == warnModulo * moduloFactor )
439  {
440  warnModulo *= moduloFactor;
441  }
442  }
443  }
444  } // end of if (fVerboseLevel)
445 
446  // Correct the energy for fields that conserve it
447  // This - hides the integration error
448  // - but gives a better physical answer
450  }
451 
452  fTransportEndSpin = aFieldTrack.GetSpin();
454  fEndPointDistance = (fTransportEndPosition - startPosition).mag() ;
455  }
456 
457  // If we are asked to go a step length of 0, and we are on a boundary
458  // then a boundary will also limit the step -> we must flag this.
459  //
460  if( currentMinimumStep == 0.0 )
461  {
462  if( currentSafety == 0.0 ) { fGeometryLimitedStep = true; }
463  }
464 
465  // Update the safety starting from the end-point,
466  // if it will become negative at the end-point.
467  //
468  if( currentSafety < fEndPointDistance )
469  {
470  if( particleCharge != 0.0 )
471  {
472  G4double endSafety =
474  currentSafety = endSafety ;
475  fPreviousSftOrigin = fTransportEndPosition ;
476  fPreviousSafety = currentSafety ;
478 
479  // Because the Stepping Manager assumes it is from the start point,
480  // add the StepLength
481  //
482  currentSafety += fEndPointDistance ;
483 
484 #ifdef G4DEBUG_TRANSPORT
485  G4cout.precision(12) ;
486  G4cout << "***G4Transportation::AlongStepGPIL ** " << G4endl ;
487  G4cout << " Called Navigator->ComputeSafety at " << fTransportEndPosition
488  << " and it returned safety= " << endSafety << G4endl ;
489  G4cout << " Adding endpoint distance " << fEndPointDistance
490  << " to obtain pseudo-safety= " << currentSafety << G4endl ;
491  }
492  else
493  {
494  G4cout << "***G4Transportation::AlongStepGPIL ** " << G4endl ;
495  G4cout << " Avoiding call to ComputeSafety : " << G4endl;
496  G4cout << " charge = " << particleCharge << G4endl;
497  G4cout << " mag moment = " << magneticMoment << G4endl;
498 #endif
499  }
500  }
501 
502  fParticleChange.ProposeTrueStepLength(geometryStepLength) ;
503 
504  return geometryStepLength ;
505 }
506 
508 //
509 // Initialize ParticleChange (by setting all its members equal
510 // to corresponding members in G4Track)
511 
513  const G4Step& stepData )
514 {
515  static G4ThreadLocal G4int noCalls=0;
516  noCalls++;
517 
518  fParticleChange.Initialize(track) ;
519 
520  // Code for specific process
521  //
526 
528 
529  G4double deltaTime = 0.0 ;
530 
531  // Calculate Lab Time of Flight (ONLY if field Equations used it!)
532  // G4double endTime = fCandidateEndGlobalTime;
533  // G4double delta_time = endTime - startTime;
534 
535  G4double startTime = track.GetGlobalTime() ;
536 
538  {
539  // The time was not integrated .. make the best estimate possible
540  //
541  G4double initialVelocity = stepData.GetPreStepPoint()->GetVelocity();
542  G4double stepLength = track.GetStepLength();
543 
544  deltaTime= 0.0; // in case initialVelocity = 0
545  if ( initialVelocity > 0.0 ) { deltaTime = stepLength/initialVelocity; }
546 
547  fCandidateEndGlobalTime = startTime + deltaTime ;
548  fParticleChange.ProposeLocalTime( track.GetLocalTime() + deltaTime) ;
549  }
550  else
551  {
552  deltaTime = fCandidateEndGlobalTime - startTime ;
554  }
555 
556 
557  // Now Correct by Lorentz factor to get delta "proper" Time
558 
559  G4double restMass = track.GetDynamicParticle()->GetMass() ;
560  G4double deltaProperTime = deltaTime*( restMass/track.GetTotalEnergy() ) ;
561 
562  fParticleChange.ProposeProperTime(track.GetProperTime() + deltaProperTime) ;
563  //fParticleChange. ProposeTrueStepLength( track.GetStepLength() ) ;
564 
565  // If the particle is caught looping or is stuck (in very difficult
566  // boundaries) in a magnetic field (doing many steps)
567  // THEN this kills it ...
568  //
569  if ( fParticleIsLooping )
570  {
572 
573  if( (endEnergy < fThreshold_Important_Energy)
575  {
576  // Kill the looping particle
577  //
579 
580  // 'Bare' statistics
581  fSumEnergyKilled += endEnergy;
582  if( endEnergy > fMaxEnergyKilled) { fMaxEnergyKilled= endEnergy; }
583 
584 #ifdef G4VERBOSE
585  if( (fVerboseLevel > 1) &&
586  ( endEnergy > fThreshold_Warning_Energy ) )
587  {
588  G4cout << " G4Transportation is killing track that is looping or stuck "
589  << G4endl
590  << " This track has " << track.GetKineticEnergy() / MeV
591  << " MeV energy." << G4endl;
592  G4cout << " Number of trials = " << fNoLooperTrials
593  << " No of calls to AlongStepDoIt = " << noCalls
594  << G4endl;
595  }
596 #endif
597  fNoLooperTrials=0;
598  }
599  else
600  {
601  fNoLooperTrials ++;
602 #ifdef G4VERBOSE
603  if( (fVerboseLevel > 2) )
604  {
605  G4cout << " G4Transportation::AlongStepDoIt(): Particle looping - "
606  << " Number of trials = " << fNoLooperTrials
607  << " No of calls to = " << noCalls
608  << G4endl;
609  }
610 #endif
611  }
612  }
613  else
614  {
615  fNoLooperTrials=0;
616  }
617 
618  // Another (sometimes better way) is to use a user-limit maximum Step size
619  // to alleviate this problem ..
620 
621  // Introduce smooth curved trajectories to particle-change
622  //
625 
626  return &fParticleChange ;
627 }
628 
630 //
631 // This ensures that the PostStep action is always called,
632 // so that it can do the relocation if it is needed.
633 //
634 
637  G4double, // previousStepSize
638  G4ForceCondition* pForceCond )
639 {
640  *pForceCond = Forced ;
641  return DBL_MAX ; // was kInfinity ; but convention now is DBL_MAX
642 }
643 
645 //
646 
648  const G4Step& )
649 {
650  G4TouchableHandle retCurrentTouchable ; // The one to return
651  G4bool isLastStep= false;
652 
653  // Initialize ParticleChange (by setting all its members equal
654  // to corresponding members in G4Track)
655  // fParticleChange.Initialize(track) ; // To initialise TouchableChange
656 
658 
659  // If the Step was determined by the volume boundary,
660  // logically relocate the particle
661 
663  {
664  // fCurrentTouchable will now become the previous touchable,
665  // and what was the previous will be freed.
666  // (Needed because the preStepPoint can point to the previous touchable)
667 
670  LocateGlobalPointAndUpdateTouchableHandle( track.GetPosition(),
671  track.GetMomentumDirection(),
673  true ) ;
674  // Check whether the particle is out of the world volume
675  // If so it has exited and must be killed.
676  //
677  if( fCurrentTouchableHandle->GetVolume() == 0 )
678  {
680  }
681  retCurrentTouchable = fCurrentTouchableHandle ;
683 
684  // Update the Step flag which identifies the Last Step in a volume
685  isLastStep = fLinearNavigator->ExitedMotherVolume()
687 
688 #ifdef G4DEBUG_TRANSPORT
689  // Checking first implementation of flagging Last Step in Volume
690  //
693  if( ! (exiting || entering) )
694  {
695  G4cout << " Transport> : Proposed isLastStep= " << isLastStep
696  << " Exiting " << fLinearNavigator->ExitedMotherVolume()
697  << " Entering " << fLinearNavigator->EnteredDaughterVolume()
698  << G4endl;
699  }
700 #endif
701  }
702  else // fGeometryLimitedStep is false
703  {
704  // This serves only to move the Navigator's location
705  //
707 
708  // The value of the track's current Touchable is retained.
709  // (and it must be correct because we must use it below to
710  // overwrite the (unset) one in particle change)
711  // It must be fCurrentTouchable too ??
712  //
714  retCurrentTouchable = track.GetTouchableHandle() ;
715 
716  isLastStep= false;
717 
718 #ifdef G4DEBUG_TRANSPORT
719  // Checking first implementation of flagging Last Step in Volume
720  //
721  G4cout << " Transport> Proposed isLastStep= " << isLastStep
722  << " Geometry did not limit step. " << G4endl;
723 #endif
724  } // endif ( fGeometryLimitedStep )
725 
727 
728  const G4VPhysicalVolume* pNewVol = retCurrentTouchable->GetVolume() ;
729  const G4Material* pNewMaterial = 0 ;
730  const G4VSensitiveDetector* pNewSensitiveDetector = 0 ;
731 
732  if( pNewVol != 0 )
733  {
734  pNewMaterial= pNewVol->GetLogicalVolume()->GetMaterial();
735  pNewSensitiveDetector= pNewVol->GetLogicalVolume()->GetSensitiveDetector();
736  }
737 
738  // ( <const_cast> pNewMaterial ) ;
739  // ( <const_cast> pNewSensitiveDetector) ;
740 
743 
744  const G4MaterialCutsCouple* pNewMaterialCutsCouple = 0;
745  if( pNewVol != 0 )
746  {
747  pNewMaterialCutsCouple=pNewVol->GetLogicalVolume()->GetMaterialCutsCouple();
748  }
749 
750  if( pNewVol!=0 && pNewMaterialCutsCouple!=0 && pNewMaterialCutsCouple->GetMaterial()!=pNewMaterial )
751  {
752  // for parametrized volume
753  //
754  pNewMaterialCutsCouple =
756  ->GetMaterialCutsCouple(pNewMaterial,
757  pNewMaterialCutsCouple->GetProductionCuts());
758  }
759  fParticleChange.SetMaterialCutsCoupleInTouchable( pNewMaterialCutsCouple );
760 
761  // temporarily until Get/Set Material of ParticleChange,
762  // and StepPoint can be made const.
763  // Set the touchable in ParticleChange
764  // this must always be done because the particle change always
765  // uses this value to overwrite the current touchable pointer.
766  //
767  fParticleChange.SetTouchableHandle(retCurrentTouchable) ;
768 
769  return &fParticleChange ;
770 }
771 
772 // New method takes over the responsibility to reset the state of G4Transportation
773 // object at the start of a new track or the resumption of a suspended track.
774 
775 void
777 {
779 
780  // The actions here are those that were taken in AlongStepGPIL
781  // when track.GetCurrentStepNumber()==1
782 
783  // reset safety value and center
784  //
785  fPreviousSafety = 0.0 ;
786  fPreviousSftOrigin = G4ThreeVector(0.,0.,0.) ;
787 
788  // reset looping counter -- for motion in field
789  fNoLooperTrials= 0;
790  // Must clear this state .. else it depends on last track's value
791  // --> a better solution would set this from state of suspended track TODO ?
792  // Was if( aTrack->GetCurrentStepNumber()==1 ) { .. }
793 
794  // ChordFinder reset internal state
795  //
796  if( DoesGlobalFieldExist() )
797  {
799  // Resets all state of field propagator class (ONLY)
800  // including safety values (in case of overlaps and to wipe for first track).
801 
802  // G4ChordFinder* chordF= fFieldPropagator->GetChordFinder();
803  // if( chordF ) chordF->ResetStepEstimate();
804  }
805 
806  // Make sure to clear the chord finders of all fields (ie managers)
807  //
809  fieldMgrStore->ClearAllChordFindersState();
810 
811  // Update the current touchable handle (from the track's)
812  //
814 }
815 
818 {
819  G4bool lastValue= fUseMagneticMoment;
820  fUseMagneticMoment= useMoment;
822  return lastValue;
823 }
void SetMaterialInTouchable(G4Material *fMaterial)
void SetTouchableHandle(const G4TouchableHandle &fTouchable)
const G4ThreeVector & GetPolarization() const
G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData)
static const double MeV
Definition: G4SIunits.hh:193
G4SafetyHelper * GetSafetyHelper() const
G4double GetLocalTime() const
G4PropagatorInField * fFieldPropagator
void SetMaterialCutsCoupleInTouchable(const G4MaterialCutsCouple *fMaterialCutsCouple)
virtual void SetChargeMomentumMass(G4ChargeState particleCharge, G4double MomentumXc, G4double MassXc2)=0
G4double GetProperTime() const
std::vector< G4ThreeVector > * GimmeTrajectoryVectorAndForgetIt() const
CLHEP::Hep3Vector G4ThreeVector
G4double fTransportEndKineticEnergy
void StartTracking(G4Track *aTrack)
G4double GetKineticEnergy() const
const G4DynamicParticle * GetDynamicParticle() const
G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData)
const G4MagIntegratorStepper * GetStepper() const
G4Material * GetMaterial() const
const G4ThreeVector & GetPosition() const
const G4ThreeVector & GetMomentumDir() const
G4TrackStatus GetTrackStatus() const
G4ThreeVector GetSpin() const
virtual G4double ComputeStep(const G4ThreeVector &pGlobalPoint, const G4ThreeVector &pDirection, const G4double pCurrentProposedStepLength, G4double &pNewSafety)
Definition: G4Navigator.cc:744
G4Navigator * GetNavigatorForTracking() const
void ProposePolarization(G4double Px, G4double Py, G4double Pz)
G4ParticleDefinition * GetDefinition() const
static const double perThousand
Definition: G4SIunits.hh:297
G4ThreeVector fTransportEndPosition
G4double GetVelocity() const
#define G4ThreadLocal
Definition: tls.hh:89
G4SafetyHelper * fpSafetyHelper
int G4int
Definition: G4Types.hh:78
void ProposePosition(G4double x, G4double y, G4double z)
G4ReferenceCountedHandle< G4VTouchable > G4TouchableHandle
G4double GetTotalMomentum() const
G4StepPoint * GetPreStepPoint() const
virtual void Initialize(const G4Track &)
virtual void StartTracking(G4Track *)
Definition: G4VProcess.cc:101
virtual void ConfigureForTrack(const G4Track *)
void SetGeometricallyLimitedStep()
G4double GetKineticEnergy() const
G4ThreeVector GetPosition() const
G4GLOB_DLL std::ostream G4cout
G4double fThreshold_Important_Energy
G4double GetMass() const
bool G4bool
Definition: G4Types.hh:79
const G4ThreeVector & GetMomentumDirection() const
G4double GetCharge() const
void ProposeTrueStepLength(G4double truePathLength)
G4bool DoesGlobalFieldExist()
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
G4bool EnteredDaughterVolume() const
G4bool DoesFieldChangeEnergy() const
Definition: G4Step.hh:76
G4double GetGlobalTime() const
G4ThreeVector fPreviousSftOrigin
G4FieldManager * FindAndSetFieldManager(G4VPhysicalVolume *pCurrentPhysVol)
G4TouchableHandle fCurrentTouchableHandle
G4Transportation(G4int verbosityLevel=1)
const G4TouchableHandle & GetTouchableHandle() const
static G4TransportationManager * GetTransportationManager()
G4double PostStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4ForceCondition *pForceCond)
static const double perMillion
Definition: G4SIunits.hh:298
G4double ComputeStep(G4FieldTrack &pFieldTrack, G4double pCurrentProposedStepLength, G4double &pNewSafety, G4VPhysicalVolume *pPhysVol=0)
static G4ProductionCutsTable * GetProductionCutsTable()
const G4ThreeVector & GetMomentumDirection() const
G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
G4LogicalVolume * GetLogicalVolume() const
void ProposeProperTime(G4double finalProperTime)
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
static G4FieldManagerStore * GetInstance()
void SetSensitiveDetectorInTouchable(G4VSensitiveDetector *fSensitiveDetector)
#define fPreviousSftOrigin
G4bool IsParticleLooping() const
G4double GetLabTimeOfFlight() const
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:44
G4double fThreshold_Warning_Energy
void SetPointerToVectorOfAuxiliaryPoints(std::vector< G4ThreeVector > *theNewVectorPointer)
G4bool IsGravityActive() const
Definition: G4Field.hh:98
void ProposeGlobalTime(G4double t)
G4ChordFinder * GetChordFinder()
G4FieldManager * GetCurrentFieldManager()
#define fPreviousSafety
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
void ProposeEnergy(G4double finalEnergy)
G4VPhysicalVolume * GetVolume() const
G4double GetTotalEnergy() const
G4double GetPDGSpin() const
void SetCurrentSafety(G4double val, const G4ThreeVector &pos)
#define G4endl
Definition: G4ios.hh:61
G4bool ExitedMotherVolume() const
void ProposeLastStepInVolume(G4bool flag)
G4ThreeVector fTransportEndSpin
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
T sqr(const T &x)
Definition: templates.hh:145
G4double GetPDGMagneticMoment() const
G4Navigator * fLinearNavigator
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
const G4Field * GetDetectorField() const
virtual G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4double pProposedMaxLength=DBL_MAX, const G4bool keepState=true)
G4ForceCondition
G4MagInt_Driver * GetIntegrationDriver()
void SetMomentumChanged(G4bool b)
G4PropagatorInField * GetPropagatorInField() const
G4double fCandidateEndGlobalTime
G4ProductionCuts * GetProductionCuts() const
#define DBL_MAX
Definition: templates.hh:83
G4VSensitiveDetector * GetSensitiveDetector() const
G4double GetMagneticMoment() const
G4ParticleChangeForTransport fParticleChange
void ProposeLocalTime(G4double t)
G4double GetStepLength() const
virtual void LocateGlobalPointWithinVolume(const G4ThreeVector &position)
Definition: G4Navigator.cc:579
G4GPILSelection
const G4Material * GetMaterial() const
G4ThreeVector fTransportEndMomentumDir
G4GLOB_DLL std::ostream G4cerr
static G4bool fUseMagneticMoment
static G4bool EnableUseMagneticMoment(G4bool useMoment=true)