Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 // 28 Oct 2011, P.Gumpl./J.Ap: Detect gravity field, use magnetic moment
44 // 20 Nov 2008, J.Apostolakis: Push safety to helper - after ComputeSafety
45 // 9 Nov 2007, J.Apostolakis: Flag for short steps, push safety to helper
46 // 19 Jan 2006, P.MoraDeFreitas: Fix for suspended tracks (StartTracking)
47 // 11 Aug 2004, M.Asai: Add G4VSensitiveDetector* for updating stepPoint.
48 // 21 June 2003, J.Apostolakis: Calling field manager with
49 // track, to enable it to configure its accuracy
50 // 13 May 2003, J.Apostolakis: Zero field areas now taken into
51 // account correclty in all cases (thanks to W Pokorski).
52 // 29 June 2001, J.Apostolakis, D.Cote-Ahern, P.Gumplinger:
53 // correction for spin tracking
54 // 20 Febr 2001, J.Apostolakis: update for new FieldTrack
55 // 22 Sept 2000, V.Grichine: update of Kinetic Energy
56 // Created: 19 March 1997, J. Apostolakis
57 // =======================================================================
58 
59 #include "G4Transportation.hh"
61 
62 #include "G4PhysicalConstants.hh"
63 #include "G4SystemOfUnits.hh"
64 #include "G4ProductionCutsTable.hh"
65 #include "G4ParticleTable.hh"
66 #include "G4ChordFinder.hh"
67 #include "G4SafetyHelper.hh"
68 #include "G4FieldManagerStore.hh"
69 
71 
73 //
74 // Constructor
75 
77  : G4VProcess( G4String("Transportation"), fTransportation ),
78  fTransportEndPosition( 0.0, 0.0, 0.0 ),
79  fTransportEndMomentumDir( 0.0, 0.0, 0.0 ),
80  fTransportEndKineticEnergy( 0.0 ),
81  fTransportEndSpin( 0.0, 0.0, 0.0 ),
82  fMomentumChanged(true),
83  fEndGlobalTimeComputed(false),
84  fCandidateEndGlobalTime(0.0),
85  fParticleIsLooping( false ),
86  fGeometryLimitedStep(true),
87  fPreviousSftOrigin( 0.,0.,0. ),
88  fPreviousSafety( 0.0 ),
89  // fParticleChange(),
90  fEndPointDistance( -1.0 ),
91  fThreshold_Warning_Energy( 100 * MeV ),
92  fThreshold_Important_Energy( 250 * MeV ),
93  fThresholdTrials( 10 ),
94  fNoLooperTrials( 0 ),
95  fSumEnergyKilled( 0.0 ), fMaxEnergyKilled( 0.0 ),
96  fShortStepOptimisation( false ), // Old default: true (=fast short steps)
97  fUseMagneticMoment( false ),
98  fVerboseLevel( verbosity )
99 {
100  // set Process Sub Type
101  SetProcessSubType(static_cast<G4int>(TRANSPORTATION));
102 
103  G4TransportationManager* transportMgr ;
104 
106 
107  fLinearNavigator = transportMgr->GetNavigatorForTracking() ;
108 
109  fFieldPropagator = transportMgr->GetPropagatorInField() ;
110 
111  fpSafetyHelper = transportMgr->GetSafetyHelper(); // New
112 
113  // Cannot determine whether a field exists here, as it would
114  // depend on the relative order of creating the detector's
115  // field and this process. That order is not guaranted.
116  // Instead later the method DoesGlobalFieldExist() is called
117 
118  static G4TouchableHandle nullTouchableHandle; // Points to (G4VTouchable*) 0
119  fCurrentTouchableHandle = nullTouchableHandle;
120 
121 #ifdef G4VERBOSE
122  if( fVerboseLevel > 0)
123  {
124  G4cout << " G4Transportation constructor> set fShortStepOptimisation to ";
125  if ( fShortStepOptimisation ) G4cout << "true" << G4endl;
126  else G4cout << "false" << G4endl;
127  }
128 #endif
129 }
130 
132 
134 {
135  if( (fVerboseLevel > 0) && (fSumEnergyKilled > 0.0 ) )
136  {
137  G4cout << " G4Transportation: Statistics for looping particles " << G4endl;
138  G4cout << " Sum of energy of loopers killed: " << fSumEnergyKilled << G4endl;
139  G4cout << " Max energy of loopers killed: " << fMaxEnergyKilled << G4endl;
140  }
141 }
142 
144 //
145 // Responsibilities:
146 // Find whether the geometry limits the Step, and to what length
147 // Calculate the new value of the safety and return it.
148 // Store the final time, position and momentum.
149 
152  G4double, // previousStepSize
153  G4double currentMinimumStep,
154  G4double& currentSafety,
155  G4GPILSelection* selection )
156 {
157  G4double geometryStepLength= -1.0, newSafety= -1.0;
158  fParticleIsLooping = false ;
159 
160  // Initial actions moved to StartTrack()
161  // --------------------------------------
162  // Note: in case another process changes touchable handle
163  // it will be necessary to add here (for all steps)
164  // fCurrentTouchableHandle = aTrack->GetTouchableHandle();
165 
166  // GPILSelection is set to defaule value of CandidateForSelection
167  // It is a return value
168  //
169  *selection = CandidateForSelection ;
170 
171  // Get initial Energy/Momentum of the track
172  //
173  const G4DynamicParticle* pParticle = track.GetDynamicParticle() ;
174  const G4ParticleDefinition* pParticleDef = pParticle->GetDefinition() ;
175  G4ThreeVector startMomentumDir = pParticle->GetMomentumDirection() ;
176  G4ThreeVector startPosition = track.GetPosition() ;
177 
178  // G4double theTime = track.GetGlobalTime() ;
179 
180  // The Step Point safety can be limited by other geometries and/or the
181  // assumptions of any process - it's not always the geometrical safety.
182  // We calculate the starting point's isotropic safety here.
183  //
184  G4ThreeVector OriginShift = startPosition - fPreviousSftOrigin ;
185  G4double MagSqShift = OriginShift.mag2() ;
186  if( MagSqShift >= sqr(fPreviousSafety) )
187  {
188  currentSafety = 0.0 ;
189  }
190  else
191  {
192  currentSafety = fPreviousSafety - std::sqrt(MagSqShift) ;
193  }
194 
195  // Is the particle charged or has it a magnetic moment?
196  //
197  G4double particleCharge = pParticle->GetCharge() ;
198  G4double magneticMoment = pParticle->GetMagneticMoment() ;
199  G4double restMass = pParticleDef->GetPDGMass() ;
200 
201  fGeometryLimitedStep = false ;
202  // fEndGlobalTimeComputed = false ;
203 
204  // There is no need to locate the current volume. It is Done elsewhere:
205  // On track construction
206  // By the tracking, after all AlongStepDoIts, in "Relocation"
207 
208  // Check if the particle has a force, EM or gravitational, exerted on it
209  //
210  G4FieldManager* fieldMgr=0;
211  G4bool fieldExertsForce = false ;
212 
213  G4bool gravityOn = false;
214  G4bool fieldExists= false; // Field is not 0 (null pointer)
215 
216  fieldMgr = fFieldPropagator->FindAndSetFieldManager( track.GetVolume() );
217  if( fieldMgr != 0 )
218  {
219  // Message the field Manager, to configure it for this track
220  fieldMgr->ConfigureForTrack( &track );
221  // Is here to allow a transition from no-field pointer
222  // to finite field (non-zero pointer).
223 
224  // If the field manager has no field ptr, the field is zero
225  // by definition ( = there is no field ! )
226  const G4Field* ptrField= fieldMgr->GetDetectorField();
227  fieldExists = (ptrField!=0) ;
228  if( fieldExists )
229  {
230  gravityOn= ptrField->IsGravityActive();
231 
232  if( (particleCharge != 0.0)
233  || (fUseMagneticMoment && (magneticMoment != 0.0) )
234  || (gravityOn && (restMass != 0.0) )
235  )
236  {
237  fieldExertsForce = fieldExists;
238  }
239  }
240  }
241  // G4cout << " G4Transport: field exerts force= " << fieldExertsForce
242  // << " fieldMgr= " << fieldMgr << G4endl;
243 
244  if( !fieldExertsForce )
245  {
246  G4double linearStepLength ;
247  if( fShortStepOptimisation && (currentMinimumStep <= currentSafety) )
248  {
249  // The Step is guaranteed to be taken
250  //
251  geometryStepLength = currentMinimumStep ;
252  fGeometryLimitedStep = false ;
253  }
254  else
255  {
256  // Find whether the straight path intersects a volume
257  //
258  linearStepLength = fLinearNavigator->ComputeStep( startPosition,
259  startMomentumDir,
260  currentMinimumStep,
261  newSafety) ;
262  // Remember last safety origin & value.
263  //
264  fPreviousSftOrigin = startPosition ;
265  fPreviousSafety = newSafety ;
266  fpSafetyHelper->SetCurrentSafety( newSafety, startPosition);
267 
268  currentSafety = newSafety ;
269 
270  fGeometryLimitedStep= (linearStepLength <= currentMinimumStep);
271  if( fGeometryLimitedStep )
272  {
273  // The geometry limits the Step size (an intersection was found.)
274  geometryStepLength = linearStepLength ;
275  }
276  else
277  {
278  // The full Step is taken.
279  geometryStepLength = currentMinimumStep ;
280  }
281  }
282  fEndPointDistance = geometryStepLength ;
283 
284  // Calculate final position
285  //
286  fTransportEndPosition = startPosition+geometryStepLength*startMomentumDir ;
287 
288  // Momentum direction, energy and polarisation are unchanged by transport
289  //
290  fTransportEndMomentumDir = startMomentumDir ;
291  fTransportEndKineticEnergy = track.GetKineticEnergy() ;
292  fTransportEndSpin = track.GetPolarization();
293  fParticleIsLooping = false ;
294  fMomentumChanged = false ;
295  fEndGlobalTimeComputed = false ;
296  }
297  else // A field exerts force
298  {
299  G4double momentumMagnitude = pParticle->GetTotalMomentum() ;
300  G4ThreeVector EndUnitMomentum ;
301  G4double lengthAlongCurve ;
302 
303  fFieldPropagator->SetChargeMomentumMass( particleCharge, // in e+ units
304  momentumMagnitude, // in Mev/c
305  restMass ) ;
306 
307  G4ThreeVector spin = track.GetPolarization() ;
308  G4FieldTrack aFieldTrack = G4FieldTrack( startPosition,
309  track.GetMomentumDirection(),
310  0.0,
311  track.GetKineticEnergy(),
312  restMass,
313  track.GetVelocity(),
314  track.GetGlobalTime(), // Lab.
315  track.GetProperTime(), // Part.
316  &spin ) ;
317  if( currentMinimumStep > 0 )
318  {
319  // Do the Transport in the field (non recti-linear)
320  //
321  lengthAlongCurve = fFieldPropagator->ComputeStep( aFieldTrack,
322  currentMinimumStep,
323  currentSafety,
324  track.GetVolume() ) ;
325  fGeometryLimitedStep= lengthAlongCurve < currentMinimumStep;
326  if( fGeometryLimitedStep )
327  {
328  geometryStepLength = lengthAlongCurve ;
329  }
330  else
331  {
332  geometryStepLength = currentMinimumStep ;
333  }
334 
335  // Remember last safety origin & value.
336  //
337  fPreviousSftOrigin = startPosition ;
338  fPreviousSafety = currentSafety ;
339  fpSafetyHelper->SetCurrentSafety( currentSafety, startPosition);
340  }
341  else
342  {
343  geometryStepLength = lengthAlongCurve= 0.0 ;
344  fGeometryLimitedStep = false ;
345  }
346 
347  // Get the End-Position and End-Momentum (Dir-ection)
348  //
349  fTransportEndPosition = aFieldTrack.GetPosition() ;
350 
351  // Momentum: Magnitude and direction can be changed too now ...
352  //
353  fMomentumChanged = true ;
354  fTransportEndMomentumDir = aFieldTrack.GetMomentumDir() ;
355 
356  fTransportEndKineticEnergy = aFieldTrack.GetKineticEnergy() ;
357 
358  if( fFieldPropagator->GetCurrentFieldManager()->DoesFieldChangeEnergy() )
359  {
360  // If the field can change energy, then the time must be integrated
361  // - so this should have been updated
362  //
363  fCandidateEndGlobalTime = aFieldTrack.GetLabTimeOfFlight();
364  fEndGlobalTimeComputed = true;
365 
366  // was ( fCandidateEndGlobalTime != track.GetGlobalTime() );
367  // a cleaner way is to have FieldTrack knowing whether time is updated.
368  }
369  else
370  {
371  // The energy should be unchanged by field transport,
372  // - so the time changed will be calculated elsewhere
373  //
374  fEndGlobalTimeComputed = false;
375 
376  // Check that the integration preserved the energy
377  // - and if not correct this!
378  G4double startEnergy= track.GetKineticEnergy();
379  G4double endEnergy= fTransportEndKineticEnergy;
380 
381  static G4int no_inexact_steps=0, no_large_ediff;
382  G4double absEdiff = std::fabs(startEnergy- endEnergy);
383  if( absEdiff > perMillion * endEnergy )
384  {
385  no_inexact_steps++;
386  // Possible statistics keeping here ...
387  }
388  if( fVerboseLevel > 1 )
389  {
390  if( std::fabs(startEnergy- endEnergy) > perThousand * endEnergy )
391  {
392  static G4int no_warnings= 0, warnModulo=1, moduloFactor= 10;
393  no_large_ediff ++;
394  if( (no_large_ediff% warnModulo) == 0 )
395  {
396  no_warnings++;
397  G4cout << "WARNING - G4Transportation::AlongStepGetPIL() "
398  << " Energy change in Step is above 1^-3 relative value. " << G4endl
399  << " Relative change in 'tracking' step = "
400  << std::setw(15) << (endEnergy-startEnergy)/startEnergy << G4endl
401  << " Starting E= " << std::setw(12) << startEnergy / MeV << " MeV " << G4endl
402  << " Ending E= " << std::setw(12) << endEnergy / MeV << " MeV " << G4endl;
403  G4cout << " Energy has been corrected -- however, review"
404  << " field propagation parameters for accuracy." << G4endl;
405  if( (fVerboseLevel > 2 ) || (no_warnings<4) || (no_large_ediff == warnModulo * moduloFactor) )
406  {
407  G4cout << " These include EpsilonStepMax(/Min) in G4FieldManager "
408  << " which determine fractional error per step for integrated quantities. " << G4endl
409  << " Note also the influence of the permitted number of integration steps."
410  << G4endl;
411  }
412  G4cerr << "ERROR - G4Transportation::AlongStepGetPIL()" << G4endl
413  << " Bad 'endpoint'. Energy change detected"
414  << " and corrected. "
415  << " Has occurred already "
416  << no_large_ediff << " times." << G4endl;
417  if( no_large_ediff == warnModulo * moduloFactor )
418  {
419  warnModulo *= moduloFactor;
420  }
421  }
422  }
423  } // end of if (fVerboseLevel)
424 
425  // Correct the energy for fields that conserve it
426  // This - hides the integration error
427  // - but gives a better physical answer
428  fTransportEndKineticEnergy= track.GetKineticEnergy();
429  }
430 
431  fTransportEndSpin = aFieldTrack.GetSpin();
432  fParticleIsLooping = fFieldPropagator->IsParticleLooping() ;
433  fEndPointDistance = (fTransportEndPosition - startPosition).mag() ;
434  }
435 
436  // If we are asked to go a step length of 0, and we are on a boundary
437  // then a boundary will also limit the step -> we must flag this.
438  //
439  if( currentMinimumStep == 0.0 )
440  {
441  if( currentSafety == 0.0 ) { fGeometryLimitedStep = true; }
442  }
443 
444  // Update the safety starting from the end-point,
445  // if it will become negative at the end-point.
446  //
447  if( currentSafety < fEndPointDistance )
448  {
449  if( particleCharge != 0.0 )
450  {
451  G4double endSafety =
452  fLinearNavigator->ComputeSafety( fTransportEndPosition) ;
453  currentSafety = endSafety ;
454  fPreviousSftOrigin = fTransportEndPosition ;
455  fPreviousSafety = currentSafety ;
456  fpSafetyHelper->SetCurrentSafety( currentSafety, fTransportEndPosition);
457 
458  // Because the Stepping Manager assumes it is from the start point,
459  // add the StepLength
460  //
461  currentSafety += fEndPointDistance ;
462 
463 #ifdef G4DEBUG_TRANSPORT
464  G4cout.precision(12) ;
465  G4cout << "***G4Transportation::AlongStepGPIL ** " << G4endl ;
466  G4cout << " Called Navigator->ComputeSafety at " << fTransportEndPosition
467  << " and it returned safety= " << endSafety << G4endl ;
468  G4cout << " Adding endpoint distance " << fEndPointDistance
469  << " to obtain pseudo-safety= " << currentSafety << G4endl ;
470  }
471  else
472  {
473  G4cout << "***G4Transportation::AlongStepGPIL ** " << G4endl ;
474  G4cout << " Avoiding call to ComputeSafety : " << G4endl;
475  G4cout << " charge = " << particleCharge << G4endl;
476  G4cout << " mag moment = " << magneticMoment << G4endl;
477 #endif
478  }
479  }
480 
481  fParticleChange.ProposeTrueStepLength(geometryStepLength) ;
482 
483  return geometryStepLength ;
484 }
485 
487 //
488 // Initialize ParticleChange (by setting all its members equal
489 // to corresponding members in G4Track)
490 
492  const G4Step& stepData )
493 {
494  static G4int noCalls=0;
495  noCalls++;
496 
497  fParticleChange.Initialize(track) ;
498 
499  // Code for specific process
500  //
501  fParticleChange.ProposePosition(fTransportEndPosition) ;
502  fParticleChange.ProposeMomentumDirection(fTransportEndMomentumDir) ;
503  fParticleChange.ProposeEnergy(fTransportEndKineticEnergy) ;
504  fParticleChange.SetMomentumChanged(fMomentumChanged) ;
505 
506  fParticleChange.ProposePolarization(fTransportEndSpin);
507 
508  G4double deltaTime = 0.0 ;
509 
510  // Calculate Lab Time of Flight (ONLY if field Equations used it!)
511  // G4double endTime = fCandidateEndGlobalTime;
512  // G4double delta_time = endTime - startTime;
513 
514  G4double startTime = track.GetGlobalTime() ;
515 
516  if (!fEndGlobalTimeComputed)
517  {
518  // The time was not integrated .. make the best estimate possible
519  //
520  G4double initialVelocity = stepData.GetPreStepPoint()->GetVelocity();
521  G4double stepLength = track.GetStepLength();
522 
523  deltaTime= 0.0; // in case initialVelocity = 0
524  if ( initialVelocity > 0.0 ) { deltaTime = stepLength/initialVelocity; }
525 
526  fCandidateEndGlobalTime = startTime + deltaTime ;
527  fParticleChange.ProposeLocalTime( track.GetLocalTime() + deltaTime) ;
528  }
529  else
530  {
531  deltaTime = fCandidateEndGlobalTime - startTime ;
532  fParticleChange.ProposeGlobalTime( fCandidateEndGlobalTime ) ;
533  }
534 
535 
536  // Now Correct by Lorentz factor to get delta "proper" Time
537 
538  G4double restMass = track.GetDynamicParticle()->GetMass() ;
539  G4double deltaProperTime = deltaTime*( restMass/track.GetTotalEnergy() ) ;
540 
541  fParticleChange.ProposeProperTime(track.GetProperTime() + deltaProperTime) ;
542  //fParticleChange. ProposeTrueStepLength( track.GetStepLength() ) ;
543 
544  // If the particle is caught looping or is stuck (in very difficult
545  // boundaries) in a magnetic field (doing many steps)
546  // THEN this kills it ...
547  //
548  if ( fParticleIsLooping )
549  {
550  G4double endEnergy= fTransportEndKineticEnergy;
551 
552  if( (endEnergy < fThreshold_Important_Energy)
553  || (fNoLooperTrials >= fThresholdTrials ) )
554  {
555  // Kill the looping particle
556  //
557  fParticleChange.ProposeTrackStatus( fStopAndKill ) ;
558 
559  // 'Bare' statistics
560  fSumEnergyKilled += endEnergy;
561  if( endEnergy > fMaxEnergyKilled) { fMaxEnergyKilled= endEnergy; }
562 
563 #ifdef G4VERBOSE
564  if( (fVerboseLevel > 1) &&
565  ( endEnergy > fThreshold_Warning_Energy ) )
566  {
567  G4cout << " G4Transportation 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  {
580  fNoLooperTrials ++;
581 #ifdef G4VERBOSE
582  if( (fVerboseLevel > 2) )
583  {
584  G4cout << " G4Transportation::AlongStepDoIt(): Particle looping - "
585  << " Number of trials = " << fNoLooperTrials
586  << " No of calls to = " << noCalls
587  << G4endl;
588  }
589 #endif
590  }
591  }
592  else
593  {
594  fNoLooperTrials=0;
595  }
596 
597  // Another (sometimes better way) is to use a user-limit maximum Step size
598  // to alleviate this problem ..
599 
600  // Introduce smooth curved trajectories to particle-change
601  //
603  (fFieldPropagator->GimmeTrajectoryVectorAndForgetIt() );
604 
605  return &fParticleChange ;
606 }
607 
609 //
610 // This ensures that the PostStep action is always called,
611 // so that it can do the relocation if it is needed.
612 //
613 
616  G4double, // previousStepSize
617  G4ForceCondition* pForceCond )
618 {
619  *pForceCond = Forced ;
620  return DBL_MAX ; // was kInfinity ; but convention now is DBL_MAX
621 }
622 
624 //
625 
627  const G4Step& )
628 {
629  G4TouchableHandle retCurrentTouchable ; // The one to return
630  G4bool isLastStep= false;
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  // Update the Step flag which identifies the Last Step in a volume
664  isLastStep = fLinearNavigator->ExitedMotherVolume()
665  | fLinearNavigator->EnteredDaughterVolume() ;
666 
667 #ifdef G4DEBUG_TRANSPORT
668  // Checking first implementation of flagging Last Step in Volume
669  //
670  G4bool exiting = fLinearNavigator->ExitedMotherVolume();
671  G4bool entering = fLinearNavigator->EnteredDaughterVolume();
672  if( ! (exiting || entering) )
673  {
674  G4cout << " Transport> : Proposed isLastStep= " << isLastStep
675  << " Exiting " << fLinearNavigator->ExitedMotherVolume()
676  << " Entering " << fLinearNavigator->EnteredDaughterVolume()
677  << G4endl;
678  }
679 #endif
680  }
681  else // fGeometryLimitedStep is false
682  {
683  // This serves only to move the Navigator's location
684  //
685  fLinearNavigator->LocateGlobalPointWithinVolume( track.GetPosition() ) ;
686 
687  // The value of the track's current Touchable is retained.
688  // (and it must be correct because we must use it below to
689  // overwrite the (unset) one in particle change)
690  // It must be fCurrentTouchable too ??
691  //
692  fParticleChange.SetTouchableHandle( track.GetTouchableHandle() ) ;
693  retCurrentTouchable = track.GetTouchableHandle() ;
694 
695  isLastStep= false;
696 
697 #ifdef G4DEBUG_TRANSPORT
698  // Checking first implementation of flagging Last Step in Volume
699  //
700  G4cout << " Transport> Proposed isLastStep= " << isLastStep
701  << " Geometry did not limit step. " << G4endl;
702 #endif
703  } // endif ( fGeometryLimitedStep )
704 
705  fParticleChange.ProposeLastStepInVolume(isLastStep);
706 
707  const G4VPhysicalVolume* pNewVol = retCurrentTouchable->GetVolume() ;
708  const G4Material* pNewMaterial = 0 ;
709  const G4VSensitiveDetector* pNewSensitiveDetector = 0 ;
710 
711  if( pNewVol != 0 )
712  {
713  pNewMaterial= pNewVol->GetLogicalVolume()->GetMaterial();
714  pNewSensitiveDetector= pNewVol->GetLogicalVolume()->GetSensitiveDetector();
715  }
716 
717  // ( <const_cast> pNewMaterial ) ;
718  // ( <const_cast> pNewSensitiveDetector) ;
719 
720  fParticleChange.SetMaterialInTouchable( (G4Material *) pNewMaterial ) ;
721  fParticleChange.SetSensitiveDetectorInTouchable( (G4VSensitiveDetector *) pNewSensitiveDetector ) ;
722 
723  const G4MaterialCutsCouple* pNewMaterialCutsCouple = 0;
724  if( pNewVol != 0 )
725  {
726  pNewMaterialCutsCouple=pNewVol->GetLogicalVolume()->GetMaterialCutsCouple();
727  }
728 
729  if( pNewVol!=0 && pNewMaterialCutsCouple!=0 && pNewMaterialCutsCouple->GetMaterial()!=pNewMaterial )
730  {
731  // for parametrized volume
732  //
733  pNewMaterialCutsCouple =
735  ->GetMaterialCutsCouple(pNewMaterial,
736  pNewMaterialCutsCouple->GetProductionCuts());
737  }
738  fParticleChange.SetMaterialCutsCoupleInTouchable( pNewMaterialCutsCouple );
739 
740  // temporarily until Get/Set Material of ParticleChange,
741  // and StepPoint can be made const.
742  // Set the touchable in ParticleChange
743  // this must always be done because the particle change always
744  // uses this value to overwrite the current touchable pointer.
745  //
746  fParticleChange.SetTouchableHandle(retCurrentTouchable) ;
747 
748  return &fParticleChange ;
749 }
750 
751 // New method takes over the responsibility to reset the state of G4Transportation
752 // object at the start of a new track or the resumption of a suspended track.
753 
754 void
756 {
758 
759  // The actions here are those that were taken in AlongStepGPIL
760  // when track.GetCurrentStepNumber()==1
761 
762  // reset safety value and center
763  //
764  fPreviousSafety = 0.0 ;
765  fPreviousSftOrigin = G4ThreeVector(0.,0.,0.) ;
766 
767  // reset looping counter -- for motion in field
768  fNoLooperTrials= 0;
769  // Must clear this state .. else it depends on last track's value
770  // --> a better solution would set this from state of suspended track TODO ?
771  // Was if( aTrack->GetCurrentStepNumber()==1 ) { .. }
772 
773  // ChordFinder reset internal state
774  //
775  if( DoesGlobalFieldExist() )
776  {
777  fFieldPropagator->ClearPropagatorState();
778  // Resets all state of field propagator class (ONLY)
779  // including safety values (in case of overlaps and to wipe for first track).
780 
781  // G4ChordFinder* chordF= fFieldPropagator->GetChordFinder();
782  // if( chordF ) chordF->ResetStepEstimate();
783  }
784 
785  // Make sure to clear the chord finders of all fields (ie managers)
786  //
788  fieldMgrStore->ClearAllChordFindersState();
789 
790  // Update the current touchable handle (from the track's)
791  //
792  fCurrentTouchableHandle = aTrack->GetTouchableHandle();
793 }