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