Geant4  10.01
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 87375 2014-12-02 08:17:28Z 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 "G4ITSafetyHelper.hh"
69 #include "G4FieldManagerStore.hh"
70 
71 #include "G4UnitsTable.hh"
72 #include "G4ReferenceCast.hh"
73 
75 
76 #ifndef State
77 #define State(theXInfo) (GetState<G4ITTransportationState>()->theXInfo)
78 #endif
79 
80 //#define DEBUG_MEM
81 
82 #ifdef DEBUG_MEM
83 #include "G4MemStat.hh"
84 using namespace G4MemStat;
85 using G4MemStat::MemStat;
86 #endif
87 
88 //#define G4DEBUG_TRANSPORT 1
89 
92  fThreshold_Warning_Energy(100 * MeV),
93  fThreshold_Important_Energy(250 * MeV),
94  fThresholdTrials(10),
95  fUnimportant_Energy(1 * MeV), // Not used
96  fSumEnergyKilled(0.0),
97  fMaxEnergyKilled(0.0),
98  fShortStepOptimisation(false), // Old default: true (=fast short steps)
99  fVerboseLevel(verbose)
100 {
102  G4TransportationManager* transportMgr;
104  G4ITTransportationManager* ITtransportMgr;
106  fLinearNavigator = ITtransportMgr->GetNavigatorForTracking();
107  fFieldPropagator = transportMgr->GetPropagatorInField();
108  fpSafetyHelper = ITtransportMgr->GetSafetyHelper(); // New
109 
110  // Cannot determine whether a field exists here, as it would
111  // depend on the relative order of creating the detector's
112  // field and this process. That order is not guaranted.
113  // Instead later the method DoesGlobalFieldExist() is called
114 
115  enableAtRestDoIt = false;
116  enableAlongStepDoIt = true;
117  enablePostStepDoIt = true;
118  SetProcessSubType(60);
122 
124  /*
125  if(fTransportationState == 0)
126  {
127  G4cout << "KILL in G4ITTransportation::G4ITTransportation" << G4endl;
128  abort();
129  }
130  */
131 }
132 
134  G4VITProcess(right)
135 {
136  // Copy attributes
145 
146  // Setup Navigators
147  G4TransportationManager* transportMgr;
149  G4ITTransportationManager* ITtransportMgr;
151  fLinearNavigator = ITtransportMgr->GetNavigatorForTracking();
152  fFieldPropagator = transportMgr->GetPropagatorInField();
153  fpSafetyHelper = ITtransportMgr->GetSafetyHelper(); // New
154 
155  // Cannot determine whether a field exists here, as it would
156  // depend on the relative order of creating the detector's
157  // field and this process. That order is not guaranted.
158  // Instead later the method DoesGlobalFieldExist() is called
159 
160  enableAtRestDoIt = false;
161  enableAlongStepDoIt = true;
162  enablePostStepDoIt = true;
163 
168 }
169 
171 {
172  if (this == &right) return *this;
173  return *this;
174 }
175 
180  G4ProcessState(), fCurrentTouchableHandle(0)
181 {
185  fTransportEndSpin = G4ThreeVector(0, 0, 0);
186  fMomentumChanged = false;
187  fEnergyChanged = false;
188  fEndGlobalTimeComputed = false;
190  fParticleIsLooping = false;
191 
192  static G4ThreadLocal G4TouchableHandle *nullTouchableHandle = 0;
193  if (!nullTouchableHandle) nullTouchableHandle = new G4TouchableHandle;
194  // Points to (G4VTouchable*) 0
195 
196  fCurrentTouchableHandle = *nullTouchableHandle;
197  fGeometryLimitedStep = false;
199  fPreviousSafety = 0.0;
200  fNoLooperTrials = false;
201  fEndPointDistance = -1;
202 }
203 
205 {
206  ;
207 }
208 
210 {
211 #ifdef G4VERBOSE
212  if ((fVerboseLevel > 0) && (fSumEnergyKilled > 0.0))
213  {
214  G4cout << " G4ITTransportation: Statistics for looping particles "
215  << G4endl;
216  G4cout << " Sum of energy of loopers killed: "
217  << fSumEnergyKilled << G4endl;
218  G4cout << " Max energy of loopers killed: "
219  << fMaxEnergyKilled << G4endl;
220  }
221 #endif
222 }
223 
225 {
227 }
228 
230 {
231  G4TransportationManager* transportMgr;
233 
234  // fFieldExists= transportMgr->GetFieldManager()->DoesFieldExist();
235  // return fFieldExists;
236  return transportMgr->GetFieldManager()->DoesFieldExist();
237 }
238 
240 //
241 // Responsibilities:
242 // Find whether the geometry limits the Step, and to what length
243 // Calculate the new value of the safety and return it.
244 // Store the final time, position and momentum.
245 G4double
248  G4double,
249  G4double currentMinimumStep,
250  G4double& currentSafety,
251  G4GPILSelection* selection)
252 {
253  G4double geometryStepLength(-1.0), newSafety(-1.0);
254 
255  State(fParticleIsLooping) = false;
256  State(fEndGlobalTimeComputed) = false;
257  State(fGeometryLimitedStep) = false;
258 
259  // Initial actions moved to StartTrack()
260  // --------------------------------------
261  // Note: in case another process changes touchable handle
262  // it will be necessary to add here (for all steps)
263  // State(fCurrentTouchableHandle) = track.GetTouchableHandle();
264 
265  // GPILSelection is set to defaule value of CandidateForSelection
266  // It is a return value
267  //
268  *selection = CandidateForSelection;
269 
270  // Get initial Energy/Momentum of the track
271  //
272  const G4DynamicParticle* pParticle = track.GetDynamicParticle();
273 // const G4ParticleDefinition* pParticleDef = pParticle->GetDefinition() ;
274  G4ThreeVector startMomentumDir = pParticle->GetMomentumDirection();
275  G4ThreeVector startPosition = track.GetPosition();
276 
277  // G4double theTime = track.GetGlobalTime() ;
278 
279  // The Step Point safety can be limited by other geometries and/or the
280  // assumptions of any process - it's not always the geometrical safety.
281  // We calculate the starting point's isotropic safety here.
282  //
283  G4ThreeVector OriginShift = startPosition - State(fPreviousSftOrigin);
284  G4double MagSqShift = OriginShift.mag2();
285  if (MagSqShift >= sqr(State(fPreviousSafety)))
286  {
287  currentSafety = 0.0;
288  }
289  else
290  {
291  currentSafety = State(fPreviousSafety) - std::sqrt(MagSqShift);
292  }
293 
294  // Is the particle charged ?
295  //
296  G4double particleCharge = pParticle->GetCharge();
297 
298  // There is no need to locate the current volume. It is Done elsewhere:
299  // On track construction
300  // By the tracking, after all AlongStepDoIts, in "Relocation"
301 
302  // Check whether the particle have an (EM) field force exerting upon it
303  //
304  G4FieldManager* fieldMgr = 0;
305  G4bool fieldExertsForce = false;
306  if ((particleCharge != 0.0))
307  {
309  if (fieldMgr != 0)
310  {
311  // Message the field Manager, to configure it for this track
312  fieldMgr->ConfigureForTrack(&track);
313  // Moved here, in order to allow a transition
314  // from a zero-field status (with fieldMgr->(field)0
315  // to a finite field status
316 
317  // If the field manager has no field, there is no field !
318  fieldExertsForce = (fieldMgr->GetDetectorField() != 0);
319  }
320  }
321 
322  // G4cout << " G4Transport: field exerts force= " << fieldExertsForce
323  // << " fieldMgr= " << fieldMgr << G4endl;
324 
325  // Choose the calculation of the transportation: Field or not
326  //
327  if (!fieldExertsForce)
328  {
329  G4double linearStepLength;
330  if (fShortStepOptimisation && (currentMinimumStep <= currentSafety))
331  {
332  // The Step is guaranteed to be taken
333  //
334  geometryStepLength = currentMinimumStep;
335  State(fGeometryLimitedStep) = false;
336  }
337  else
338  {
339  // Find whether the straight path intersects a volume
340  //
341  // fLinearNavigator->SetNavigatorState(GetIT(track)->GetTrackingInfo()->GetNavigatorState());
342  linearStepLength = fLinearNavigator->ComputeStep(startPosition,
343  startMomentumDir,
344  currentMinimumStep,
345  newSafety);
346 
347 // G4cout << "linearStepLength : " << G4BestUnit(linearStepLength,"Length")
348 // << " | currentMinimumStep: " << currentMinimumStep
349 // << " | trackID: " << track.GetTrackID() << G4endl;
350 
351  // Remember last safety origin & value.
352  //
353  State(fPreviousSftOrigin) = startPosition;
354  State(fPreviousSafety) = newSafety;
355 
356  G4TrackStateManager& trackStateMan = GetIT(track)->GetTrackingInfo()
358  fpSafetyHelper->LoadTrackState(trackStateMan);
359  // fpSafetyHelper->SetTrackState(state);
360  fpSafetyHelper->SetCurrentSafety(newSafety,
361  State(fTransportEndPosition));
363 
364  // The safety at the initial point has been re-calculated:
365  //
366  currentSafety = newSafety;
367 
368  State(fGeometryLimitedStep) = (linearStepLength <= currentMinimumStep);
369  if (State(fGeometryLimitedStep))
370  {
371  // The geometry limits the Step size (an intersection was found.)
372  geometryStepLength = linearStepLength;
373  }
374  else
375  {
376  // The full Step is taken.
377  geometryStepLength = currentMinimumStep;
378  }
379  }
380  State(fEndPointDistance) = geometryStepLength;
381 
382  // Calculate final position
383  //
384  State(fTransportEndPosition) = startPosition
385  + geometryStepLength * startMomentumDir;
386 
387  // Momentum direction, energy and polarisation are unchanged by transport
388  //
389  State(fTransportEndMomentumDir) = startMomentumDir;
390  State(fTransportEndKineticEnergy) = track.GetKineticEnergy();
391  State(fTransportEndSpin) = track.GetPolarization();
392  State(fParticleIsLooping) = false;
393  State(fMomentumChanged) = false;
394  State(fEndGlobalTimeComputed) = true;
395  State(theInteractionTimeLeft) = State(fEndPointDistance)
396  / track.GetVelocity();
397  State(fCandidateEndGlobalTime) = State(theInteractionTimeLeft)
398  + track.GetGlobalTime();
399 /*
400  G4cout << "track.GetVelocity() : "
401  << track.GetVelocity() << G4endl;
402  G4cout << "State(endpointDistance) : "
403  << G4BestUnit(State(endpointDistance),"Length") << G4endl;
404  G4cout << "State(theInteractionTimeLeft) : "
405  << G4BestUnit(State(theInteractionTimeLeft),"Time") << G4endl;
406  G4cout << "track.GetGlobalTime() : "
407  << G4BestUnit(track.GetGlobalTime(),"Time") << G4endl;
408 */
409  }
410  else // A field exerts force
411  {
412 
413  G4ExceptionDescription exceptionDescription;
414  exceptionDescription
415  << "ITTransportation does not support external fields.";
416  exceptionDescription
417  << " If you are dealing with a tradiational MC simulation, ";
418  exceptionDescription << "please use G4Transportation.";
419 
420  G4Exception("G4ITTransportation::AlongStepGetPhysicalInteractionLength",
421  "NoExternalFieldSupport", FatalException, exceptionDescription);
422  /*
423  G4double momentumMagnitude = pParticle->GetTotalMomentum() ;
424  // G4ThreeVector EndUnitMomentum ;
425  G4double lengthAlongCurve ;
426  G4double restMass = pParticleDef->GetPDGMass() ;
427 
428  fFieldPropagator->SetChargeMomentumMass( particleCharge, // in e+ units
429  momentumMagnitude, // in Mev/c
430  restMass ) ;
431 
432  G4ThreeVector spin = track.GetPolarization() ;
433  G4FieldTrack aFieldTrack = G4FieldTrack( startPosition,
434  track.GetMomentumDirection(),
435  0.0,
436  track.GetKineticEnergy(),
437  restMass,
438  track.GetVelocity(),
439  track.GetGlobalTime(), // Lab.
440  track.GetProperTime(), // Part.
441  &spin ) ;
442  if( currentMinimumStep > 0 )
443  {
444  // Do the Transport in the field (non recti-linear)
445  //
446  lengthAlongCurve = fFieldPropagator->ComputeStep( aFieldTrack,
447  currentMinimumStep,
448  currentSafety,
449  track.GetVolume() ) ;
450  State(fGeometryLimitedStep)= lengthAlongCurve < currentMinimumStep;
451  if( State(fGeometryLimitedStep) )
452  {
453  geometryStepLength = lengthAlongCurve ;
454  }
455  else
456  {
457  geometryStepLength = currentMinimumStep ;
458  }
459 
460  // Remember last safety origin & value.
461  //
462  State(fPreviousSftOrigin) = startPosition ;
463  State(fPreviousSafety) = currentSafety ;
464  fpSafetyHelper->SetCurrentSafety( newSafety, startPosition);
465  }
466  else
467  {
468  geometryStepLength = lengthAlongCurve= 0.0 ;
469  State(fGeometryLimitedStep) = false ;
470  }
471 
472  // Get the End-Position and End-Momentum (Dir-ection)
473  //
474  State(fTransportEndPosition) = aFieldTrack.GetPosition() ;
475 
476  // Momentum: Magnitude and direction can be changed too now ...
477  //
478  State(fMomentumChanged) = true ;
479  State(fTransportEndMomentumDir) = aFieldTrack.GetMomentumDir() ;
480 
481  State(fTransportEndKineticEnergy) = aFieldTrack.GetKineticEnergy() ;
482 
483  if( fFieldPropagator->GetCurrentFieldManager()->DoesFieldChangeEnergy() )
484  {
485  // If the field can change energy, then the time must be integrated
486  // - so this should have been updated
487  //
488  State(fCandidateEndGlobalTime) = aFieldTrack.GetLabTimeOfFlight();
489  State(fEndGlobalTimeComputed) = true;
490 
491  State(theInteractionTimeLeft) = State(fCandidateEndGlobalTime) -
492  track.GetGlobalTime() ;
493 
494  // was ( State(fCandidateEndGlobalTime) != track.GetGlobalTime() );
495  // a cleaner way is to have FieldTrack knowing whether time is updated.
496  }
497  else
498  {
499  // The energy should be unchanged by field transport,
500  // - so the time changed will be calculated elsewhere
501  //
502  State(fEndGlobalTimeComputed) = false;
503 
504  // Check that the integration preserved the energy
505  // - and if not correct this!
506  G4double startEnergy= track.GetKineticEnergy();
507  G4double endEnergy= State(fTransportEndKineticEnergy);
508 
509  static G4int no_inexact_steps=0, no_large_ediff;
510  G4double absEdiff = std::fabs(startEnergy- endEnergy);
511  if( absEdiff > perMillion * endEnergy )
512  {
513  no_inexact_steps++;
514  // Possible statistics keeping here ...
515  }
516  #ifdef G4VERBOSE
517  if( fVerboseLevel > 1 )
518  {
519  if( std::fabs(startEnergy- endEnergy) > perThousand * endEnergy )
520  {
521  static G4int no_warnings= 0, warnModulo=1, moduloFactor= 10;
522  no_large_ediff ++;
523  if( (no_large_ediff% warnModulo) == 0 )
524  {
525  no_warnings++;
526  G4cout << "WARNING - G4Transportation::AlongStepGetPIL() "
527  << " Energy change in Step is above 1^-3 relative value. " << G4endl
528  << " Relative change in 'tracking' step = "
529  << std::setw(15) << (endEnergy-startEnergy)/startEnergy << G4endl
530  << " Starting E= " << std::setw(12) << startEnergy / MeV << " MeV "
531  << G4endl
532  << " Ending E= " << std::setw(12) << endEnergy / MeV << " MeV "
533  << G4endl;
534  G4cout << " Energy has been corrected -- however, review"
535  << " field propagation parameters for accuracy." << G4endl;
536  if( (fVerboseLevel > 2 ) || (no_warnings<4) ||
537  (no_large_ediff == warnModulo * moduloFactor) )
538  {
539  G4cout << " These include EpsilonStepMax(/Min) in G4FieldManager "
540  << " which determine fractional error per step for integrated quantities. "
541  << G4endl
542  << " Note also the influence of the permitted number of integration steps."
543  << G4endl;
544  }
545  G4cerr << "ERROR - G4Transportation::AlongStepGetPIL()" << G4endl
546  << " Bad 'endpoint'. Energy change detected"
547  << " and corrected. "
548  << " Has occurred already "
549  << no_large_ediff << " times." << G4endl;
550  if( no_large_ediff == warnModulo * moduloFactor )
551  {
552  warnModulo *= moduloFactor;
553  }
554  }
555  }
556  } // end of if (fVerboseLevel)
557  #endif
558  // Correct the energy for fields that conserve it
559  // This - hides the integration error
560  // - but gives a better physical answer
561  State(fTransportEndKineticEnergy)= track.GetKineticEnergy();
562  }
563 
564  State(fTransportEndSpin) = aFieldTrack.GetSpin();
565  State(fParticleIsLooping) = fFieldPropagator->IsParticleLooping() ;
566  State(endpointDistance) = (State(fTransportEndPosition) -
567  startPosition).mag() ;
568  // State(theInteractionTimeLeft) =
569  track.GetVelocity()/State(endpointDistance) ;
570  */
571  }
572 
573  // If we are asked to go a step length of 0, and we are on a boundary
574  // then a boundary will also limit the step -> we must flag this.
575  //
576  if (currentMinimumStep == 0.0)
577  {
578  if (currentSafety == 0.0)
579  {
580  State(fGeometryLimitedStep) = true;
581 // G4cout << "!!!! Safety is NULL, on the Boundary !!!!!" << G4endl;
582 // G4cout << " Track position : " << track.GetPosition() /nanometer
583 // << G4endl;
584  }
585  }
586 
587  // Update the safety starting from the end-point,
588  // if it will become negative at the end-point.
589  //
590  if (currentSafety < State(fEndPointDistance))
591  {
592  // if( particleCharge == 0.0 )
593  // G4cout << " Avoiding call to ComputeSafety : charge = 0.0 " << G4endl;
594 
595  if (particleCharge != 0.0)
596  {
597 
598  G4double endSafety = fLinearNavigator->ComputeSafety(
599  State(fTransportEndPosition));
600  currentSafety = endSafety;
601  State(fPreviousSftOrigin) = State(fTransportEndPosition);
602  State(fPreviousSafety) = currentSafety;
603 
604  /*
605  G4VTrackStateHandle state =
606  GetIT(track)->GetTrackingInfo()->GetTrackState(fpSafetyHelper);
607  */
608  G4TrackStateManager& trackStateMan = GetIT(track)->GetTrackingInfo()
610  fpSafetyHelper->LoadTrackState(trackStateMan);
611  // fpSafetyHelper->SetTrackState(state);
612  fpSafetyHelper->SetCurrentSafety(currentSafety,
613  State(fTransportEndPosition));
615 
616  // Because the Stepping Manager assumes it is from the start point,
617  // add the StepLength
618  //
619  currentSafety += State(fEndPointDistance);
620 
621 #ifdef G4DEBUG_TRANSPORT
622  G4cout.precision(12);
623  G4cout << "***G4Transportation::AlongStepGPIL ** " << G4endl;
624  G4cout << " Called Navigator->ComputeSafety at "
625  << State(fTransportEndPosition)
626  << " and it returned safety= " << endSafety << G4endl;
627  G4cout << " Adding endpoint distance " << State(fEndPointDistance)
628  << " to obtain pseudo-safety= " << currentSafety << G4endl;
629 #endif
630  }
631  }
632 
633  // fParticleChange.ProposeTrueStepLength(geometryStepLength) ;
634 
635 // G4cout << "G4ITTransportation::AlongStepGetPhysicalInteractionLength = "
636 // << G4BestUnit(geometryStepLength,"Length") << G4endl;
637 
638  return geometryStepLength;
639 }
640 
642  const G4Step& /*step*/,
643  const double timeStep,
644  double& oPhysicalStep)
645 {
646  const G4DynamicParticle* pParticle = track.GetDynamicParticle();
647  G4ThreeVector startMomentumDir = pParticle->GetMomentumDirection();
648  G4ThreeVector startPosition = track.GetPosition();
649 
650  track.CalculateVelocity();
651  G4double initialVelocity = track.GetVelocity();
652 
653  State(fGeometryLimitedStep) = false;
654 
656  // !!! CASE NO FIELD !!!
658  State(fCandidateEndGlobalTime) = timeStep + track.GetGlobalTime();
659  State(fEndGlobalTimeComputed) = true;
660 
661  // Choose the calculation of the transportation: Field or not
662  //
663  if (!State(fMomentumChanged))
664  {
665  // G4cout << "Momentum has not changed" << G4endl;
666  fParticleChange.ProposeVelocity(initialVelocity);
667  oPhysicalStep = initialVelocity * timeStep;
668 
669  // Calculate final position
670  //
671  State(fTransportEndPosition) = startPosition
672  + oPhysicalStep * startMomentumDir;
673  }
674 }
675 
677 //
678 // Initialize ParticleChange (by setting all its members equal
679 // to corresponding members in G4Track)
680 #include "G4ParticleTable.hh"
682  const G4Step& stepData)
683 {
684 
685 #if defined (DEBUG_MEM)
686  MemStat mem_first, mem_second, mem_diff;
687 #endif
688 
689 #if defined (DEBUG_MEM)
690  mem_first = MemoryUsage();
691 #endif
692 
693  // G4cout << "G4ITTransportation::AlongStepDoIt" << G4endl;
694  // set pdefOpticalPhoton
695  // Andrea Dotti: the following statement should be in a single line:
696  // G4-MT transformation tools get confused if statement spans two lines
697  // If needed contact: adotti@slac.stanford.edu
698  static G4ThreadLocal G4ParticleDefinition* pdefOpticalPhoton = 0;
699  if (!pdefOpticalPhoton) pdefOpticalPhoton =
701 
702  static G4ThreadLocal G4int noCalls = 0;
703  noCalls++;
704 
706 
707  // Code for specific process
708  //
709  fParticleChange.ProposePosition(State(fTransportEndPosition));
710  fParticleChange.ProposeMomentumDirection(State(fTransportEndMomentumDir));
711  fParticleChange.ProposeEnergy(State(fTransportEndKineticEnergy));
712  fParticleChange.SetMomentumChanged(State(fMomentumChanged));
713 
714  fParticleChange.ProposePolarization(State(fTransportEndSpin));
715 
716  G4double deltaTime = 0.0;
717 
718  // Calculate Lab Time of Flight (ONLY if field Equations used it!)
719  // G4double endTime = State(fCandidateEndGlobalTime);
720  // G4double delta_time = endTime - startTime;
721 
722  G4double startTime = track.GetGlobalTime();
726  if (State(fEndGlobalTimeComputed) == false)
727  {
728  // The time was not integrated .. make the best estimate possible
729  //
730  G4double initialVelocity = stepData.GetPreStepPoint()->GetVelocity();
731  G4double stepLength = track.GetStepLength();
732 
733  deltaTime = 0.0; // in case initialVelocity = 0
734  if (track.GetParticleDefinition() == pdefOpticalPhoton)
735  {
736  // For only Optical Photon, final velocity is used
737  double finalVelocity = track.CalculateVelocityForOpticalPhoton();
738  fParticleChange.ProposeVelocity(finalVelocity);
739  deltaTime = stepLength / finalVelocity;
740  }
741  else if (initialVelocity > 0.0)
742  {
743  deltaTime = stepLength / initialVelocity;
744  }
745 
746  State(fCandidateEndGlobalTime) = startTime + deltaTime;
747  }
748  else
749  {
750  deltaTime = State(fCandidateEndGlobalTime) - startTime;
751  }
752 
753  fParticleChange.ProposeGlobalTime(State(fCandidateEndGlobalTime));
754  fParticleChange.ProposeLocalTime(track.GetLocalTime() + deltaTime);
755  /*
756  // Now Correct by Lorentz factor to get delta "proper" Time
757 
758  G4double restMass = track.GetDynamicParticle()->GetMass() ;
759  G4double deltaProperTime = deltaTime*( restMass/track.GetTotalEnergy() ) ;
760 
761  fParticleChange.ProposeProperTime(track.GetProperTime() + deltaProperTime) ;
762  */
763 
765 
768 
769  // If the particle is caught looping or is stuck (in very difficult
770  // boundaries) in a magnetic field (doing many steps)
771  // THEN this kills it ...
772  //
773  if (State(fParticleIsLooping))
774  {
775  G4double endEnergy = State(fTransportEndKineticEnergy);
776 
777  if ((endEnergy < fThreshold_Important_Energy) || (State(fNoLooperTrials)
778  >= fThresholdTrials))
779  {
780  // Kill the looping particle
781  //
782  // G4cout << "G4ITTransportation will killed the molecule"<< G4endl;
784 
785  // 'Bare' statistics
786  fSumEnergyKilled += endEnergy;
787  if (endEnergy > fMaxEnergyKilled)
788  {
789  fMaxEnergyKilled = endEnergy;
790  }
791 
792 #ifdef G4VERBOSE
793  if ((fVerboseLevel > 1) || (endEnergy > fThreshold_Warning_Energy))
794  {
795  G4cout
796  << " G4ITTransportation is killing track that is looping or stuck "
797  << G4endl<< " This track has " << track.GetKineticEnergy() / MeV
798  << " MeV energy." << G4endl;
799  G4cout << " Number of trials = " << State(fNoLooperTrials)
800  << " No of calls to AlongStepDoIt = " << noCalls
801  << G4endl;
802  }
803 #endif
804  State(fNoLooperTrials) = 0;
805  }
806  else
807  {
808  State(fNoLooperTrials)++;
809 #ifdef G4VERBOSE
810  if ((fVerboseLevel > 2))
811  {
812  G4cout << " G4ITTransportation::AlongStepDoIt(): Particle looping - "
813  << " Number of trials = " << State(fNoLooperTrials)
814  << " No of calls to = " << noCalls << G4endl;
815  }
816 #endif
817  }
818  }
819  else
820  {
821  State(fNoLooperTrials)=0;
822  }
823 
824  // Another (sometimes better way) is to use a user-limit maximum Step size
825  // to alleviate this problem ..
826 
827  // Introduce smooth curved trajectories to particle-change
828  //
831 
832 #if defined (DEBUG_MEM)
833  mem_second = MemoryUsage();
834  mem_diff = mem_second-mem_first;
835  G4cout << "\t || MEM || End of G4ITTransportation::AlongStepDoIt, diff is: "
836  << mem_diff << G4endl;
837 #endif
838 
839  return &fParticleChange;
840 }
841 
843 //
844 // This ensures that the PostStep action is always called,
845 // so that it can do the relocation if it is needed.
846 //
847 
848 G4double
851  G4double, // previousStepSize
852  G4ForceCondition* pForceCond)
853 {
854  *pForceCond = Forced;
855  return DBL_MAX; // was kInfinity ; but convention now is DBL_MAX
856 }
857 
859 //
860 
862  const G4Step&)
863 {
864  // G4cout << "G4ITTransportation::PostStepDoIt" << G4endl;
865  G4TouchableHandle retCurrentTouchable; // The one to return
866  G4bool isLastStep = false;
867 
868  // Initialize ParticleChange (by setting all its members equal
869  // to corresponding members in G4Track)
870  fParticleChange.Initialize(track); // To initialise TouchableChange
871 
873 
874  // If the Step was determined by the volume boundary,
875  // logically relocate the particle
876 
877  if (State(fGeometryLimitedStep))
878  {
879 
880  if(fVerboseLevel)
881  {
882  G4cout << "Step is limited by geometry "
883  << "track ID : " << track.GetTrackID() << G4endl;
884  }
885 
886  // fCurrentTouchable will now become the previous touchable,
887  // and what was the previous will be freed.
888  // (Needed because the preStepPoint can point to the previous touchable)
889 
890  if ( State(fCurrentTouchableHandle)->GetVolume() == 0)
891  {
892  G4ExceptionDescription exceptionDescription;
893  exceptionDescription << "No current touchable found ";
894  G4Exception(" G4ITTransportation::PostStepDoIt", "G4ITTransportation001",
895  FatalErrorInArgument, exceptionDescription);
896  }
897 
898  fLinearNavigator->SetGeometricallyLimitedStep();
899  fLinearNavigator->LocateGlobalPointAndUpdateTouchableHandle(
900  track.GetPosition(), track.GetMomentumDirection(),
901  State(fCurrentTouchableHandle), true);
902  // Check whether the particle is out of the world volume
903  // If so it has exited and must be killed.
904  //
905  if ( State(fCurrentTouchableHandle)->GetVolume() == 0)
906  {
907  // abort();
908 #ifdef G4VERBOSE
909  if (fVerboseLevel > 0)
910  {
911  G4cout << "Track position : " << track.GetPosition() / nanometer
912  << " [nm]" << " Track ID : " << track.GetTrackID() << G4endl;
913  G4cout << "G4ITTransportation will killed the track because "
914  "State(fCurrentTouchableHandle)->GetVolume() == 0"<< G4endl;
915  }
916 #endif
918  }
919 
920  retCurrentTouchable = State(fCurrentTouchableHandle);
921 
922 // G4cout << "Current volume : " << track.GetVolume()->GetName()
923 // << " Next volume : "
924 // << (State(fCurrentTouchableHandle)->GetVolume() ?
925 // State(fCurrentTouchableHandle)->GetVolume()->GetName():"OutWorld")
926 // << " Position : " << track.GetPosition() / nanometer
927 // << " track ID : " << track.GetTrackID()
928 // << G4endl;
929 
930  fParticleChange.SetTouchableHandle(State(fCurrentTouchableHandle));
931 
932  // Update the Step flag which identifies the Last Step in a volume
933  isLastStep = fLinearNavigator->ExitedMotherVolume()
934  | fLinearNavigator->EnteredDaughterVolume();
935 
936 #ifdef G4DEBUG_TRANSPORT
937  // Checking first implementation of flagging Last Step in Volume
938  G4bool exiting = fLinearNavigator->ExitedMotherVolume();
939  G4bool entering = fLinearNavigator->EnteredDaughterVolume();
940 
941  if( ! (exiting || entering) )
942  {
943  G4cout << " Transport> : Proposed isLastStep= " << isLastStep
944  << " Exiting " << fLinearNavigator->ExitedMotherVolume()
945  << " Entering " << fLinearNavigator->EnteredDaughterVolume()
946  << " Track position : " << track.GetPosition() /nanometer << " [nm]"
947  << G4endl;
948  G4cout << " Track position : " << track.GetPosition() /nanometer
949  << G4endl;
950  }
951 #endif
952  }
953  else // fGeometryLimitedStep is false
954  {
955  // This serves only to move the Navigator's location
956  //
957 // abort();
958  fLinearNavigator->LocateGlobalPointWithinVolume(track.GetPosition());
959 
960  // The value of the track's current Touchable is retained.
961  // (and it must be correct because we must use it below to
962  // overwrite the (unset) one in particle change)
963  // It must be fCurrentTouchable too ??
964  //
966  retCurrentTouchable = track.GetTouchableHandle();
967 
968  isLastStep = false;
969 #ifdef G4DEBUG_TRANSPORT
970  // Checking first implementation of flagging Last Step in Volume
971  G4cout << " Transport> Proposed isLastStep= " << isLastStep
972  << " Geometry did not limit step. Position : "
973  << track.GetPosition()/ nanometer << G4endl;
974 #endif
975  } // endif ( fGeometryLimitedStep )
976 
978 
979  const G4VPhysicalVolume* pNewVol = retCurrentTouchable->GetVolume();
980  const G4Material* pNewMaterial = 0;
981  const G4VSensitiveDetector* pNewSensitiveDetector = 0;
982 
983  if (pNewVol != 0)
984  {
985  pNewMaterial = pNewVol->GetLogicalVolume()->GetMaterial();
986  pNewSensitiveDetector = pNewVol->GetLogicalVolume()->GetSensitiveDetector();
987  }
988 
989  // ( <const_cast> pNewMaterial ) ;
990  // ( <const_cast> pNewSensitiveDetector) ;
991 
994  (G4VSensitiveDetector *) pNewSensitiveDetector);
995 
996  const G4MaterialCutsCouple* pNewMaterialCutsCouple = 0;
997  if (pNewVol != 0)
998  {
999  pNewMaterialCutsCouple =
1001  }
1002 
1003  if (pNewVol != 0 && pNewMaterialCutsCouple != 0
1004  && pNewMaterialCutsCouple->GetMaterial() != pNewMaterial)
1005  {
1006  // for parametrized volume
1007  //
1008  pNewMaterialCutsCouple = G4ProductionCutsTable::GetProductionCutsTable()
1009  ->GetMaterialCutsCouple(pNewMaterial,
1010  pNewMaterialCutsCouple->GetProductionCuts());
1011  }
1012  fParticleChange.SetMaterialCutsCoupleInTouchable(pNewMaterialCutsCouple);
1013 
1014  // temporarily until Get/Set Material of ParticleChange,
1015  // and StepPoint can be made const.
1016  // Set the touchable in ParticleChange
1017  // this must always be done because the particle change always
1018  // uses this value to overwrite the current touchable pointer.
1019  //
1020  fParticleChange.SetTouchableHandle(retCurrentTouchable);
1021 
1022  return &fParticleChange;
1023 }
1024 
1025 // New method takes over the responsibility to reset the state of
1026 // G4Transportation object at the start of a new track or the resumption of
1027 // a suspended track.
1028 
1030 {
1033  {
1034 // G4VITProcess::fpState = new G4ITTransportationState();
1036  // Will set in the same time fTransportationState
1037  }
1038 
1041  GetIT(track)->GetTrackingInfo()->GetTrackStateManager());
1042 
1043  // The actions here are those that were taken in AlongStepGPIL
1044  // when track.GetCurrentStepNumber()==1
1045 
1046  // reset safety value and center
1047  //
1048  // State(fPreviousSafety) = 0.0 ;
1049  // State(fPreviousSftOrigin) = G4ThreeVector(0.,0.,0.) ;
1050 
1051  // reset looping counter -- for motion in field
1052  // State(fNoLooperTrials)= 0;
1053  // Must clear this state .. else it depends on last track's value
1054  // --> a better solution would set this from state of suspended track TODO ?
1055  // Was if( aTrack->GetCurrentStepNumber()==1 ) { .. }
1056 
1057  // ChordFinder reset internal state
1058  //
1059  if (DoesGlobalFieldExist())
1060  {
1062  // Resets all state of field propagator class (ONLY)
1063  // including safety values (in case of overlaps and to wipe for first track).
1064 
1065  // G4ChordFinder* chordF= fFieldPropagator->GetChordFinder();
1066  // if( chordF ) chordF->ResetStepEstimate();
1067  }
1068 
1069  // Make sure to clear the chord finders of all fields (ie managers)
1070  static G4ThreadLocal G4FieldManagerStore* fieldMgrStore = 0;
1071  if (!fieldMgrStore) fieldMgrStore = G4FieldManagerStore::GetInstance();
1072  fieldMgrStore->ClearAllChordFindersState();
1073 
1074  // Update the current touchable handle (from the track's)
1075  //
1076  State(fCurrentTouchableHandle) = track->GetTouchableHandle();
1077 
1079 }
1080 
1081 #undef State
void SetMaterialInTouchable(G4Material *fMaterial)
G4ITSafetyHelper * GetSafetyHelper() const
void SetTouchableHandle(const G4TouchableHandle &fTouchable)
const G4ThreeVector & GetPolarization() const
static const double MeV
Definition: G4SIunits.hh:193
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
G4double GetLocalTime() const
G4ParticleChangeForTransport fParticleChange
void SetMaterialCutsCoupleInTouchable(const G4MaterialCutsCouple *fMaterialCutsCouple)
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
std::vector< G4ThreeVector > * GimmeTrajectoryVectorAndForgetIt() const
CLHEP::Hep3Vector G4ThreeVector
G4double GetVelocity() const
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
const G4DynamicParticle * GetDynamicParticle() const
G4Material * GetMaterial() const
static const double nanometer
Definition: G4SIunits.hh:91
const G4ThreeVector & GetPosition() const
G4TrackStatus GetTrackStatus() const
virtual void LoadTrackState(G4TrackStateManager &manager)
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &)
void ProposePolarization(G4double Px, G4double Py, G4double Pz)
virtual G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData)
G4double GetVelocity() const
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
#define G4ThreadLocal
Definition: tls.hh:84
G4bool enableAtRestDoIt
Definition: G4VProcess.hh:350
#define State(theXInfo)
int G4int
Definition: G4Types.hh:78
MemStat MemoryUsage()
Definition: G4MemStat.cc:55
void SetInstantiateProcessState(G4bool flag)
void ProposePosition(G4double x, G4double y, G4double z)
G4ReferenceCountedHandle< G4VTouchable > G4TouchableHandle
G4ITNavigator * fLinearNavigator
G4StepPoint * GetPreStepPoint() const
virtual void Initialize(const G4Track &)
virtual void StartTracking(G4Track *)
Definition: G4VProcess.cc:101
G4double CalculateVelocityForOpticalPhoton() const
Definition: G4Track.cc:247
virtual void ConfigureForTrack(const G4Track *)
G4ITTransportation & operator=(const G4ITTransportation &)
virtual void StartTracking(G4Track *aTrack)
G4IT * GetIT(const G4Track *track)
Definition: G4IT.cc:49
G4double GetKineticEnergy() const
static G4ITTransportationManager * GetTransportationManager()
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &, G4double, G4ForceCondition *pForceCond)
G4GLOB_DLL std::ostream G4cout
virtual void ResetTrackState()
virtual void NewTrackState()
virtual void ComputeStep(const G4Track &, const G4Step &, const double timeStep, double &spaceStep)
bool G4bool
Definition: G4Types.hh:79
const G4ThreeVector & GetMomentumDirection() const
G4PropagatorInField * fFieldPropagator
G4double GetCharge() const
void ProposeTrueStepLength(G4double truePathLength)
G4bool DoesFieldExist() const
virtual void StartTracking(G4Track *)
Definition: G4VITProcess.cc:84
void SetInstantiateProcessState(G4bool flag)
const G4ParticleDefinition * GetParticleDefinition() const
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
Definition: G4Step.hh:76
G4bool enablePostStepDoIt
Definition: G4VProcess.hh:352
G4int GetTrackID() const
G4double GetGlobalTime() const
G4double CalculateVelocity() const
Definition: G4Track.cc:214
G4FieldManager * FindAndSetFieldManager(G4VPhysicalVolume *pCurrentPhysVol)
const G4TouchableHandle & GetTouchableHandle() const
G4ITSafetyHelper * fpSafetyHelper
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4TransportationManager * GetTransportationManager()
G4FieldManager * GetFieldManager() const
G4TrackingInformation * GetTrackingInfo()
Definition: G4IT.hh:152
static G4ProductionCutsTable * GetProductionCutsTable()
const G4ThreeVector & GetMomentumDirection() const
G4LogicalVolume * GetLogicalVolume() const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
static G4FieldManagerStore * GetInstance()
void SetSensitiveDetectorInTouchable(G4VSensitiveDetector *fSensitiveDetector)
static G4ParticleTable * GetParticleTable()
#define fPreviousSftOrigin
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:44
void SetPointerToVectorOfAuxiliaryPoints(std::vector< G4ThreeVector > *theNewVectorPointer)
void ProposeGlobalTime(G4double t)
G4ITTransportation(const G4String &aName="ITTransportation", G4int verbosityLevel=0)
#define fPreviousSafety
G4ITNavigator * GetNavigatorForTracking() const
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
void ProposeEnergy(G4double finalEnergy)
G4::shared_ptr< G4ProcessState > fpState
G4VPhysicalVolume * GetVolume() const
#define G4endl
Definition: G4ios.hh:61
G4bool enableAlongStepDoIt
Definition: G4VProcess.hh:351
void ProposeLastStepInVolume(G4bool flag)
virtual void SaveTrackState(G4TrackStateManager &manager)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
T sqr(const T &x)
Definition: templates.hh:145
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
G4double fThreshold_Important_Energy
const G4Field * GetDetectorField() const
G4TrackStateManager & GetTrackStateManager()
G4ForceCondition
G4VITProcess inherits from G4VProcess.
Definition: G4VITProcess.hh:99
void SetMomentumChanged(G4bool b)
G4PropagatorInField * GetPropagatorInField() const
void ProposeVelocity(G4double finalVelocity)
G4ProductionCuts * GetProductionCuts() const
#define DBL_MAX
Definition: templates.hh:83
G4VSensitiveDetector * GetSensitiveDetector() const
G4double * theInteractionTimeLeft
void ProposeLocalTime(G4double t)
void SetCurrentSafety(G4double val, const G4ThreeVector &pos)
G4double GetStepLength() const
G4GPILSelection
const G4Material * GetMaterial() const