Geant4  10.01
G4CoupledTransportation.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 //
27 // $Id: G4CoupledTransportation.cc 86964 2014-11-21 11:47:44Z gcosmo $
28 //
29 // ------------------------------------------------------------
30 // GEANT 4 class implementation
31 // =======================================================================
32 // Modified:
33 // 13 May 2006, J. Apostolakis: Revised for parallel navigation (PathFinder)
34 // 19 Jan 2006, P.MoraDeFreitas: Fix for suspended tracks (StartTracking)
35 // 11 Aug 2004, M.Asai: Add G4VSensitiveDetector* for updating stepPoint.
36 // 21 June 2003, J.Apostolakis: Calling field manager with
37 // track, to enable it to configure its accuracy
38 // 13 May 2003, J.Apostolakis: Zero field areas now taken into
39 // account correclty in all cases (thanks to W Pokorski).
40 // 29 June 2001, J.Apostolakis, D.Cote-Ahern, P.Gumplinger:
41 // correction for spin tracking
42 // 20 Febr 2001, J.Apostolakis: update for new FieldTrack
43 // 22 Sept 2000, V.Grichine: update of Kinetic Energy
44 // Created: 19 March 1997, J. Apostolakis
45 // =======================================================================
46 
48 
49 #include "G4PhysicalConstants.hh"
50 #include "G4SystemOfUnits.hh"
52 #include "G4ProductionCutsTable.hh"
53 #include "G4ParticleTable.hh"
54 #include "G4ChordFinder.hh"
55 #include "G4Field.hh"
56 #include "G4FieldTrack.hh"
57 #include "G4FieldManagerStore.hh"
58 
60 
63 //
64 // Constructor
65 
67  : G4VProcess( G4String("CoupledTransportation"), fTransportation ),
68  fTransportEndPosition(0.0, 0.0, 0.0),
69  fTransportEndMomentumDir(0.0, 0.0, 0.0),
70  fTransportEndKineticEnergy(0.0),
71  fTransportEndSpin(0.0, 0.0, 0.0), // fTransportEndPolarization(0.0, 0.0, 0.0),
72  fMomentumChanged(false),
73  fEndGlobalTimeComputed(false),
74  fCandidateEndGlobalTime(0.0),
75  fParticleIsLooping( false ),
76  fPreviousSftOrigin( 0.,0.,0. ),
77  fPreviousMassSafety( 0.0 ),
78  fPreviousFullSafety( 0.0 ),
79  fMassGeometryLimitedStep( false ),
80  fAnyGeometryLimitedStep( false ),
81  fEndpointDistance( -1.0 ),
82  fThreshold_Warning_Energy( 100 * MeV ),
83  fThreshold_Important_Energy( 250 * MeV ),
84  fThresholdTrials( 10 ),
85  fNoLooperTrials( 0 ),
86  fSumEnergyKilled( 0.0 ), fMaxEnergyKilled( 0.0 ),
87  fVerboseLevel( verbosity )
88 {
89  // set Process Sub Type
90  SetProcessSubType(static_cast<G4int>(COUPLED_TRANSPORTATION));
91 
92  G4TransportationManager* transportMgr ;
93 
95 
96  fMassNavigator = transportMgr->GetNavigatorForTracking() ;
97  fFieldPropagator = transportMgr->GetPropagatorInField() ;
98  // fGlobalFieldMgr = transportMgr->GetFieldManager() ;
100  if( fVerboseLevel > 0 )
101  {
102  G4cout << " G4CoupledTransportation constructor: ----- " << G4endl;
103  G4cout << " Verbose level is " << fVerboseLevel << G4endl;
104  G4cout << " Navigator Id obtained in G4CoupledTransportation constructor "
105  << fNavigatorId << G4endl;
106  }
108  fpSafetyHelper = transportMgr->GetSafetyHelper(); // New
109 
110  // Following assignment is to fix small memory leak from simple use of 'new'
111  static G4ThreadLocal G4TouchableHandle* pNullTouchableHandle = 0;
112  if ( !pNullTouchableHandle) { pNullTouchableHandle = new G4TouchableHandle; }
113  fCurrentTouchableHandle = *pNullTouchableHandle;
114  // Points to (G4VTouchable*) 0
115 
116  G4FieldManager *globalFieldMgr= transportMgr->GetFieldManager();
117  fGlobalFieldExists= globalFieldMgr ? globalFieldMgr->GetDetectorField() : 0 ;
118 }
119 
121 
123 {
124  // fCurrentTouchableHandle is a data member - no deletion required
125 
126  if( (fVerboseLevel > 0) || (fSumEnergyKilled > 0.0 ) )
127  {
128  G4cout << " G4CoupledTransportation: Statistics for looping particles " << G4endl;
129  G4cout << " Sum of energy of loopers killed: " << fSumEnergyKilled << G4endl;
130  G4cout << " Max energy of loopers killed: " << fMaxEnergyKilled << G4endl;
131  }
132 }
133 
135 //
136 // Responsibilities:
137 // Find whether the geometry limits the Step, and to what length
138 // Calculate the new value of the safety and return it.
139 // Store the final time, position and momentum.
140 
143  G4double, // previousStepSize
144  G4double currentMinimumStep,
145  G4double& proposedSafetyForStart,
146  G4GPILSelection* selection )
147 {
148  G4double geometryStepLength;
149  G4double startMassSafety= 0.0; // estimated safety for start point (mass geometry)
150  G4double startFullSafety= 0.0; // estimated safety for start point (all geometries)
151  G4double safetyProposal= -1.0; // local copy of proposal
152 
153  G4ThreeVector EndUnitMomentum ;
154  G4double lengthAlongCurve=0.0 ;
155 
156  fParticleIsLooping = false ;
157 
158  // Initial actions moved to StartTrack()
159  // --------------------------------------
160  // Note: in case another process changes touchable handle
161  // it will be necessary to add here (for all steps)
162  // fCurrentTouchableHandle = aTrack->GetTouchableHandle();
163 
164  // GPILSelection is set to defaule value of CandidateForSelection
165  // It is a return value
166  //
167  *selection = CandidateForSelection ;
168 
169  // Get initial Energy/Momentum of the track
170  //
171  const G4DynamicParticle* pParticle = track.GetDynamicParticle() ;
172  const G4ParticleDefinition* pParticleDef = pParticle->GetDefinition() ;
173  G4ThreeVector startMomentumDir = pParticle->GetMomentumDirection() ;
174  G4ThreeVector startPosition = track.GetPosition() ;
175  G4VPhysicalVolume* currentVolume= track.GetVolume();
176 
177 #ifdef G4DEBUG_TRANSPORT
178  if( fVerboseLevel > 1 )
179  {
180  G4cout << "G4CoupledTransportation::AlongStepGPIL> called in volume "
181  << currentVolume->GetName() << G4endl;
182  }
183 #endif
184  // G4double theTime = track.GetGlobalTime() ;
185 
186  // The Step Point safety can be limited by other geometries and/or the
187  // assumptions of any process - it's not always the geometrical safety.
188  // We calculate the starting point's isotropic safety here.
189  //
190  G4ThreeVector OriginShift = startPosition - fPreviousSftOrigin ;
191  G4double MagSqShift = OriginShift.mag2() ;
192  startMassSafety = 0.0;
193  startFullSafety= 0.0;
194 
195  // Recall that FullSafety <= MassSafety
196  // Original: if( MagSqShift < sqr(fPreviousMassSafety) ) {
197  if( MagSqShift < sqr(fPreviousFullSafety) ) // Revision proposed by Alex H, 2 Oct 07
198  {
199  G4double mag_shift= std::sqrt(MagSqShift);
200  startMassSafety = std::max( (fPreviousMassSafety - mag_shift), 0.0);
201  startFullSafety = std::max( (fPreviousFullSafety - mag_shift), 0.0);
202  // Need to be consistent between full safety with Mass safety
203  // in order reproduce results in simple case --> use same calculation method
204 
205  // Only compute full safety if massSafety > 0. Else it remains 0
206  // startFullSafety = fPathFinder->ComputeSafety( startPosition );
207  }
208 
209  // Is the particle charged or has it a magnetic moment?
210  //
211  G4double particleCharge = pParticle->GetCharge() ;
212  G4double magneticMoment = pParticle->GetMagneticMoment() ;
213  G4double restMass = pParticleDef->GetPDGMass() ;
214 
215  fMassGeometryLimitedStep = false ; // Set default - alex
216  fAnyGeometryLimitedStep = false;
217 
218  // fEndGlobalTimeComputed = false ;
219 
220  // There is no need to locate the current volume. It is Done elsewhere:
221  // On track construction
222  // By the tracking, after all AlongStepDoIts, in "Relocation"
223 
224  // Check if the particle has a force, EM or gravitational, exerted on it
225  //
226  G4FieldManager* fieldMgr=0;
227  G4bool fieldExertsForce = false ;
228 
229  G4bool gravityOn = false;
230  const G4Field* ptrField= 0;
231 
232  fieldMgr = fFieldPropagator->FindAndSetFieldManager( track.GetVolume() );
233  if( fieldMgr != 0 )
234  {
235  // Message the field Manager, to configure it for this track
236  fieldMgr->ConfigureForTrack( &track );
237  // Here it can transition from a null field-ptr to a finite field
238 
239  // If the field manager has no field ptr, the field is zero
240  // by definition ( = there is no field ! )
241  ptrField= fieldMgr->GetDetectorField();
242 
243  if( ptrField != 0)
244  {
245  gravityOn= ptrField->IsGravityActive();
246  if( (particleCharge != 0.0)
247  || (fUseMagneticMoment && (magneticMoment != 0.0) )
248  || (gravityOn && (restMass != 0.0))
249  )
250  {
251  fieldExertsForce = true;
252  }
253  }
254  }
255  G4double momentumMagnitude = pParticle->GetTotalMomentum() ;
256 
257  if( fieldExertsForce )
258  {
259  G4EquationOfMotion* equationOfMotion = 0;
260 
261  // equationOfMotion =
262  // (fFieldPropagator->GetChordFinder()->GetIntegrationDriver()->GetStepper())
263  // ->GetEquationOfMotion();
264 
265  // Consolidate into auxiliary method G4EquationOfMotion* GetEquationOfMotion()
266  // equationOfMotion= fFieldPropagator->GetCurrentEquationOfMotion();
267  G4MagIntegratorStepper* pStepper= 0;
268 
270  if( pChordFinder )
271  {
272  G4MagInt_Driver* pIntDriver= 0;
273 
274  pIntDriver= pChordFinder->GetIntegrationDriver();
275  if( pIntDriver )
276  {
277  pStepper= pIntDriver->GetStepper();
278  }
279  if( pStepper )
280  {
281  equationOfMotion= pStepper->GetEquationOfMotion();
282  }
283  }
284  // End of proto GetEquationOfMotion()
285 
286  G4ChargeState chargeState(particleCharge, // The charge can change (dynamic)
287  magneticMoment,
288  pParticleDef->GetPDGSpin() );
289  // For insurance, could set it again
290  // chargeState.SetPDGSpin( pParticleDef->GetPDGSpin() ); // Newly/provisionally in same object
291 
292  if( equationOfMotion )
293  {
294  equationOfMotion->SetChargeMomentumMass( chargeState,
295  momentumMagnitude,
296  restMass );
297  }
298  }
299 
300  G4ThreeVector polarizationVec = track.GetPolarization() ;
301  G4FieldTrack aFieldTrack = G4FieldTrack( startPosition,
302  track.GetGlobalTime(), // Lab.
303  // track.GetProperTime(), // Particle rest frame
304  track.GetMomentumDirection(),
305  track.GetKineticEnergy(),
306  restMass,
307  particleCharge,
308  polarizationVec,
309  pParticleDef->GetPDGMagneticMoment(),
310  0.0, // Length along track
311  pParticleDef->GetPDGSpin()
312  ) ;
313  G4int stepNo= track.GetCurrentStepNumber();
314 
315  ELimited limitedStep;
316  G4FieldTrack endTrackState('a'); // Default values
317 
318  fMassGeometryLimitedStep = false ; // default
319  fAnyGeometryLimitedStep = false ;
320  if( currentMinimumStep > 0 )
321  {
322  G4double newMassSafety= 0.0; // temp. for recalculation
323 
324  // Do the Transport in the field (non recti-linear)
325  //
326  lengthAlongCurve = fPathFinder->ComputeStep( aFieldTrack,
327  currentMinimumStep,
328  fNavigatorId,
329  stepNo,
330  newMassSafety,
331  limitedStep,
332  endTrackState,
333  currentVolume ) ;
334  // G4cout << " PathFinder ComputeStep returns " << lengthAlongCurve << G4endl;
335 
336  G4double newFullSafety= fPathFinder->GetCurrentSafety();
337  // this was estimated already in step above
338  // G4double newFullStep= fPathFinder->GetMinimumStep();
339 
340  if( limitedStep == kUnique || limitedStep == kSharedTransport )
341  {
342  fMassGeometryLimitedStep = true ;
343  }
344 
346 
347 //#ifdef G4DEBUG_TRANSPORT
349  {
350  G4cerr << " Error in determining geometries limiting the step" << G4endl;
351  G4cerr << " Limiting: mass=" << fMassGeometryLimitedStep
352  << " any= " << fAnyGeometryLimitedStep << G4endl;
353  G4Exception("G4CoupledTransportation::AlongStepGetPhysicalInteractionLength()",
354  "PathFinderConfused", FatalException,
355  "Incompatible conditions - was limited by a geometry?");
356  }
357 //#endif
358 
359  // Other potential
360  // fAnyGeometryLimitedStep = newFullStep < currentMinimumStep;
361  // ^^^ Not good enough;
362  // Must compare with maximum requested step size
363  // (eg in case another process requested bigger, got this!)
364 
365  geometryStepLength = std::min( lengthAlongCurve, currentMinimumStep);
366 
367  // Momentum: Magnitude and direction can be changed too now ...
368  //
369  fMomentumChanged = true ;
370  fTransportEndMomentumDir = endTrackState.GetMomentumDir() ;
371 
372  // Remember last safety origin & value.
373  fPreviousSftOrigin = startPosition ;
374  fPreviousMassSafety = newMassSafety ;
375  fPreviousFullSafety = newFullSafety ;
376  // fpSafetyHelper->SetCurrentSafety( newFullSafety, startPosition);
377 
378 #ifdef G4DEBUG_TRANSPORT
379  if( fVerboseLevel > 1 )
380  {
381  G4cout << "G4Transport:CompStep> "
382  << " called the pathfinder for a new step at " << startPosition
383  << " and obtained step = " << lengthAlongCurve << G4endl;
384  G4cout << " New safety (preStep) = " << newMassSafety
385  << " versus precalculated = " << startMassSafety << G4endl;
386  }
387 #endif
388 
389  // Store as best estimate value
390  startMassSafety = newMassSafety ;
391  startFullSafety = newFullSafety ;
392 
393  // Get the End-Position and End-Momentum (Dir-ection)
394  fTransportEndPosition = endTrackState.GetPosition() ;
395  fTransportEndKineticEnergy = endTrackState.GetKineticEnergy() ;
396  }
397  else
398  {
399  geometryStepLength = lengthAlongCurve= 0.0 ;
400  fMomentumChanged = false ;
401  // fMassGeometryLimitedStep = false ; // --- ???
402  // fAnyGeometryLimitedStep = true;
405 
406  fTransportEndPosition = startPosition;
407 
408  endTrackState= aFieldTrack; // Ensures that time is updated
409 
410  // If the step length requested is 0, and we are on a boundary
411  // then a boundary will also limit the step.
412  if( startMassSafety == 0.0 )
413  {
414  fMassGeometryLimitedStep = true ;
416  }
417  // TODO: Add explicit logical status for being at a boundary
418  }
419  // G4FieldTrack aTrackState(endTrackState);
420 
421  if( !fieldExertsForce )
422  {
423  fParticleIsLooping = false ;
424  fMomentumChanged = false ;
425  fEndGlobalTimeComputed = false ;
426  // G4cout << " global time is false " << G4endl;
427  }
428  else
429  {
430 
431 #ifdef G4DEBUG_TRANSPORT
432  if( fVerboseLevel > 1 )
433  {
434  G4cout << " G4CT::CS End Position = " << fTransportEndPosition << G4endl;
435  G4cout << " G4CT::CS End Direction = " << fTransportEndMomentumDir << G4endl;
436  }
437 #endif
439  {
440  // If the field can change energy, then the time must be integrated
441  // - so this should have been updated
442  //
444  fEndGlobalTimeComputed = true;
445 
446  // was ( fCandidateEndGlobalTime != track.GetGlobalTime() );
447  // a cleaner way is to have FieldTrack knowing whether time is updated.
448  }
449  else
450  {
451  // The energy should be unchanged by field transport,
452  // - so the time changed will be calculated elsewhere
453  //
454  fEndGlobalTimeComputed = false;
455 
456  // Check that the integration preserved the energy
457  // - and if not correct this!
458  G4double startEnergy= track.GetKineticEnergy();
460 
461  static G4ThreadLocal G4int no_inexact_steps=0; // , no_large_ediff;
462  G4double absEdiff = std::fabs(startEnergy- endEnergy);
463  if( absEdiff > perMillion * endEnergy )
464  {
465  no_inexact_steps++;
466  // Possible statistics keeping here ...
467  }
468 #ifdef G4VERBOSE
469  if( (fVerboseLevel > 1) && ( absEdiff > perThousand * endEnergy) )
470  {
471  ReportInexactEnergy(startEnergy, endEnergy);
472  } // end of if (fVerboseLevel)
473 #endif
474  // Correct the energy for fields that conserve it
475  // This - hides the integration error
476  // - but gives a better physical answer
478  }
479  }
480 
481  fEndpointDistance = (fTransportEndPosition - startPosition).mag() ;
483 
484  fTransportEndSpin = endTrackState.GetSpin();
485 
486  // Calculate the safety
487  safetyProposal= startFullSafety; // used to be startMassSafety
488  // Changed to accomodate processes that cannot update the safety -- JA 22 Nov 06
489 
490  // Update safety for the end-point, if becomes negative at the end-point.
491 
492  if( (startFullSafety < fEndpointDistance )
493  && ( particleCharge != 0.0 ) ) // Only needed to prepare for Mult Scat.
494  // && !fAnyGeometryLimitedStep ) // To-Try: No safety update if at a boundary
495  {
496  G4double endFullSafety =
498  // Expected mission -- only mass geometry's safety
499  // fMassNavigator->ComputeSafety( fTransportEndPosition) ;
500  // Yet discrete processes only have poststep -- and this cannot
501  // currently revise the safety
502  // ==> so we use the all-geometry safety as a precaution
503 
505  // Pushing safety to Helper avoids recalculation at this point
506 
507  G4ThreeVector centerPt= G4ThreeVector(0.0, 0.0, 0.0); // Used for return value
508  G4double endMassSafety= fPathFinder->ObtainSafety( fNavigatorId, centerPt);
509  // Retrieves the mass value from PathFinder (it calculated it)
510 
511  fPreviousMassSafety = endMassSafety ;
512  fPreviousFullSafety = endFullSafety;
513  fPreviousSftOrigin = fTransportEndPosition ;
514 
515  // The convention (Stepping Manager's) is safety from the start point
516  //
517  safetyProposal = endFullSafety + fEndpointDistance;
518  // --> was endMassSafety
519  // Changed to accomodate processes that cannot update the safety -- JA 22 Nov 06
520 
521  // #define G4DEBUG_TRANSPORT 1
522 
523 #ifdef G4DEBUG_TRANSPORT
524  G4int prec= G4cout.precision(12) ;
525  G4cout << "***Transportation::AlongStepGPIL ** " << G4endl ;
526  G4cout << " Revised Safety at endpoint " << fTransportEndPosition
527  << " give safety values: Mass= " << endMassSafety
528  << " All= " << endFullSafety << G4endl ;
529  G4cout << " Adding endpoint distance " << fEndpointDistance
530  << " to obtain pseudo-safety= " << safetyProposal << G4endl ;
531  G4cout.precision(prec);
532  }
533  else
534  {
535  G4int prec= G4cout.precision(12) ;
536  G4cout << "***Transportation::AlongStepGPIL ** " << G4endl ;
537  G4cout << " Quick Safety estimate at endpoint " << fTransportEndPosition
538  << " gives safety endpoint value = " << startFullSafety - fEndpointDistance
539  << " using start-point value " << startFullSafety
540  << " and endpointDistance " << fEndpointDistance << G4endl;
541  G4cout.precision(prec);
542 #endif
543  }
544 
545  proposedSafetyForStart= safetyProposal;
546  fParticleChange.ProposeTrueStepLength(geometryStepLength) ;
547 
548  return geometryStepLength ;
549 }
550 
552 
555  const G4Step& stepData )
556 {
557  static G4ThreadLocal G4int noCalls=0;
558  noCalls++;
559 
560  fParticleChange.Initialize(track) ;
561  // sets all its members to the value of corresponding members in G4Track
562 
563  // Code specific for Transport
564  //
566  // G4cout << " G4CoupledTransportation::AlongStepDoIt"
567  // << " proposes position = " << fTransportEndPosition
568  // << " and end momentum direction = " << fTransportEndMomentumDir << G4endl;
572 
574 
575  G4double deltaTime = 0.0 ;
576 
577  // Calculate Lab Time of Flight (ONLY if field Equations used it!)
578  // G4double endTime = fCandidateEndGlobalTime;
579  // G4double delta_time = endTime - startTime;
580 
581  G4double startTime = track.GetGlobalTime() ;
582 
584  {
585  G4double finalInverseVel= DBL_MAX, initialInverseVel=DBL_MAX;
586 
587  // The time was not integrated .. make the best estimate possible
588  //
589  G4double finalVelocity = track.GetVelocity() ;
590  if( finalVelocity > 0.0 ) { finalInverseVel= 1.0 / finalVelocity; }
591  G4double initialVelocity = stepData.GetPreStepPoint()->GetVelocity() ;
592  if( initialVelocity > 0.0 ) { initialInverseVel= 1.0 / initialVelocity; }
593  G4double stepLength = track.GetStepLength() ;
594 
595  if (finalVelocity > 0.0)
596  {
597  // deltaTime = stepLength/finalVelocity ;
598  G4double meanInverseVelocity = 0.5 * ( initialInverseVel + finalInverseVel );
599  deltaTime = stepLength * meanInverseVelocity ;
600  // G4cout << " dt = s * mean(1/v) , with " << " s = " << stepLength
601  // << " mean(1/v)= " << meanInverseVelocity << G4endl;
602  }
603  else
604  {
605  deltaTime = stepLength * initialInverseVel ;
606  // G4cout << " dt = s / initV " << " s = " << stepLength
607  // << " 1 / initV= " << initialInverseVel << G4endl;
608  } // Could do with better estimate for final step (finalVelocity = 0) ?
609 
610  fCandidateEndGlobalTime = startTime + deltaTime ;
611  fParticleChange.ProposeLocalTime( track.GetLocalTime() + deltaTime) ;
612 
613  // G4cout << " Calculated global time from start = " << startTime << " and "
614  // << " delta time = " << deltaTime << G4endl;
615  }
616  else
617  {
618  deltaTime = fCandidateEndGlobalTime - startTime ;
620  // G4cout << " Calculated global time from candidate end time = "
621  // << fCandidateEndGlobalTime << " and start time = " << startTime << G4endl;
622  }
623 
624  // G4cout << " G4CoupledTransportation::AlongStepDoIt "
625  // << " flag whether computed time = " << fEndGlobalTimeComputed << " and "
626  // << " is proposes end time " << fCandidateEndGlobalTime << G4endl;
627 
628  // Now Correct by Lorentz factor to get "proper" deltaTime
629 
630  G4double restMass = track.GetDynamicParticle()->GetMass() ;
631  G4double deltaProperTime = deltaTime*( restMass/track.GetTotalEnergy() ) ;
632 
633  fParticleChange.ProposeProperTime(track.GetProperTime() + deltaProperTime) ;
634  //fParticleChange. ProposeTrueStepLength( track.GetStepLength() ) ;
635 
636  // If the particle is caught looping or is stuck (in very difficult
637  // boundaries) in a magnetic field (doing many steps)
638  // THEN this kills it ...
639  //
640  if ( fParticleIsLooping )
641  {
643 
644  if( (endEnergy < fThreshold_Important_Energy)
646  {
647  // Kill the looping particle
648  //
650 
651  // 'Bare' statistics
652  fSumEnergyKilled += endEnergy;
653  if( endEnergy > fMaxEnergyKilled) { fMaxEnergyKilled= endEnergy; }
654 
655 #ifdef G4VERBOSE
656  if((fVerboseLevel > 1) && ( endEnergy > fThreshold_Warning_Energy ))
657  {
658  G4cout << " G4CoupledTransportation is killing track that is looping or stuck " << G4endl
659  << " This track has " << track.GetKineticEnergy() / MeV
660  << " MeV energy." << G4endl;
661  }
662  if( fVerboseLevel > 0 )
663  {
664  G4cout << " Steps by this track: " << track.GetCurrentStepNumber() << G4endl;
665  }
666 #endif
667  fNoLooperTrials=0;
668  }
669  else
670  {
671  fNoLooperTrials ++;
672 #ifdef G4VERBOSE
673  if( (fVerboseLevel > 2) )
674  {
675  G4cout << " ** G4CoupledTransportation::AlongStepDoIt(): Particle looping - " << G4endl
676  << " Number of consecutive problem step (this track) = " << fNoLooperTrials << G4endl
677  << " Steps by this track: " << track.GetCurrentStepNumber() << G4endl
678  << " Total no of calls to this method (all tracks) = " << noCalls << G4endl;
679  }
680 #endif
681  }
682  }
683  else
684  {
685  fNoLooperTrials=0;
686  }
687 
688  // Another (sometimes better way) is to use a user-limit maximum Step size
689  // to alleviate this problem ..
690 
691  // Add smooth curved trajectories to particle-change
692  //
693  // fParticleChange.SetPointerToVectorOfAuxiliaryPoints
694  // (fFieldPropagator->GimmeTrajectoryVectorAndForgetIt() );
695 
696  return &fParticleChange ;
697 }
698 
700 //
701 // This ensures that the PostStep action is always called,
702 // so that it can do the relocation if it is needed.
703 //
704 
707  G4double, // previousStepSize
708  G4ForceCondition* pForceCond )
709 {
710  // Must act as PostStep action -- to relocate particle
711  *pForceCond = Forced ;
712  return DBL_MAX ;
713 }
714 
716 ReportMove( G4ThreeVector OldVector, G4ThreeVector NewVector, const G4String& Quantity )
717 {
718  G4ThreeVector moveVec = ( NewVector - OldVector );
719 
720  G4cerr << G4endl
721  << "**************************************************************" << G4endl;
722  G4cerr << "Endpoint has moved between value expected from TransportEndPosition "
723  << " and value from Track in PostStepDoIt. " << G4endl
724  << "Change of " << Quantity << " is " << moveVec.mag() / mm << " mm long, "
725  << " and its vector is " << (1.0/mm) * moveVec << " mm " << G4endl
726  << "Endpoint of ComputeStep was " << OldVector
727  << " and current position to locate is " << NewVector << G4endl;
728 }
729 
731 
733  const G4Step& )
734 {
735  G4TouchableHandle retCurrentTouchable ; // The one to return
736 
737  // Initialize ParticleChange (by setting all its members equal
738  // to corresponding members in G4Track)
739  // fParticleChange.Initialize(track) ; // To initialise TouchableChange
740 
742 
743  // Check that the end position and direction are preserved
744  // since call to AlongStepDoIt
745 
746 #ifdef G4DEBUG_TRANSPORT
747  if( ( fVerboseLevel > 0 )
748  && ((fTransportEndPosition - track.GetPosition()).mag2() >= 1.0e-16) )
749  {
750  ReportMove( track.GetPosition(), fTransportEndPosition, "End of Step Position" );
751  G4cerr << " Problem in G4CoupledTransportation::PostStepDoIt " << G4endl;
752  }
753 
754  // If the Step was determined by the volume boundary, relocate the particle
755  // The pathFinder will know that the geometry limited the step (!?)
756 
757  if( fVerboseLevel > 0 )
758  {
759  G4cout << " Calling PathFinder::Locate() from "
760  << " G4CoupledTransportation::PostStepDoIt " << G4endl;
761  G4cout << " fAnyGeometryLimitedStep is " << fAnyGeometryLimitedStep << G4endl;
762 
763  }
764 #endif
765 
767  {
768  fPathFinder->Locate( track.GetPosition(),
769  track.GetMomentumDirection(),
770  true);
771 
772  // fCurrentTouchable will now become the previous touchable,
773  // and what was the previous will be freed.
774  // (Needed because the preStepPoint can point to the previous touchable)
775 
778 
779 #ifdef G4DEBUG_TRANSPORT
780  if( fVerboseLevel > 0 )
781  {
782  G4cout << "G4CoupledTransportation::PostStepDoIt --- fNavigatorId = "
783  << fNavigatorId << G4endl;
784  }
785  if( fVerboseLevel > 1 )
786  {
788  G4cout << "CHECK !!!!!!!!!!! fCurrentTouchableHandle->GetVolume() = " << vol;
789  if( vol ) { G4cout << "Name=" << vol->GetName(); }
790  G4cout << G4endl;
791  }
792 #endif
793 
794  // Check whether the particle is out of the world volume
795  // If so it has exited and must be killed.
796  //
797  if( fCurrentTouchableHandle->GetVolume() == 0 )
798  {
800  }
801  retCurrentTouchable = fCurrentTouchableHandle ;
802  // fParticleChange.SetTouchableHandle( fCurrentTouchableHandle ) ;
803 
804  // Notify particle change that this is last step in volume
806  }
807  else // fAnyGeometryLimitedStep is false
808  {
809 #ifdef G4DEBUG_TRANSPORT
810  if( fVerboseLevel > 1 )
811  {
812  G4cout << "G4CoupledTransportation::PostStepDoIt -- "
813  << " fAnyGeometryLimitedStep = " << fAnyGeometryLimitedStep
814  << " must be false " << G4endl;
815  }
816 #endif
817  // This serves only to move each of the Navigator's location
818  //
819  // fLinearNavigator->LocateGlobalPointWithinVolume( track.GetPosition() ) ;
820 
821  // G4cout << "G4CoupledTransportation calling PathFinder::ReLocate() " << G4endl;
822  fPathFinder->ReLocate( track.GetPosition() );
823  // track.GetMomentumDirection() );
824 
825  // Keep the value of the track's current Touchable is retained,
826  // and use it to overwrite the (unset) one in particle change.
827  // Expect this must be fCurrentTouchable too
828  // - could it be different, eg at the start of a step ?
829  //
830  retCurrentTouchable = track.GetTouchableHandle() ;
831  // fParticleChange.SetTouchableHandle( track.GetTouchableHandle() ) ;
832 
833  // Have not reached a boundary
835  } // endif ( fAnyGeometryLimitedStep )
836 
837  const G4VPhysicalVolume* pNewVol = retCurrentTouchable->GetVolume() ;
838  const G4Material* pNewMaterial = 0 ;
839  const G4VSensitiveDetector* pNewSensitiveDetector = 0 ;
840 
841  if( pNewVol != 0 )
842  {
843  pNewMaterial= pNewVol->GetLogicalVolume()->GetMaterial();
844  pNewSensitiveDetector= pNewVol->GetLogicalVolume()->GetSensitiveDetector();
845  }
846 
847  // ( const_cast<G4Material *> pNewMaterial ) ;
848  // ( const_cast<G4VSensitiveDetetor *> pNewSensitiveDetector) ;
849 
852  // "temporarily" until Get/Set Material of ParticleChange,
853  // and StepPoint can be made const.
854 
855  const G4MaterialCutsCouple* pNewMaterialCutsCouple = 0;
856  if( pNewVol != 0 )
857  {
858  pNewMaterialCutsCouple=pNewVol->GetLogicalVolume()->GetMaterialCutsCouple();
859  if( pNewMaterialCutsCouple!=0
860  && pNewMaterialCutsCouple->GetMaterial()!=pNewMaterial )
861  {
862  // for parametrized volume
863  //
864  pNewMaterialCutsCouple =
866  ->GetMaterialCutsCouple(pNewMaterial,
867  pNewMaterialCutsCouple->GetProductionCuts());
868  }
869  }
870  fParticleChange.SetMaterialCutsCoupleInTouchable( pNewMaterialCutsCouple );
871 
872  // Must always set the touchable in ParticleChange, whether relocated or not
873  fParticleChange.SetTouchableHandle(retCurrentTouchable) ;
874 
875  return &fParticleChange ;
876 }
877 
878 // New method takes over the responsibility to reset the state of
879 // G4CoupledTransportation object:
880 // - at the start of a new track, and
881 // - on the resumption of a suspended track.
882 
883 void
885 {
886 
887  G4TransportationManager* transportMgr =
889 
890  // G4VProcess::StartTracking(aTrack);
891 
892  // The 'initialising' actions
893  // once taken in AlongStepGPIL -- if ( track.GetCurrentStepNumber()==1 )
894 
895  // fStartedNewTrack= true;
896 
897  fMassNavigator = transportMgr->GetNavigatorForTracking() ;
898  fNavigatorId= transportMgr->ActivateNavigator( fMassNavigator ); // Confirm it!
899 
900  // if( fVerboseLevel > 1 ){
901  // G4cout << " Navigator Id obtained in StartTracking " << fNavigatorId << G4endl;
902  // }
903  G4ThreeVector position = aTrack->GetPosition();
904  G4ThreeVector direction = aTrack->GetMomentumDirection();
905 
906  // if( fVerboseLevel > 1 ){
907  // G4cout << " Calling PathFinder::PrepareNewTrack from "
908  // << " G4CoupledTransportation::StartTracking -- which calls Locate()" << G4endl;
909  // }
910  fPathFinder->PrepareNewTrack( position, direction);
911  // This implies a call to fPathFinder->Locate( position, direction );
912 
913  // Global field, if any, must exist before tracking is started
915  // reset safety value and center
916  //
917  fPreviousMassSafety = 0.0 ;
918  fPreviousFullSafety = 0.0 ;
919  fPreviousSftOrigin = G4ThreeVector(0.,0.,0.) ;
920 
921  // reset looping counter -- for motion in field
922  fNoLooperTrials= 0;
923  // Must clear this state .. else it depends on last track's value
924  // --> a better solution would set this from state of suspended track TODO ?
925  // Was if( aTrack->GetCurrentStepNumber()==1 ) { .. }
926 
927  // ChordFinder reset internal state
928  //
929  if( fGlobalFieldExists )
930  {
932  // Resets safety values, in case of overlaps.
933 
935  if( chordF ) { chordF->ResetStepEstimate(); }
936  }
937 
938  // Clear the chord finders of all fields (ie managers) derived objects
939  //
941  fieldMgrStore->ClearAllChordFindersState();
942 
943 #ifdef G4DEBUG_TRANSPORT
944  if( fVerboseLevel > 1 )
945  {
946  G4cout << " Returning touchable handle " << fCurrentTouchableHandle << G4endl;
947  }
948 #endif
949 
950  // Update the current touchable handle (from the track's)
951  //
953 }
954 
955 void
957 {
959  fPathFinder->EndTrack(); // Resets TransportationManager to use ordinary Navigator
960 }
961 
962 void
964 ReportInexactEnergy(G4double startEnergy, G4double endEnergy)
965 {
966  static G4ThreadLocal G4int no_warnings= 0, warnModulo=1, moduloFactor= 10, no_large_ediff= 0;
967 
968  if( std::fabs(startEnergy- endEnergy) > perThousand * endEnergy )
969  {
970  no_large_ediff ++;
971  if( (no_large_ediff% warnModulo) == 0 )
972  {
973  no_warnings++;
974  G4cout << "WARNING - G4CoupledTransportation::AlongStepGetPIL() "
975  << " Energy change in Step is above 1^-3 relative value. " << G4endl
976  << " Relative change in 'tracking' step = "
977  << std::setw(15) << (endEnergy-startEnergy)/startEnergy << G4endl
978  << " Starting E= " << std::setw(12) << startEnergy / MeV << " MeV " << G4endl
979  << " Ending E= " << std::setw(12) << endEnergy / MeV << " MeV " << G4endl;
980  G4cout << " Energy has been corrected -- however, review"
981  << " field propagation parameters for accuracy." << G4endl;
982  if( (fVerboseLevel > 2 ) || (no_warnings<4) || (no_large_ediff == warnModulo * moduloFactor) )
983  {
984  G4cout << " These include EpsilonStepMax(/Min) in G4FieldManager "
985  << " which determine fractional error per step for integrated quantities. " << G4endl
986  << " Note also the influence of the permitted number of integration steps."
987  << G4endl;
988  }
989  G4cerr << "ERROR - G4CoupledTransportation::AlongStepGetPIL()" << G4endl
990  << " Bad 'endpoint'. Energy change detected"
991  << " and corrected. "
992  << " Has occurred already "
993  << no_large_ediff << " times." << G4endl;
994  if( no_large_ediff == warnModulo * moduloFactor )
995  {
996  warnModulo *= moduloFactor;
997  }
998  }
999  }
1000 }
1001 
1002 #include "G4Transportation.hh"
1004 {
1005  G4bool lastValue= fUseMagneticMoment;
1006  fUseMagneticMoment= useMoment;
1008  return lastValue;
1009 }
void PrepareNewTrack(const G4ThreeVector &position, const G4ThreeVector &direction, G4VPhysicalVolume *massStartVol=0)
static G4PathFinder * GetInstance()
Definition: G4PathFinder.cc:57
void SetMaterialInTouchable(G4Material *fMaterial)
void SetTouchableHandle(const G4TouchableHandle &fTouchable)
const G4ThreeVector & GetPolarization() const
G4double ObtainSafety(G4int navId, G4ThreeVector &globalCenterPoint)
static const double MeV
Definition: G4SIunits.hh:193
void Locate(const G4ThreeVector &position, const G4ThreeVector &direction, G4bool relativeSearch=true)
G4SafetyHelper * GetSafetyHelper() const
G4double GetLocalTime() const
void ReportMove(G4ThreeVector OldVector, G4ThreeVector NewVector, const G4String &Quantity)
void SetMaterialCutsCoupleInTouchable(const G4MaterialCutsCouple *fMaterialCutsCouple)
void StartTracking(G4Track *aTrack)
virtual void SetChargeMomentumMass(G4ChargeState particleCharge, G4double MomentumXc, G4double MassXc2)=0
G4double GetProperTime() const
CLHEP::Hep3Vector G4ThreeVector
G4double GetVelocity() const
G4double GetCurrentSafety() const
G4double GetKineticEnergy() const
G4TouchableHandle fCurrentTouchableHandle
const G4DynamicParticle * GetDynamicParticle() const
const G4MagIntegratorStepper * GetStepper() const
G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData)
G4Material * GetMaterial() const
void ReLocate(const G4ThreeVector &position)
ELimited
const G4ThreeVector & GetPosition() const
const G4ThreeVector & GetMomentumDir() const
G4TrackStatus GetTrackStatus() const
G4ThreeVector GetSpin() const
G4Navigator * GetNavigatorForTracking() const
G4TouchableHandle CreateTouchableHandle(G4int navId) const
void ProposePolarization(G4double Px, G4double Py, G4double Pz)
G4ParticleDefinition * GetDefinition() const
static const double perThousand
Definition: G4SIunits.hh:297
G4double GetVelocity() const
#define G4ThreadLocal
Definition: tls.hh:84
int G4int
Definition: G4Types.hh:78
void ProposePosition(G4double x, G4double y, G4double z)
G4ReferenceCountedHandle< G4VTouchable > G4TouchableHandle
unsigned int GetNumberGeometriesLimitingStep() const
G4double GetTotalMomentum() const
static G4bool EnableUseMagneticMoment(G4bool useMoment=true)
G4StepPoint * GetPreStepPoint() const
virtual void Initialize(const G4Track &)
virtual void ConfigureForTrack(const G4Track *)
static const double prec
Definition: RanecuEngine.cc:58
G4double GetKineticEnergy() const
G4EquationOfMotion * GetEquationOfMotion()
#define position
Definition: xmlparse.cc:605
G4ThreeVector GetPosition() const
G4GLOB_DLL std::ostream G4cout
void ReportInexactEnergy(G4double startEnergy, G4double endEnergy)
G4int GetCurrentStepNumber() const
const G4String & GetName() const
G4double GetMass() const
bool G4bool
Definition: G4Types.hh:79
const G4ThreeVector & GetMomentumDirection() const
G4double GetCharge() const
void ProposeTrueStepLength(G4double truePathLength)
G4double ComputeStep(const G4FieldTrack &pFieldTrack, G4double pCurrentProposedStepLength, G4int navigatorId, G4int stepNo, G4double &pNewSafety, ELimited &limitedStep, G4FieldTrack &EndState, G4VPhysicalVolume *currentVolume)
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
G4bool DoesFieldChangeEnergy() const
Definition: G4Step.hh:76
G4double ComputeSafety(const G4ThreeVector &globalPoint)
G4double GetGlobalTime() const
G4FieldManager * FindAndSetFieldManager(G4VPhysicalVolume *pCurrentPhysVol)
const G4TouchableHandle & GetTouchableHandle() const
G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4TransportationManager * GetTransportationManager()
G4FieldManager * GetFieldManager() const
static const double perMillion
Definition: G4SIunits.hh:298
static G4ProductionCutsTable * GetProductionCutsTable()
G4int ActivateNavigator(G4Navigator *aNavigator)
const G4ThreeVector & GetMomentumDirection() const
G4LogicalVolume * GetLogicalVolume() const
void ProposeProperTime(G4double finalProperTime)
G4double GetPDGMass() const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
static G4FieldManagerStore * GetInstance()
void SetSensitiveDetectorInTouchable(G4VSensitiveDetector *fSensitiveDetector)
#define fPreviousSftOrigin
G4bool IsParticleLooping() const
G4double GetLabTimeOfFlight() const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:44
G4bool IsGravityActive() const
Definition: G4Field.hh:98
void ProposeGlobalTime(G4double t)
G4ChordFinder * GetChordFinder()
G4FieldManager * GetCurrentFieldManager()
void ProposeEnergy(G4double finalEnergy)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4VPhysicalVolume * GetVolume() const
G4double GetTotalEnergy() const
G4double GetPDGSpin() const
void SetCurrentSafety(G4double val, const G4ThreeVector &pos)
#define G4endl
Definition: G4ios.hh:61
void ProposeLastStepInVolume(G4bool flag)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
T sqr(const T &x)
Definition: templates.hh:145
G4double GetPDGMagneticMoment() const
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
const G4Field * GetDetectorField() const
G4ForceCondition
G4MagInt_Driver * GetIntegrationDriver()
void SetMomentumChanged(G4bool b)
G4PropagatorInField * GetPropagatorInField() const
G4ProductionCuts * GetProductionCuts() const
#define DBL_MAX
Definition: templates.hh:83
G4VSensitiveDetector * GetSensitiveDetector() const
static const double mm
Definition: G4SIunits.hh:102
G4double GetMagneticMoment() const
G4double PostStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4ForceCondition *pForceCond)
void ProposeLocalTime(G4double t)
G4PropagatorInField * fFieldPropagator
G4double GetStepLength() const
void ResetStepEstimate()
G4GPILSelection
const G4Material * GetMaterial() const
G4GLOB_DLL std::ostream G4cerr
G4ParticleChangeForTransport fParticleChange
static G4bool fUseMagneticMoment
G4CoupledTransportation(G4int verbosityLevel=0)