Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4ITTransportation.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 // $Id: G4ITTransportation.cc 64374 2012-10-31 16:37:23Z gcosmo $
27 //
31 //
32 // Contact : Mathieu Karamitros (kara (AT) cenbg . in2p3 . fr)
33 //
34 // History :
35 // -----------
36 // =======================================================================
37 // Modified:
38 // 28 Oct 2011, P.Gumpl./J.Ap: Detect gravity field, use magnetic moment
39 // 20 Nov 2008, J.Apostolakis: Push safety to helper - after ComputeSafety
40 // 9 Nov 2007, J.Apostolakis: Flag for short steps, push safety to helper
41 // 19 Jan 2006, P.MoraDeFreitas: Fix for suspended tracks (StartTracking)
42 // 11 Aug 2004, M.Asai: Add G4VSensitiveDetector* for updating stepPoint.
43 // 21 June 2003, J.Apostolakis: Calling field manager with
44 // track, to enable it to configure its accuracy
45 // 13 May 2003, J.Apostolakis: Zero field areas now taken into
46 // account correclty in all cases (thanks to W Pokorski).
47 // 29 June 2001, J.Apostolakis, D.Cote-Ahern, P.Gumplinger:
48 // correction for spin tracking
49 // 20 Febr 2001, J.Apostolakis: update for new FieldTrack
50 // 22 Sept 2000, V.Grichine: update of Kinetic Energy
51 // ---------------------------------------------------
52 // 10 Oct 2011, M.Karamitros: G4ITTransportation created
53 // Created: 19 March 1997, J. Apostolakis
54 // =======================================================================
55 //
56 // -------------------------------------------------------------------
57 
58 #include "G4ITTransportation.hh"
59 #include "G4SystemOfUnits.hh"
62 #include "G4ProductionCutsTable.hh"
63 #include "G4ParticleTable.hh"
64 #include "G4ITNavigator.hh"
65 #include "G4PropagatorInField.hh"
66 #include "G4FieldManager.hh"
67 #include "G4ChordFinder.hh"
68 #include "G4SafetyHelper.hh"
69 #include "G4FieldManagerStore.hh"
70 
71 #include "G4UnitsTable.hh"
72 
74 
75 #ifndef State
76 #define State(theXInfo) (fTransportationState->theXInfo)
77 #endif
78 
79 //#define G4DEBUG_TRANSPORT 1
80 
83  InitProcessState(fTransportationState, fpState),
84  fThreshold_Warning_Energy( 100 * MeV ),
85  fThreshold_Important_Energy( 250 * MeV ),
86  fThresholdTrials( 10 ),
87  fUnimportant_Energy( 1 * MeV ), // Not used
88  fSumEnergyKilled( 0.0 ), fMaxEnergyKilled( 0.0 ),
89  fShortStepOptimisation(false), // Old default: true (=fast short steps)
90  fVerboseLevel( verbose )
91 {
93  G4TransportationManager* transportMgr ;
95  G4ITTransportationManager* ITtransportMgr ;
97  fLinearNavigator = ITtransportMgr->GetNavigatorForTracking() ;
98  fFieldPropagator = transportMgr->GetPropagatorInField() ;
99  fpSafetyHelper = 0; // transportMgr->GetSafetyHelper(); // New
100 
101  // Cannot determine whether a field exists here, as it would
102  // depend on the relative order of creating the detector's
103  // field and this process. That order is not guaranted.
104  // Instead later the method DoesGlobalFieldExist() is called
105 
106  enableAtRestDoIt = false;
107  enableAlongStepDoIt = true;
108  enablePostStepDoIt = true;
109  SetProcessSubType(60);
112  fInstantiateProcessState = true;
113 }
114 
115 
117  G4VITProcess(right),
118  InitProcessState(fTransportationState, fpState)
119 {
120  // Copy attributes
121  fVerboseLevel = right.fVerboseLevel ;
129 
130  // Setup Navigators
131  G4TransportationManager* transportMgr ;
133  G4ITTransportationManager* ITtransportMgr ;
135  fLinearNavigator = ITtransportMgr->GetNavigatorForTracking() ;
136  fFieldPropagator = transportMgr->GetPropagatorInField() ;
137  fpSafetyHelper = 0; //transportMgr->GetSafetyHelper(); // New
138 
139  // Cannot determine whether a field exists here, as it would
140  // depend on the relative order of creating the detector's
141  // field and this process. That order is not guaranted.
142  // Instead later the method DoesGlobalFieldExist() is called
143 
144  enableAtRestDoIt = false;
145  enableAlongStepDoIt = true;
146  enablePostStepDoIt = true;
147 
151  fInstantiateProcessState = right.fInstantiateProcessState;
152 }
153 
154 G4ITTransportation& G4ITTransportation::operator=(const G4ITTransportation& right)
155 {
156  if(this == &right) return *this;
157  return *this;
158 }
159 
164  fCurrentTouchableHandle(0)
165 {
170  fMomentumChanged = false;
171  fEnergyChanged = false;
172  fEndGlobalTimeComputed = false;
174  fParticleIsLooping = false;
175  static G4TouchableHandle nullTouchableHandle; // Points to (G4VTouchable*) 0
176  fCurrentTouchableHandle = nullTouchableHandle;
177  fGeometryLimitedStep = false;
179  fPreviousSafety = 0.0;
180  fNoLooperTrials = false;
181  endpointDistance= -1;
182 }
183 
185 {
186  ;
187 }
188 
190 {
191 #ifdef G4VERBOSE
192  if( (fVerboseLevel > 0) && (fSumEnergyKilled > 0.0 ) )
193  {
194  G4cout << " G4ITTransportation: Statistics for looping particles " << G4endl;
195  G4cout << " Sum of energy of loopers killed: " << fSumEnergyKilled << G4endl;
196  G4cout << " Max energy of loopers killed: " << fMaxEnergyKilled << G4endl;
197  }
198 #endif
199 }
200 
202 {
203  G4TransportationManager* transportMgr;
205 
206  // fFieldExists= transportMgr->GetFieldManager()->DoesFieldExist();
207  // return fFieldExists;
208  return transportMgr->GetFieldManager()->DoesFieldExist();
209 }
210 
211 
213 //
214 // Responsibilities:
215 // Find whether the geometry limits the Step, and to what length
216 // Calculate the new value of the safety and return it.
217 // Store the final time, position and momentum.
219  const G4Track& track,
220  G4double ,
221  G4double currentMinimumStep,
222  G4double& currentSafety,
223  G4GPILSelection* selection)
224 {
225  G4double geometryStepLength(-1.0), newSafety(-1.0) ;
226 
227  State(fParticleIsLooping) = false ;
228  State(fEndGlobalTimeComputed) = false ;
229  State(fGeometryLimitedStep) = false ;
230 
231  // Initial actions moved to StartTrack()
232  // --------------------------------------
233  // Note: in case another process changes touchable handle
234  // it will be necessary to add here (for all steps)
235  State(fCurrentTouchableHandle) = track.GetTouchableHandle();
236 
237  // GPILSelection is set to defaule value of CandidateForSelection
238  // It is a return value
239  //
240  *selection = CandidateForSelection ;
241 
242  // Get initial Energy/Momentum of the track
243  //
244  const G4DynamicParticle* pParticle = track.GetDynamicParticle() ;
245 // const G4ParticleDefinition* pParticleDef = pParticle->GetDefinition() ;
246  G4ThreeVector startMomentumDir = pParticle->GetMomentumDirection() ;
247  G4ThreeVector startPosition = track.GetPosition() ;
248 
249  // G4double theTime = track.GetGlobalTime() ;
250 
251  // The Step Point safety can be limited by other geometries and/or the
252  // assumptions of any process - it's not always the geometrical safety.
253  // We calculate the starting point's isotropic safety here.
254  //
255  G4ThreeVector OriginShift = startPosition - State(fPreviousSftOrigin) ;
256  G4double MagSqShift = OriginShift.mag2() ;
257  if( MagSqShift >= sqr(State(fPreviousSafety)) )
258  {
259  currentSafety = 0.0 ;
260  }
261  else
262  {
263  currentSafety = State(fPreviousSafety) - std::sqrt(MagSqShift) ;
264  }
265 
266  // Is the particle charged ?
267  //
268  G4double particleCharge = pParticle->GetCharge() ;
269 
270 
271  // There is no need to locate the current volume. It is Done elsewhere:
272  // On track construction
273  // By the tracking, after all AlongStepDoIts, in "Relocation"
274 
275  // Check whether the particle have an (EM) field force exerting upon it
276  //
277  G4FieldManager* fieldMgr=0;
278  G4bool fieldExertsForce = false ;
279  if( (particleCharge != 0.0) )
280  {
281  fieldMgr= fFieldPropagator->FindAndSetFieldManager( track.GetVolume() );
282  if (fieldMgr != 0)
283  {
284  // Message the field Manager, to configure it for this track
285  fieldMgr->ConfigureForTrack( &track );
286  // Moved here, in order to allow a transition
287  // from a zero-field status (with fieldMgr->(field)0
288  // to a finite field status
289 
290  // If the field manager has no field, there is no field !
291  fieldExertsForce = (fieldMgr->GetDetectorField() != 0);
292  }
293  }
294 
295  // G4cout << " G4Transport: field exerts force= " << fieldExertsForce
296  // << " fieldMgr= " << fieldMgr << G4endl;
297 
298  // Choose the calculation of the transportation: Field or not
299  //
300  if( !fieldExertsForce )
301  {
302  G4double linearStepLength ;
303  if( fShortStepOptimisation && (currentMinimumStep <= currentSafety) )
304  {
305  // The Step is guaranteed to be taken
306  //
307  geometryStepLength = currentMinimumStep ;
308  State(fGeometryLimitedStep) = false ;
309  }
310  else
311  {
312  // Find whether the straight path intersects a volume
313  //
314  linearStepLength = fLinearNavigator->ComputeStep( startPosition,
315  startMomentumDir,
316  currentMinimumStep,
317  newSafety) ;
318  // Remember last safety origin & value.
319  //
320  State(fPreviousSftOrigin) = startPosition ;
321  State(fPreviousSafety) = newSafety ;
322  // fpSafetyHelper->SetCurrentSafety( newSafety, startPosition);
323 
324  // The safety at the initial point has been re-calculated:
325  //
326  currentSafety = newSafety ;
327 
328  State(fGeometryLimitedStep)= (linearStepLength <= currentMinimumStep);
329  if( State(fGeometryLimitedStep) )
330  {
331  // The geometry limits the Step size (an intersection was found.)
332  geometryStepLength = linearStepLength ;
333  }
334  else
335  {
336  // The full Step is taken.
337  geometryStepLength = currentMinimumStep ;
338  }
339  }
340  State(endpointDistance) = geometryStepLength ;
341 
342  // Calculate final position
343  //
344  State(fTransportEndPosition) = startPosition+geometryStepLength*startMomentumDir ;
345 
346  // Momentum direction, energy and polarisation are unchanged by transport
347  //
348  State(fTransportEndMomentumDir) = startMomentumDir ;
349  State(fTransportEndKineticEnergy) = track.GetKineticEnergy() ;
350  State(fTransportEndSpin) = track.GetPolarization();
351  State(fParticleIsLooping) = false ;
352  State(fMomentumChanged) = false ;
353  State(fEndGlobalTimeComputed) = true ;
354  State(theInteractionTimeLeft) = State(endpointDistance)/track.GetVelocity();
355  State(fCandidateEndGlobalTime) = State(theInteractionTimeLeft)+track.GetGlobalTime();
356 
357  // G4cout << "track.GetVelocity() : " << track.GetVelocity() << G4endl;
358  // G4cout << "State(endpointDistance) : " << G4BestUnit(State(endpointDistance),"Length") << G4endl;
359  // G4cout << "State(theInteractionTimeLeft) : " << G4BestUnit(State(theInteractionTimeLeft),"Time") << G4endl;
360  // G4cout << "track.GetGlobalTime() : " << G4BestUnit(track.GetGlobalTime(),"Time") << G4endl;
361 
362  }
363  else // A field exerts force
364  {
365 
366  G4ExceptionDescription exceptionDescription;
367  exceptionDescription << "ITTransportation does not support external fields.";
368  exceptionDescription << " If you are dealing with a tradiational MC simulation, ";
369  exceptionDescription << "please use G4Transportation.";
370 
371  G4Exception("G4ITTransportation::AlongStepGetPhysicalInteractionLength","NoExternalFieldSupport",
372  FatalException,exceptionDescription);
373  /*
374  G4double momentumMagnitude = pParticle->GetTotalMomentum() ;
375  // G4ThreeVector EndUnitMomentum ;
376  G4double lengthAlongCurve ;
377  G4double restMass = pParticleDef->GetPDGMass() ;
378 
379  fFieldPropagator->SetChargeMomentumMass( particleCharge, // in e+ units
380  momentumMagnitude, // in Mev/c
381  restMass ) ;
382 
383  G4ThreeVector spin = track.GetPolarization() ;
384  G4FieldTrack aFieldTrack = G4FieldTrack( startPosition,
385  track.GetMomentumDirection(),
386  0.0,
387  track.GetKineticEnergy(),
388  restMass,
389  track.GetVelocity(),
390  track.GetGlobalTime(), // Lab.
391  track.GetProperTime(), // Part.
392  &spin ) ;
393  if( currentMinimumStep > 0 )
394  {
395  // Do the Transport in the field (non recti-linear)
396  //
397  lengthAlongCurve = fFieldPropagator->ComputeStep( aFieldTrack,
398  currentMinimumStep,
399  currentSafety,
400  track.GetVolume() ) ;
401  State(fGeometryLimitedStep)= lengthAlongCurve < currentMinimumStep;
402  if( State(fGeometryLimitedStep) )
403  {
404  geometryStepLength = lengthAlongCurve ;
405  }
406  else
407  {
408  geometryStepLength = currentMinimumStep ;
409  }
410 
411  // Remember last safety origin & value.
412  //
413  State(fPreviousSftOrigin) = startPosition ;
414  State(fPreviousSafety) = currentSafety ;
415  fpSafetyHelper->SetCurrentSafety( newSafety, startPosition);
416  }
417  else
418  {
419  geometryStepLength = lengthAlongCurve= 0.0 ;
420  State(fGeometryLimitedStep) = false ;
421  }
422 
423  // Get the End-Position and End-Momentum (Dir-ection)
424  //
425  State(fTransportEndPosition) = aFieldTrack.GetPosition() ;
426 
427  // Momentum: Magnitude and direction can be changed too now ...
428  //
429  State(fMomentumChanged) = true ;
430  State(fTransportEndMomentumDir) = aFieldTrack.GetMomentumDir() ;
431 
432  State(fTransportEndKineticEnergy) = aFieldTrack.GetKineticEnergy() ;
433 
434  if( fFieldPropagator->GetCurrentFieldManager()->DoesFieldChangeEnergy() )
435  {
436  // If the field can change energy, then the time must be integrated
437  // - so this should have been updated
438  //
439  State(fCandidateEndGlobalTime) = aFieldTrack.GetLabTimeOfFlight();
440  State(fEndGlobalTimeComputed) = true;
441 
442  State(theInteractionTimeLeft) = State(fCandidateEndGlobalTime) - track.GetGlobalTime() ;
443 
444  // was ( State(fCandidateEndGlobalTime) != track.GetGlobalTime() );
445  // a cleaner way is to have FieldTrack knowing whether time is updated.
446  }
447  else
448  {
449  // The energy should be unchanged by field transport,
450  // - so the time changed will be calculated elsewhere
451  //
452  State(fEndGlobalTimeComputed) = false;
453 
454  // Check that the integration preserved the energy
455  // - and if not correct this!
456  G4double startEnergy= track.GetKineticEnergy();
457  G4double endEnergy= State(fTransportEndKineticEnergy);
458 
459  static G4int no_inexact_steps=0, no_large_ediff;
460  G4double absEdiff = std::fabs(startEnergy- endEnergy);
461  if( absEdiff > perMillion * endEnergy )
462  {
463  no_inexact_steps++;
464  // Possible statistics keeping here ...
465  }
466 #ifdef G4VERBOSE
467  if( fVerboseLevel > 1 )
468  {
469  if( std::fabs(startEnergy- endEnergy) > perThousand * endEnergy )
470  {
471  static G4int no_warnings= 0, warnModulo=1, moduloFactor= 10;
472  no_large_ediff ++;
473  if( (no_large_ediff% warnModulo) == 0 )
474  {
475  no_warnings++;
476  G4cout << "WARNING - G4Transportation::AlongStepGetPIL() "
477  << " Energy change in Step is above 1^-3 relative value. " << G4endl
478  << " Relative change in 'tracking' step = "
479  << std::setw(15) << (endEnergy-startEnergy)/startEnergy << G4endl
480  << " Starting E= " << std::setw(12) << startEnergy / MeV << " MeV " << G4endl
481  << " Ending E= " << std::setw(12) << endEnergy / MeV << " MeV " << G4endl;
482  G4cout << " Energy has been corrected -- however, review"
483  << " field propagation parameters for accuracy." << G4endl;
484  if( (fVerboseLevel > 2 ) || (no_warnings<4) || (no_large_ediff == warnModulo * moduloFactor) )
485  {
486  G4cout << " These include EpsilonStepMax(/Min) in G4FieldManager "
487  << " which determine fractional error per step for integrated quantities. " << G4endl
488  << " Note also the influence of the permitted number of integration steps."
489  << G4endl;
490  }
491  G4cerr << "ERROR - G4Transportation::AlongStepGetPIL()" << G4endl
492  << " Bad 'endpoint'. Energy change detected"
493  << " and corrected. "
494  << " Has occurred already "
495  << no_large_ediff << " times." << G4endl;
496  if( no_large_ediff == warnModulo * moduloFactor )
497  {
498  warnModulo *= moduloFactor;
499  }
500  }
501  }
502  } // end of if (fVerboseLevel)
503 #endif
504  // Correct the energy for fields that conserve it
505  // This - hides the integration error
506  // - but gives a better physical answer
507  State(fTransportEndKineticEnergy)= track.GetKineticEnergy();
508  }
509 
510  State(fTransportEndSpin) = aFieldTrack.GetSpin();
511  State(fParticleIsLooping) = fFieldPropagator->IsParticleLooping() ;
512  State(endpointDistance) = (State(fTransportEndPosition) - startPosition).mag() ;
513  // State(theInteractionTimeLeft) = track.GetVelocity()/State(endpointDistance) ;
514 */
515  }
516 
517  // If we are asked to go a step length of 0, and we are on a boundary
518  // then a boundary will also limit the step -> we must flag this.
519  //
520  if( currentMinimumStep == 0.0 )
521  {
522  if( currentSafety == 0.0 )
523  {
524  State(fGeometryLimitedStep) = true ;
525 // G4cout << "!!!! Safety is NULL, on the Boundary !!!!!" << G4endl;
526 // G4cout << " Track position : " << track.GetPosition() /nanometer << G4endl;
527  }
528  }
529 
530  // Update the safety starting from the end-point,
531  // if it will become negative at the end-point.
532  //
533  if( currentSafety < State(endpointDistance) )
534  {
535  // if( particleCharge == 0.0 )
536  // G4cout << " Avoiding call to ComputeSafety : charge = 0.0 " << G4endl;
537 
538  if( particleCharge != 0.0 )
539  {
540 
541  G4double endSafety =
542  fLinearNavigator->ComputeSafety( State(fTransportEndPosition)) ;
543  currentSafety = endSafety ;
544  State(fPreviousSftOrigin) = State(fTransportEndPosition) ;
545  State(fPreviousSafety) = currentSafety ;
546  fpSafetyHelper->SetCurrentSafety( currentSafety, State(fTransportEndPosition));
547 
548  // Because the Stepping Manager assumes it is from the start point,
549  // add the StepLength
550  //
551  currentSafety += State(endpointDistance) ;
552 
553 #ifdef G4DEBUG_TRANSPORT
554  G4cout.precision(12) ;
555  G4cout << "***G4Transportation::AlongStepGPIL ** " << G4endl ;
556  G4cout << " Called Navigator->ComputeSafety at " << State(fTransportEndPosition)
557  << " and it returned safety= " << endSafety << G4endl ;
558  G4cout << " Adding endpoint distance " << State(endpointDistance)
559  << " to obtain pseudo-safety= " << currentSafety << G4endl ;
560 #endif
561  }
562  }
563 
564  // fParticleChange.ProposeTrueStepLength(geometryStepLength) ;
565 
566  return geometryStepLength ;
567 }
568 
569 
571  const G4Step& /*step*/,
572  const double timeStep,
573  double& oPhysicalStep)
574 {
575  const G4DynamicParticle* pParticle = track.GetDynamicParticle() ;
576  G4ThreeVector startMomentumDir = pParticle->GetMomentumDirection() ;
577  G4ThreeVector startPosition = track.GetPosition() ;
578 
579  track.CalculateVelocity();
580  G4double initialVelocity = track.GetVelocity() ;
581 
582  State(fGeometryLimitedStep) = false;
583 
585  // !!! CASE NO FIELD !!!
587  State(fCandidateEndGlobalTime) = timeStep + track.GetGlobalTime();
588  State(fEndGlobalTimeComputed) = true ;
589 
590  // Choose the calculation of the transportation: Field or not
591  //
592  if( !State(fMomentumChanged) )
593  {
594  // G4cout << "Momentum has not changed" << G4endl;
595  fParticleChange.ProposeVelocity(initialVelocity);
596  oPhysicalStep = initialVelocity*timeStep ;
597 
598  // Calculate final position
599  //
600  State(fTransportEndPosition) = startPosition + oPhysicalStep*startMomentumDir ;
601  }
602 }
603 
604 
606 //
607 // Initialize ParticleChange (by setting all its members equal
608 // to corresponding members in G4Track)
609 #include "G4ParticleTable.hh"
611  const G4Step& stepData )
612 {
613 
614  // G4cout << "G4ITTransportation::AlongStepDoIt" << G4endl;
615  // set pdefOpticalPhoton
616  static G4ParticleDefinition* pdefOpticalPhoton =
618 
619  static G4int noCalls=0;
620  noCalls++;
621 
622  fParticleChange.Initialize(track) ;
623 
624  // Code for specific process
625  //
626  fParticleChange.ProposePosition(State(fTransportEndPosition)) ;
627  fParticleChange.ProposeMomentumDirection(State(fTransportEndMomentumDir)) ;
628  fParticleChange.ProposeEnergy(State(fTransportEndKineticEnergy)) ;
629  fParticleChange.SetMomentumChanged(State(fMomentumChanged)) ;
630 
631  fParticleChange.ProposePolarization(State(fTransportEndSpin));
632 
633  G4double deltaTime = 0.0 ;
634 
635  // Calculate Lab Time of Flight (ONLY if field Equations used it!)
636  // G4double endTime = State(fCandidateEndGlobalTime);
637  // G4double delta_time = endTime - startTime;
638 
639  G4double startTime = track.GetGlobalTime() ;
643  if (State(fEndGlobalTimeComputed) == false)
644  {
645  // The time was not integrated .. make the best estimate possible
646  //
647  G4double initialVelocity = stepData.GetPreStepPoint()->GetVelocity() ;
648  G4double stepLength = track.GetStepLength() ;
649 
650  deltaTime= 0.0; // in case initialVelocity = 0
651  if (track.GetParticleDefinition() == pdefOpticalPhoton)
652  {
653  // For only Optical Photon, final velocity is used
654  double finalVelocity = track.CalculateVelocityForOpticalPhoton();
655  fParticleChange.ProposeVelocity(finalVelocity);
656  deltaTime = stepLength/finalVelocity ;
657  }
658  else if( initialVelocity > 0.0 )
659  {
660  deltaTime = stepLength/initialVelocity ;
661  }
662 
663  State(fCandidateEndGlobalTime) = startTime + deltaTime ;
664  }
665  else
666  {
667  deltaTime = State(fCandidateEndGlobalTime) - startTime ;
668  }
669 
670  fParticleChange.ProposeGlobalTime( State(fCandidateEndGlobalTime) ) ;
671  fParticleChange.ProposeLocalTime( track.GetLocalTime() + deltaTime) ;
672  /*
673  // Now Correct by Lorentz factor to get delta "proper" Time
674 
675  G4double restMass = track.GetDynamicParticle()->GetMass() ;
676  G4double deltaProperTime = deltaTime*( restMass/track.GetTotalEnergy() ) ;
677 
678  fParticleChange.ProposeProperTime(track.GetProperTime() + deltaProperTime) ;
679 */
680 
681  fParticleChange. ProposeTrueStepLength( track.GetStepLength() ) ;
682 
685 
686  // If the particle is caught looping or is stuck (in very difficult
687  // boundaries) in a magnetic field (doing many steps)
688  // THEN this kills it ...
689  //
690  if ( State(fParticleIsLooping) )
691  {
692  G4double endEnergy= State(fTransportEndKineticEnergy);
693 
694  if( (endEnergy < fThreshold_Important_Energy)
695  || (State(fNoLooperTrials) >= fThresholdTrials ) )
696  {
697  // Kill the looping particle
698  //
699  // G4cout << "G4ITTransportation will killed the molecule"<< G4endl;
701 
702  // 'Bare' statistics
703  fSumEnergyKilled += endEnergy;
704  if( endEnergy > fMaxEnergyKilled)
705  {
706  fMaxEnergyKilled= endEnergy;
707  }
708 
709 #ifdef G4VERBOSE
710  if( (fVerboseLevel > 1) ||
711  ( endEnergy > fThreshold_Warning_Energy ) )
712  {
713  G4cout << " G4ITTransportation is killing track that is looping or stuck "
714  << G4endl
715  << " This track has " << track.GetKineticEnergy() / MeV
716  << " MeV energy." << G4endl;
717  G4cout << " Number of trials = " << State(fNoLooperTrials)
718  << " No of calls to AlongStepDoIt = " << noCalls
719  << G4endl;
720  }
721 #endif
722  State(fNoLooperTrials)=0;
723  }
724  else
725  {
726  State(fNoLooperTrials) ++;
727 #ifdef G4VERBOSE
728  if( (fVerboseLevel > 2) )
729  {
730  G4cout << " G4ITTransportation::AlongStepDoIt(): Particle looping - "
731  << " Number of trials = " << State(fNoLooperTrials)
732  << " No of calls to = " << noCalls
733  << G4endl;
734  }
735 #endif
736  }
737  }
738  else
739  {
740  State(fNoLooperTrials)=0;
741  }
742 
743  // Another (sometimes better way) is to use a user-limit maximum Step size
744  // to alleviate this problem ..
745 
746  // Introduce smooth curved trajectories to particle-change
747  //
750 
751  return &fParticleChange ;
752 }
753 
755 //
756 // This ensures that the PostStep action is always called,
757 // so that it can do the relocation if it is needed.
758 //
759 
762  const G4Track& , // track
763  G4double , // previousStepSize
764  G4ForceCondition* pForceCond
765  )
766 {
767  *pForceCond = Forced ;
768  return DBL_MAX ; // was kInfinity ; but convention now is DBL_MAX
769 }
770 
772 //
773 
775  const G4Step& )
776 {
777  // G4cout << "G4ITTransportation::PostStepDoIt" << G4endl;
778  G4TouchableHandle retCurrentTouchable ; // The one to return
779  G4bool isLastStep= false;
780 
781  // Initialize ParticleChange (by setting all its members equal
782  // to corresponding members in G4Track)
783  fParticleChange.Initialize(track) ; // To initialise TouchableChange
784 
786 
787  // If the Step was determined by the volume boundary,
788  // logically relocate the particle
789 
790  if(State(fGeometryLimitedStep))
791  {
792 
793  // G4cout << "Step is limited by geometry " << "track ID : " << track.GetTrackID() << G4endl;
794 
795  // fCurrentTouchable will now become the previous touchable,
796  // and what was the previous will be freed.
797  // (Needed because the preStepPoint can point to the previous touchable)
798 
799  if( State(fCurrentTouchableHandle)->GetVolume() == 0 )
800  {
801  G4ExceptionDescription exceptionDescription ;
802  exceptionDescription << "No current touchable found " ;
803  G4Exception(" G4ITTransportation::PostStepDoIt","G4ITTransportation001",
804  FatalErrorInArgument,exceptionDescription);
805  }
806 
809  LocateGlobalPointAndUpdateTouchableHandle( track.GetPosition(),
810  track.GetMomentumDirection(),
811  State(fCurrentTouchableHandle),
812  true ) ;
813  // Check whether the particle is out of the world volume
814  // If so it has exited and must be killed.
815  //
816  if( State(fCurrentTouchableHandle)->GetVolume() == 0 )
817  {
818 #ifdef G4VERBOSE
819  if(fVerboseLevel > 0)
820  {
821  G4cout << "Track position : " << track.GetPosition() / nanometer << " [nm]"
822  << " Track ID : " << track.GetTrackID()<< G4endl;
823  G4cout << "G4ITTransportation will killed the track because State(fCurrentTouchableHandle)->GetVolume() == 0"<< G4endl;
824  }
825 #endif
827  }
828 
829  retCurrentTouchable = State(fCurrentTouchableHandle) ;
830 
831 // G4cout << "Current volume : " << track.GetVolume()->GetName()
832 // << " Next volume : "
833 // << (State(fCurrentTouchableHandle)->GetVolume() ? State(fCurrentTouchableHandle)->GetVolume()->GetName():"OutWorld")
834 // << " Position : " << track.GetPosition() / nanometer
835 // << " track ID : " << track.GetTrackID()
836 // << G4endl;
837 
838  fParticleChange.SetTouchableHandle( State(fCurrentTouchableHandle) ) ;
839 
840  // Update the Step flag which identifies the Last Step in a volume
841  isLastStep = fLinearNavigator->ExitedMotherVolume()
843 
844 #ifdef G4DEBUG_TRANSPORT
845  // Checking first implementation of flagging Last Step in Volume
848 
849  if( ! (exiting || entering) )
850  {
851  G4cout << " Transport> : Proposed isLastStep= " << isLastStep
852  << " Exiting " << fLinearNavigator->ExitedMotherVolume()
853  << " Entering " << fLinearNavigator->EnteredDaughterVolume()
854  << " Track position : " << track.GetPosition() /nanometer << " [nm]"
855  << G4endl;
856  G4cout << " Track position : " << track.GetPosition() /nanometer << G4endl;
857  }
858 #endif
859  }
860  else // fGeometryLimitedStep is false
861  {
862  // This serves only to move the Navigator's location
863  //
865 
866  // The value of the track's current Touchable is retained.
867  // (and it must be correct because we must use it below to
868  // overwrite the (unset) one in particle change)
869  // It must be fCurrentTouchable too ??
870  //
872  retCurrentTouchable = track.GetTouchableHandle() ;
873 
874  isLastStep= false;
875 #ifdef G4DEBUG_TRANSPORT
876  // Checking first implementation of flagging Last Step in Volume
877  G4cout << " Transport> Proposed isLastStep= " << isLastStep
878  << " Geometry did not limit step. Position : "
879  << track.GetPosition()/ nanometer << G4endl;
880 #endif
881  } // endif ( fGeometryLimitedStep )
882 
884 
885  const G4VPhysicalVolume* pNewVol = retCurrentTouchable->GetVolume() ;
886  const G4Material* pNewMaterial = 0 ;
887  const G4VSensitiveDetector* pNewSensitiveDetector = 0 ;
888 
889  if( pNewVol != 0 )
890  {
891  pNewMaterial= pNewVol->GetLogicalVolume()->GetMaterial();
892  pNewSensitiveDetector= pNewVol->GetLogicalVolume()->GetSensitiveDetector();
893  }
894 
895  // ( <const_cast> pNewMaterial ) ;
896  // ( <const_cast> pNewSensitiveDetector) ;
897 
900 
901  const G4MaterialCutsCouple* pNewMaterialCutsCouple = 0;
902  if( pNewVol != 0 )
903  {
904  pNewMaterialCutsCouple=pNewVol->GetLogicalVolume()->GetMaterialCutsCouple();
905  }
906 
907  if( pNewVol!=0 && pNewMaterialCutsCouple!=0 && pNewMaterialCutsCouple->GetMaterial()!=pNewMaterial )
908  {
909  // for parametrized volume
910  //
911  pNewMaterialCutsCouple =
913  ->GetMaterialCutsCouple(pNewMaterial,
914  pNewMaterialCutsCouple->GetProductionCuts());
915  }
916  fParticleChange.SetMaterialCutsCoupleInTouchable( pNewMaterialCutsCouple );
917 
918  // temporarily until Get/Set Material of ParticleChange,
919  // and StepPoint can be made const.
920  // Set the touchable in ParticleChange
921  // this must always be done because the particle change always
922  // uses this value to overwrite the current touchable pointer.
923  //
924  fParticleChange.SetTouchableHandle(retCurrentTouchable) ;
925 
926  return &fParticleChange ;
927 }
928 
929 // New method takes over the responsibility to reset the state of G4Transportation
930 // object at the start of a new track or the resumption of a suspended track.
931 
932 void
934 {
936  if(fInstantiateProcessState)
938  // Will set in the same time fTransportationState
939 
940  // The actions here are those that were taken in AlongStepGPIL
941  // when track.GetCurrentStepNumber()==1
942 
943  // reset safety value and center
944  //
945  // State(fPreviousSafety) = 0.0 ;
946  // State(fPreviousSftOrigin) = G4ThreeVector(0.,0.,0.) ;
947 
948  // reset looping counter -- for motion in field
949  // State(fNoLooperTrials)= 0;
950  // Must clear this state .. else it depends on last track's value
951  // --> a better solution would set this from state of suspended track TODO ?
952  // Was if( aTrack->GetCurrentStepNumber()==1 ) { .. }
953 
954  // ChordFinder reset internal state
955  //
956  if( DoesGlobalFieldExist() )
957  {
959  // Resets all state of field propagator class (ONLY)
960  // including safety values (in case of overlaps and to wipe for first track).
961 
962  // G4ChordFinder* chordF= fFieldPropagator->GetChordFinder();
963  // if( chordF ) chordF->ResetStepEstimate();
964  }
965 
966  // Make sure to clear the chord finders of all fields (ie managers)
967  static G4FieldManagerStore* fieldMgrStore= G4FieldManagerStore::GetInstance();
968  fieldMgrStore->ClearAllChordFindersState();
969 
970  // Update the current touchable handle (from the track's)
971  //
972  State(fCurrentTouchableHandle) = track->GetTouchableHandle();
973 
975 }
976 
977 #undef State