Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 69887 2013-05-17 08:17:02Z 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 "G4FieldManagerStore.hh"
57 
59 
61 //
62 // Constructor
63 
65  : G4VProcess( G4String("CoupledTransportation"), fTransportation ),
66  fParticleIsLooping( false ),
67  fPreviousSftOrigin( 0.,0.,0. ),
68  fPreviousMassSafety( 0.0 ),
69  fPreviousFullSafety( 0.0 ),
70  fThreshold_Warning_Energy( 100 * MeV ),
71  fThreshold_Important_Energy( 250 * MeV ),
72  fThresholdTrials( 10 ),
73  fNoLooperTrials( 0 ),
74  fSumEnergyKilled( 0.0 ), fMaxEnergyKilled( 0.0 ),
75  fUseMagneticMoment( false ),
76  fVerboseLevel( verbosity )
77 {
78  // set Process Sub Type
79  SetProcessSubType(static_cast<G4int>(COUPLED_TRANSPORTATION));
80 
81  G4TransportationManager* transportMgr ;
82 
84 
85  fMassNavigator = transportMgr->GetNavigatorForTracking() ;
86  fFieldPropagator = transportMgr->GetPropagatorInField() ;
87  // fGlobalFieldMgr = transportMgr->GetFieldManager() ;
88  fNavigatorId= transportMgr->ActivateNavigator( fMassNavigator );
89  if( fVerboseLevel > 0 )
90  {
91  G4cout << " G4CoupledTransportation constructor: ----- " << G4endl;
92  G4cout << " Verbose level is " << fVerboseLevel << G4endl;
93  G4cout << " Navigator Id obtained in G4CoupledTransportation constructor "
94  << fNavigatorId << G4endl;
95  }
96  fPathFinder= G4PathFinder::GetInstance();
97  fpSafetyHelper = transportMgr->GetSafetyHelper(); // New
98 
99  // Following assignment is to fix small memory leak from simple use of 'new'
100  static G4TouchableHandle nullTouchableHandle; // Points to (G4VTouchable*) 0
101  fCurrentTouchableHandle = nullTouchableHandle;
102 
103  fEndGlobalTimeComputed = false;
104  fCandidateEndGlobalTime = 0;
105 }
106 
108 
110 {
111  // fCurrentTouchableHandle is a data member - no deletion required
112 
113  if( (fVerboseLevel > 0) || (fSumEnergyKilled > 0.0 ) )
114  {
115  G4cout << " G4CoupledTransportation: Statistics for looping particles " << G4endl;
116  G4cout << " Sum of energy of loopers killed: " << fSumEnergyKilled << G4endl;
117  G4cout << " Max energy of loopers killed: " << fMaxEnergyKilled << G4endl;
118  }
119 }
120 
122 //
123 // Responsibilities:
124 // Find whether the geometry limits the Step, and to what length
125 // Calculate the new value of the safety and return it.
126 // Store the final time, position and momentum.
127 
130  G4double, // previousStepSize
131  G4double currentMinimumStep,
132  G4double& proposedSafetyForStart,
133  G4GPILSelection* selection )
134 {
135  G4double geometryStepLength;
136  G4double startMassSafety= 0.0; // estimated safety for start point (mass geometry)
137  G4double startFullSafety= 0.0; // estimated safety for start point (all geometries)
138  G4double safetyProposal= -1.0; // local copy of proposal
139 
140  G4ThreeVector EndUnitMomentum ;
141  G4double lengthAlongCurve=0.0 ;
142 
143  fParticleIsLooping = false ;
144 
145  // Initial actions moved to StartTrack()
146  // --------------------------------------
147  // Note: in case another process changes touchable handle
148  // it will be necessary to add here (for all steps)
149  // fCurrentTouchableHandle = aTrack->GetTouchableHandle();
150 
151  // GPILSelection is set to defaule value of CandidateForSelection
152  // It is a return value
153  //
154  *selection = CandidateForSelection ;
155 
156  // Get initial Energy/Momentum of the track
157  //
158  const G4DynamicParticle* pParticle = track.GetDynamicParticle() ;
159  const G4ParticleDefinition* pParticleDef = pParticle->GetDefinition() ;
160  G4ThreeVector startMomentumDir = pParticle->GetMomentumDirection() ;
161  G4ThreeVector startPosition = track.GetPosition() ;
162  G4VPhysicalVolume* currentVolume= track.GetVolume();
163 
164 #ifdef G4DEBUG_TRANSPORT
165  if( fVerboseLevel > 1 )
166  {
167  G4cout << "G4CoupledTransportation::AlongStepGPIL> called in volume "
168  << currentVolume->GetName() << G4endl;
169  }
170 #endif
171  // G4double theTime = track.GetGlobalTime() ;
172 
173  // The Step Point safety can be limited by other geometries and/or the
174  // assumptions of any process - it's not always the geometrical safety.
175  // We calculate the starting point's isotropic safety here.
176  //
177  G4ThreeVector OriginShift = startPosition - fPreviousSftOrigin ;
178  G4double MagSqShift = OriginShift.mag2() ;
179  startMassSafety = 0.0;
180  startFullSafety= 0.0;
181 
182  // Recall that FullSafety <= MassSafety
183  // Original: if( MagSqShift < sqr(fPreviousMassSafety) ) {
184  if( MagSqShift < sqr(fPreviousFullSafety) ) // Revision proposed by Alex H, 2 Oct 07
185  {
186  G4double mag_shift= std::sqrt(MagSqShift);
187  startMassSafety = std::max( (fPreviousMassSafety - mag_shift), 0.0);
188  startFullSafety = std::max( (fPreviousFullSafety - mag_shift), 0.0);
189  // Need to be consistent between full safety with Mass safety
190  // in order reproduce results in simple case --> use same calculation method
191 
192  // Only compute full safety if massSafety > 0. Else it remains 0
193  // startFullSafety = fPathFinder->ComputeSafety( startPosition );
194  }
195 
196  // Is the particle charged or has it a magnetic moment?
197  //
198  G4double particleCharge = pParticle->GetCharge() ;
199  G4double magneticMoment = pParticle->GetMagneticMoment() ;
200  G4double restMass = pParticleDef->GetPDGMass() ;
201 
202  fMassGeometryLimitedStep = false ; // Set default - alex
203  fAnyGeometryLimitedStep = false;
204 
205  // fEndGlobalTimeComputed = false ;
206 
207  // There is no need to locate the current volume. It is Done elsewhere:
208  // On track construction
209  // By the tracking, after all AlongStepDoIts, in "Relocation"
210 
211  // Check if the particle has a force, EM or gravitational, exerted on it
212  //
213  G4FieldManager* fieldMgr=0;
214  G4bool fieldExertsForce = false ;
215 
216  G4bool gravityOn = false;
217  const G4Field* ptrField= 0;
218 
219  fieldMgr = fFieldPropagator->FindAndSetFieldManager( track.GetVolume() );
220  if( fieldMgr != 0 )
221  {
222  // Message the field Manager, to configure it for this track
223  fieldMgr->ConfigureForTrack( &track );
224  // Here it can transition from a null field-ptr to a finite field
225 
226  // If the field manager has no field ptr, the field is zero
227  // by definition ( = there is no field ! )
228  ptrField= fieldMgr->GetDetectorField();
229 
230  if( ptrField != 0)
231  {
232  gravityOn= ptrField->IsGravityActive();
233  if( (particleCharge != 0.0)
234  || (fUseMagneticMoment && (magneticMoment != 0.0) )
235  || (gravityOn && (restMass != 0.0))
236  )
237  {
238  fieldExertsForce = true;
239  }
240  }
241  }
242  G4double momentumMagnitude = pParticle->GetTotalMomentum() ;
243 
244  fFieldPropagator->SetChargeMomentumMass( particleCharge, // in e+ units
245  momentumMagnitude, // in Mev/c
246  restMass ) ;
247  // This should be obsolete - the call to SetChargeAndMoments below should do the work
248 
249  G4ThreeVector spin = track.GetPolarization() ;
250  G4FieldTrack theFieldTrack = G4FieldTrack( startPosition,
251  track.GetMomentumDirection(),
252  0.0,
253  track.GetKineticEnergy(),
254  restMass,
255  0.0, // UNUSED: track.GetVelocity(),
256  track.GetGlobalTime(), // Lab.
257  track.GetProperTime(), // Part.
258  &spin ) ;
259  theFieldTrack.SetChargeAndMoments( particleCharge ); // EM moments -- future extension
260 
261  G4int stepNo= track.GetCurrentStepNumber();
262 
263  ELimited limitedStep;
264  G4FieldTrack endTrackState('a'); // Default values
265 
266  fMassGeometryLimitedStep = false ; // default
267  fAnyGeometryLimitedStep = false ;
268  if( currentMinimumStep > 0 )
269  {
270  G4double newMassSafety= 0.0; // temp. for recalculation
271 
272  // Do the Transport in the field (non recti-linear)
273  //
274  lengthAlongCurve = fPathFinder->ComputeStep( theFieldTrack,
275  currentMinimumStep,
276  fNavigatorId,
277  stepNo,
278  newMassSafety,
279  limitedStep,
280  endTrackState,
281  currentVolume ) ;
282  // G4cout << " PathFinder ComputeStep returns " << lengthAlongCurve << G4endl;
283 
284  G4double newFullSafety= fPathFinder->GetCurrentSafety();
285  // this was estimated already in step above
286  // G4double newFullStep= fPathFinder->GetMinimumStep();
287 
288  if( limitedStep == kUnique || limitedStep == kSharedTransport )
289  {
290  fMassGeometryLimitedStep = true ;
291  }
292 
293  fAnyGeometryLimitedStep = (fPathFinder->GetNumberGeometriesLimitingStep() != 0);
294 
295 //#ifdef G4DEBUG_TRANSPORT
296  if( fMassGeometryLimitedStep && !fAnyGeometryLimitedStep )
297  {
298  G4cerr << " Error in determining geometries limiting the step" << G4endl;
299  G4cerr << " Limiting: mass=" << fMassGeometryLimitedStep
300  << " any= " << fAnyGeometryLimitedStep << G4endl;
301  G4Exception("G4CoupledTransportation::AlongStepGetPhysicalInteractionLength()",
302  "PathFinderConfused", FatalException,
303  "Incompatible conditions - was limited by a geometry?");
304  }
305 //#endif
306 
307  // Other potential
308  // fAnyGeometryLimitedStep = newFullStep < currentMinimumStep;
309  // ^^^ Not good enough;
310  // Must compare with maximum requested step size
311  // (eg in case another process requested bigger, got this!)
312 
313  geometryStepLength = std::min( lengthAlongCurve, currentMinimumStep);
314 
315  // Momentum: Magnitude and direction can be changed too now ...
316  //
317  fMomentumChanged = true ;
318  fTransportEndMomentumDir = endTrackState.GetMomentumDir() ;
319 
320  // Remember last safety origin & value.
321  fPreviousSftOrigin = startPosition ;
322  fPreviousMassSafety = newMassSafety ;
323  fPreviousFullSafety = newFullSafety ;
324  // fpSafetyHelper->SetCurrentSafety( newFullSafety, startPosition);
325 
326 #ifdef G4DEBUG_TRANSPORT
327  if( fVerboseLevel > 1 )
328  {
329  G4cout << "G4Transport:CompStep> "
330  << " called the pathfinder for a new step at " << startPosition
331  << " and obtained step = " << lengthAlongCurve << G4endl;
332  G4cout << " New safety (preStep) = " << newMassSafety
333  << " versus precalculated = " << startMassSafety << G4endl;
334  }
335 #endif
336 
337  // Store as best estimate value
338  startMassSafety = newMassSafety ;
339  startFullSafety = newFullSafety ;
340 
341  // Get the End-Position and End-Momentum (Dir-ection)
342  fTransportEndPosition = endTrackState.GetPosition() ;
343  fTransportEndKineticEnergy = endTrackState.GetKineticEnergy() ;
344  }
345  else
346  {
347  geometryStepLength = lengthAlongCurve= 0.0 ;
348  fMomentumChanged = false ;
349  // fMassGeometryLimitedStep = false ; // --- ???
350  // fAnyGeometryLimitedStep = true;
351  fTransportEndMomentumDir = track.GetMomentumDirection();
352  fTransportEndKineticEnergy = track.GetKineticEnergy();
353 
354  fTransportEndPosition = startPosition;
355  // If the step length requested is 0, and we are on a boundary
356  // then a boundary will also limit the step.
357  if( startMassSafety == 0.0 )
358  {
359  fMassGeometryLimitedStep = true ;
360  fAnyGeometryLimitedStep = true;
361  }
362  // TODO: Add explicit logical status for being at a boundary
363  }
364  // G4FieldTrack aTrackState(endTrackState);
365 
366  if( !fieldExertsForce )
367  {
368  fParticleIsLooping = false ;
369  fMomentumChanged = false ;
370  fEndGlobalTimeComputed = false ;
371  // G4cout << " global time is false " << G4endl;
372  }
373  else
374  {
375 
376 #ifdef G4DEBUG_TRANSPORT
377  if( fVerboseLevel > 1 )
378  {
379  G4cout << " G4CT::CS End Position = " << fTransportEndPosition << G4endl;
380  G4cout << " G4CT::CS End Direction = " << fTransportEndMomentumDir << G4endl;
381  }
382 #endif
383  if( fFieldPropagator->GetCurrentFieldManager()->DoesFieldChangeEnergy() )
384  {
385  // If the field can change energy, then the time must be integrated
386  // - so this should have been updated
387  //
388  fCandidateEndGlobalTime = endTrackState.GetLabTimeOfFlight();
389  fEndGlobalTimeComputed = true;
390 
391  // was ( fCandidateEndGlobalTime != track.GetGlobalTime() );
392  // a cleaner way is to have FieldTrack knowing whether time is updated.
393  }
394  else
395  {
396  // The energy should be unchanged by field transport,
397  // - so the time changed will be calculated elsewhere
398  //
399  fEndGlobalTimeComputed = false;
400 
401  // Check that the integration preserved the energy
402  // - and if not correct this!
403  G4double startEnergy= track.GetKineticEnergy();
404  G4double endEnergy= fTransportEndKineticEnergy;
405 
406  static G4int no_inexact_steps=0; // , no_large_ediff;
407  G4double absEdiff = std::fabs(startEnergy- endEnergy);
408  if( absEdiff > perMillion * endEnergy )
409  {
410  no_inexact_steps++;
411  // Possible statistics keeping here ...
412  }
413 #ifdef G4VERBOSE
414  if( (fVerboseLevel > 1) && ( absEdiff > perThousand * endEnergy) )
415  {
416  ReportInexactEnergy(startEnergy, endEnergy);
417  } // end of if (fVerboseLevel)
418 #endif
419  // Correct the energy for fields that conserve it
420  // This - hides the integration error
421  // - but gives a better physical answer
422  fTransportEndKineticEnergy= track.GetKineticEnergy();
423  }
424  }
425 
426  endpointDistance = (fTransportEndPosition - startPosition).mag() ;
427  fParticleIsLooping = fFieldPropagator->IsParticleLooping() ;
428 
429  fTransportEndSpin = endTrackState.GetSpin();
430 
431  // Calculate the safety
432  safetyProposal= startFullSafety; // used to be startMassSafety
433  // Changed to accomodate processes that cannot update the safety -- JA 22 Nov 06
434 
435  // Update safety for the end-point, if becomes negative at the end-point.
436 
437  if( (startFullSafety < endpointDistance )
438  && ( particleCharge != 0.0 ) ) // Only needed to prepare for Mult Scat.
439  // && !fAnyGeometryLimitedStep ) // To-Try: No safety update if at a boundary
440  {
441  G4double endFullSafety =
442  fPathFinder->ComputeSafety( fTransportEndPosition);
443  // Expected mission -- only mass geometry's safety
444  // fMassNavigator->ComputeSafety( fTransportEndPosition) ;
445  // Yet discrete processes only have poststep -- and this cannot
446  // currently revise the safety
447  // ==> so we use the all-geometry safety as a precaution
448 
449  fpSafetyHelper->SetCurrentSafety( endFullSafety, fTransportEndPosition);
450  // Pushing safety to Helper avoids recalculation at this point
451 
452  G4ThreeVector centerPt= G4ThreeVector(0.0, 0.0, 0.0); // Used for return value
453  G4double endMassSafety= fPathFinder->ObtainSafety( fNavigatorId, centerPt);
454  // Retrieves the mass value from PathFinder (it calculated it)
455 
456  fPreviousMassSafety = endMassSafety ;
457  fPreviousFullSafety = endFullSafety;
458  fPreviousSftOrigin = fTransportEndPosition ;
459 
460  // The convention (Stepping Manager's) is safety from the start point
461  //
462  safetyProposal = endFullSafety + endpointDistance;
463  // --> was endMassSafety
464  // Changed to accomodate processes that cannot update the safety -- JA 22 Nov 06
465 
466  // #define G4DEBUG_TRANSPORT 1
467 
468 #ifdef G4DEBUG_TRANSPORT
469  G4int prec= G4cout.precision(12) ;
470  G4cout << "***Transportation::AlongStepGPIL ** " << G4endl ;
471  G4cout << " Revised Safety at endpoint " << fTransportEndPosition
472  << " give safety values: Mass= " << endMassSafety
473  << " All= " << endFullSafety << G4endl ;
474  G4cout << " Adding endpoint distance " << endpointDistance
475  << " to obtain pseudo-safety= " << safetyProposal << G4endl ;
476  G4cout.precision(prec);
477  }
478  else
479  {
480  G4int prec= G4cout.precision(12) ;
481  G4cout << "***Transportation::AlongStepGPIL ** " << G4endl ;
482  G4cout << " Quick Safety estimate at endpoint " << fTransportEndPosition
483  << " gives safety endpoint value = " << startFullSafety - endpointDistance
484  << " using start-point value " << startFullSafety
485  << " and endpointDistance " << endpointDistance << G4endl;
486  G4cout.precision(prec);
487 #endif
488  }
489 
490  proposedSafetyForStart= safetyProposal;
491  fParticleChange.ProposeTrueStepLength(geometryStepLength) ;
492 
493  return geometryStepLength ;
494 }
495 
497 
500  const G4Step& stepData )
501 {
502  static G4int noCalls=0;
503  noCalls++;
504 
505  fParticleChange.Initialize(track) ;
506  // sets all its members to the value of corresponding members in G4Track
507 
508  // Code specific for Transport
509  //
510  fParticleChange.ProposePosition(fTransportEndPosition) ;
511  // G4cout << " G4CoupledTransportation::AlongStepDoIt"
512  // << " proposes position = " << fTransportEndPosition
513  // << " and end momentum direction = " << fTransportEndMomentumDir << G4endl;
514  fParticleChange.ProposeMomentumDirection(fTransportEndMomentumDir) ;
515  fParticleChange.ProposeEnergy(fTransportEndKineticEnergy) ;
516  fParticleChange.SetMomentumChanged(fMomentumChanged) ;
517 
518  fParticleChange.ProposePolarization(fTransportEndSpin);
519 
520  G4double deltaTime = 0.0 ;
521 
522  // Calculate Lab Time of Flight (ONLY if field Equations used it!)
523  // G4double endTime = fCandidateEndGlobalTime;
524  // G4double delta_time = endTime - startTime;
525 
526  G4double startTime = track.GetGlobalTime() ;
527 
528  if (!fEndGlobalTimeComputed)
529  {
530  G4double finalInverseVel= DBL_MAX, initialInverseVel=DBL_MAX;
531 
532  // The time was not integrated .. make the best estimate possible
533  //
534  G4double finalVelocity = track.GetVelocity() ;
535  if( finalVelocity > 0.0 ) { finalInverseVel= 1.0 / finalVelocity; }
536  G4double initialVelocity = stepData.GetPreStepPoint()->GetVelocity() ;
537  if( initialVelocity > 0.0 ) { initialInverseVel= 1.0 / initialVelocity; }
538  G4double stepLength = track.GetStepLength() ;
539 
540  if (finalVelocity > 0.0)
541  {
542  // deltaTime = stepLength/finalVelocity ;
543  G4double meanInverseVelocity = 0.5 * ( initialInverseVel + finalInverseVel );
544  deltaTime = stepLength * meanInverseVelocity ;
545  // G4cout << " dt = s * mean(1/v) , with " << " s = " << stepLength
546  // << " mean(1/v)= " << meanInverseVelocity << G4endl;
547  }
548  else
549  {
550  deltaTime = stepLength * initialInverseVel ;
551  // G4cout << " dt = s / initV " << " s = " << stepLength
552  // << " 1 / initV= " << initialInverseVel << G4endl;
553  } // Could do with better estimate for final step (finalVelocity = 0) ?
554 
555  fCandidateEndGlobalTime = startTime + deltaTime ;
556  fParticleChange.ProposeLocalTime( track.GetLocalTime() + deltaTime) ;
557 
558  // G4cout << " Calculated global time from start = " << startTime << " and "
559  // << " delta time = " << deltaTime << G4endl;
560  }
561  else
562  {
563  deltaTime = fCandidateEndGlobalTime - startTime ;
564  fParticleChange.ProposeGlobalTime( fCandidateEndGlobalTime ) ;
565  // G4cout << " Calculated global time from candidate end time = "
566  // << fCandidateEndGlobalTime << " and start time = " << startTime << G4endl;
567  }
568 
569  // G4cout << " G4CoupledTransportation::AlongStepDoIt "
570  // << " flag whether computed time = " << fEndGlobalTimeComputed << " and "
571  // << " is proposes end time " << fCandidateEndGlobalTime << G4endl;
572 
573  // Now Correct by Lorentz factor to get "proper" deltaTime
574 
575  G4double restMass = track.GetDynamicParticle()->GetMass() ;
576  G4double deltaProperTime = deltaTime*( restMass/track.GetTotalEnergy() ) ;
577 
578  fParticleChange.ProposeProperTime(track.GetProperTime() + deltaProperTime) ;
579  //fParticleChange. ProposeTrueStepLength( track.GetStepLength() ) ;
580 
581  // If the particle is caught looping or is stuck (in very difficult
582  // boundaries) in a magnetic field (doing many steps)
583  // THEN this kills it ...
584  //
585  if ( fParticleIsLooping )
586  {
587  G4double endEnergy= fTransportEndKineticEnergy;
588 
589  if( (endEnergy < fThreshold_Important_Energy)
590  || (fNoLooperTrials >= fThresholdTrials ) )
591  {
592  // Kill the looping particle
593  //
594  fParticleChange.ProposeTrackStatus( fStopAndKill ) ;
595 
596  // 'Bare' statistics
597  fSumEnergyKilled += endEnergy;
598  if( endEnergy > fMaxEnergyKilled) { fMaxEnergyKilled= endEnergy; }
599 
600 #ifdef G4VERBOSE
601  if((fVerboseLevel > 1) && ( endEnergy > fThreshold_Warning_Energy ))
602  {
603  G4cout << " G4CoupledTransportation is killing track that is looping or stuck " << G4endl
604  << " This track has " << track.GetKineticEnergy() / MeV
605  << " MeV energy." << G4endl;
606  }
607  if( fVerboseLevel > 0 )
608  {
609  G4cout << " Steps by this track: " << track.GetCurrentStepNumber() << G4endl;
610  }
611 #endif
612  fNoLooperTrials=0;
613  }
614  else
615  {
616  fNoLooperTrials ++;
617 #ifdef G4VERBOSE
618  if( (fVerboseLevel > 2) )
619  {
620  G4cout << " ** G4CoupledTransportation::AlongStepDoIt(): Particle looping - " << G4endl
621  << " Number of consecutive problem step (this track) = " << fNoLooperTrials << G4endl
622  << " Steps by this track: " << track.GetCurrentStepNumber() << G4endl
623  << " Total no of calls to this method (all tracks) = " << noCalls << G4endl;
624  }
625 #endif
626  }
627  }
628  else
629  {
630  fNoLooperTrials=0;
631  }
632 
633  // Another (sometimes better way) is to use a user-limit maximum Step size
634  // to alleviate this problem ..
635 
636  // Add smooth curved trajectories to particle-change
637  //
638  // fParticleChange.SetPointerToVectorOfAuxiliaryPoints
639  // (fFieldPropagator->GimmeTrajectoryVectorAndForgetIt() );
640 
641  return &fParticleChange ;
642 }
643 
645 //
646 // This ensures that the PostStep action is always called,
647 // so that it can do the relocation if it is needed.
648 //
649 
652  G4double, // previousStepSize
653  G4ForceCondition* pForceCond )
654 {
655  // Must act as PostStep action -- to relocate particle
656  *pForceCond = Forced ;
657  return DBL_MAX ;
658 }
659 
661 ReportMove( G4ThreeVector OldVector, G4ThreeVector NewVector, const G4String& Quantity )
662 {
663  G4ThreeVector moveVec = ( NewVector - OldVector );
664 
665  G4cerr << G4endl
666  << "**************************************************************" << G4endl;
667  G4cerr << "Endpoint has moved between value expected from TransportEndPosition "
668  << " and value from Track in PostStepDoIt. " << G4endl
669  << "Change of " << Quantity << " is " << moveVec.mag() / mm << " mm long, "
670  << " and its vector is " << (1.0/mm) * moveVec << " mm " << G4endl
671  << "Endpoint of ComputeStep was " << OldVector
672  << " and current position to locate is " << NewVector << G4endl;
673 }
674 
676 
678  const G4Step& )
679 {
680  G4TouchableHandle retCurrentTouchable ; // The one to return
681 
682  // Initialize ParticleChange (by setting all its members equal
683  // to corresponding members in G4Track)
684  // fParticleChange.Initialize(track) ; // To initialise TouchableChange
685 
686  fParticleChange.ProposeTrackStatus(track.GetTrackStatus()) ;
687 
688  // Check that the end position and direction are preserved
689  // since call to AlongStepDoIt
690 
691 #ifdef G4DEBUG_TRANSPORT
692  if( ( fVerboseLevel > 0 ) && ((fTransportEndPosition - track.GetPosition()).mag2() >= 1.0e-16) )
693  {
694  ReportMove( track.GetPosition(), fTransportEndPosition, "End of Step Position" );
695  G4cerr << " Problem in G4CoupledTransportation::PostStepDoIt " << G4endl;
696  }
697 
698  // If the Step was determined by the volume boundary, relocate the particle
699  // The pathFinder will know that the geometry limited the step (!?)
700 
701  if( fVerboseLevel > 0 )
702  {
703  G4cout << " Calling PathFinder::Locate() from "
704  << " G4CoupledTransportation::PostStepDoIt " << G4endl;
705  G4cout << " fAnyGeometryLimitedStep is " << fAnyGeometryLimitedStep << G4endl;
706 
707  }
708 #endif
709 
710  if(fAnyGeometryLimitedStep)
711  {
712  fPathFinder->Locate( track.GetPosition(),
713  track.GetMomentumDirection(),
714  true);
715 
716  // fCurrentTouchable will now become the previous touchable,
717  // and what was the previous will be freed.
718  // (Needed because the preStepPoint can point to the previous touchable)
719 
720  fCurrentTouchableHandle=
721  fPathFinder->CreateTouchableHandle( fNavigatorId );
722 
723 #ifdef G4DEBUG_TRANSPORT
724  if( fVerboseLevel > 0 )
725  {
726  G4cout << "G4CoupledTransportation::PostStepDoIt --- fNavigatorId = "
727  << fNavigatorId << G4endl;
728  }
729  if( fVerboseLevel > 1 )
730  {
731  G4VPhysicalVolume* vol= fCurrentTouchableHandle->GetVolume();
732  G4cout << "CHECK !!!!!!!!!!! fCurrentTouchableHandle->GetVolume() = " << vol;
733  if( vol ) { G4cout << "Name=" << vol->GetName(); }
734  G4cout << G4endl;
735  }
736 #endif
737 
738  // Check whether the particle is out of the world volume
739  // If so it has exited and must be killed.
740  //
741  if( fCurrentTouchableHandle->GetVolume() == 0 )
742  {
743  fParticleChange.ProposeTrackStatus( fStopAndKill ) ;
744  }
745  retCurrentTouchable = fCurrentTouchableHandle ;
746  // fParticleChange.SetTouchableHandle( fCurrentTouchableHandle ) ;
747 
748  // Notify particle change that this is last step in volume
749  fParticleChange.ProposeLastStepInVolume(true);
750  }
751  else // fAnyGeometryLimitedStep is false
752  {
753 #ifdef G4DEBUG_TRANSPORT
754  if( fVerboseLevel > 1 )
755  {
756  G4cout << "G4CoupledTransportation::PostStepDoIt -- "
757  << " fAnyGeometryLimitedStep = " << fAnyGeometryLimitedStep
758  << " must be false " << G4endl;
759  }
760 #endif
761  // This serves only to move each of the Navigator's location
762  //
763  // fLinearNavigator->LocateGlobalPointWithinVolume( track.GetPosition() ) ;
764 
765  // G4cout << "G4CoupledTransportation calling PathFinder::ReLocate() " << G4endl;
766  fPathFinder->ReLocate( track.GetPosition() );
767  // track.GetMomentumDirection() );
768 
769  // Keep the value of the track's current Touchable is retained,
770  // and use it to overwrite the (unset) one in particle change.
771  // Expect this must be fCurrentTouchable too
772  // - could it be different, eg at the start of a step ?
773  //
774  retCurrentTouchable = track.GetTouchableHandle() ;
775  // fParticleChange.SetTouchableHandle( track.GetTouchableHandle() ) ;
776 
777  // Have not reached a boundary
778  fParticleChange.ProposeLastStepInVolume(false);
779  } // endif ( fAnyGeometryLimitedStep )
780 
781  const G4VPhysicalVolume* pNewVol = retCurrentTouchable->GetVolume() ;
782  const G4Material* pNewMaterial = 0 ;
783  const G4VSensitiveDetector* pNewSensitiveDetector = 0 ;
784 
785  if( pNewVol != 0 )
786  {
787  pNewMaterial= pNewVol->GetLogicalVolume()->GetMaterial();
788  pNewSensitiveDetector= pNewVol->GetLogicalVolume()->GetSensitiveDetector();
789  }
790 
791  // ( const_cast<G4Material *> pNewMaterial ) ;
792  // ( const_cast<G4VSensitiveDetetor *> pNewSensitiveDetector) ;
793 
794  fParticleChange.SetMaterialInTouchable( (G4Material *) pNewMaterial ) ;
795  fParticleChange.SetSensitiveDetectorInTouchable( (G4VSensitiveDetector *) pNewSensitiveDetector ) ;
796  // "temporarily" until Get/Set Material of ParticleChange,
797  // and StepPoint can be made const.
798 
799  const G4MaterialCutsCouple* pNewMaterialCutsCouple = 0;
800  if( pNewVol != 0 )
801  {
802  pNewMaterialCutsCouple=pNewVol->GetLogicalVolume()->GetMaterialCutsCouple();
803  if( pNewMaterialCutsCouple!=0
804  && pNewMaterialCutsCouple->GetMaterial()!=pNewMaterial )
805  {
806  // for parametrized volume
807  //
808  pNewMaterialCutsCouple =
810  ->GetMaterialCutsCouple(pNewMaterial,
811  pNewMaterialCutsCouple->GetProductionCuts());
812  }
813  }
814  fParticleChange.SetMaterialCutsCoupleInTouchable( pNewMaterialCutsCouple );
815 
816  // Must always set the touchable in ParticleChange, whether relocated or not
817  fParticleChange.SetTouchableHandle(retCurrentTouchable) ;
818 
819  return &fParticleChange ;
820 }
821 
822 // New method takes over the responsibility to reset the state of
823 // G4CoupledTransportation object:
824 // - at the start of a new track, and
825 // - on the resumption of a suspended track.
826 
827 void
829 {
830 
831  G4TransportationManager* transportMgr =
833 
834  // G4VProcess::StartTracking(aTrack);
835 
836  // The 'initialising' actions
837  // once taken in AlongStepGPIL -- if ( track.GetCurrentStepNumber()==1 )
838 
839  // fStartedNewTrack= true;
840 
841  fMassNavigator = transportMgr->GetNavigatorForTracking() ;
842  fNavigatorId= transportMgr->ActivateNavigator( fMassNavigator ); // Confirm it!
843 
844  // if( fVerboseLevel > 1 ){
845  // G4cout << " Navigator Id obtained in StartTracking " << fNavigatorId << G4endl;
846  // }
847  G4ThreeVector position = aTrack->GetPosition();
848  G4ThreeVector direction = aTrack->GetMomentumDirection();
849 
850  // if( fVerboseLevel > 1 ){
851  // G4cout << " Calling PathFinder::PrepareNewTrack from "
852  // << " G4CoupledTransportation::StartTracking -- which calls Locate()" << G4endl;
853  // }
854  fPathFinder->PrepareNewTrack( position, direction);
855  // This implies a call to fPathFinder->Locate( position, direction );
856 
857  // Global field, if any, must exist before tracking is started
858  fGlobalFieldExists= DoesGlobalFieldExist();
859  // reset safety value and center
860  //
861  fPreviousMassSafety = 0.0 ;
862  fPreviousFullSafety = 0.0 ;
863  fPreviousSftOrigin = G4ThreeVector(0.,0.,0.) ;
864 
865  // reset looping counter -- for motion in field
866  fNoLooperTrials= 0;
867  // Must clear this state .. else it depends on last track's value
868  // --> a better solution would set this from state of suspended track TODO ?
869  // Was if( aTrack->GetCurrentStepNumber()==1 ) { .. }
870 
871  // ChordFinder reset internal state
872  //
873  if( fGlobalFieldExists )
874  {
875  fFieldPropagator->ClearPropagatorState();
876  // Resets safety values, in case of overlaps.
877 
878  G4ChordFinder* chordF= fFieldPropagator->GetChordFinder();
879  if( chordF ) { chordF->ResetStepEstimate(); }
880  }
881 
882  // Clear the chord finders of all fields (ie managers) derived objects
883  //
885  fieldMgrStore->ClearAllChordFindersState();
886 
887 #ifdef G4DEBUG_TRANSPORT
888  if( fVerboseLevel > 1 )
889  {
890  G4cout << " Returning touchable handle " << fCurrentTouchableHandle << G4endl;
891  }
892 #endif
893 
894  // Update the current touchable handle (from the track's)
895  //
896  fCurrentTouchableHandle = aTrack->GetTouchableHandle();
897 }
898 
899 void
901 {
903 }
904 
905 void
907 ReportInexactEnergy(G4double startEnergy, G4double endEnergy)
908 {
909  static G4int no_warnings= 0, warnModulo=1, moduloFactor= 10, no_large_ediff= 0;
910 
911  if( std::fabs(startEnergy- endEnergy) > perThousand * endEnergy )
912  {
913  no_large_ediff ++;
914  if( (no_large_ediff% warnModulo) == 0 )
915  {
916  no_warnings++;
917  G4cout << "WARNING - G4CoupledTransportation::AlongStepGetPIL() "
918  << " Energy change in Step is above 1^-3 relative value. " << G4endl
919  << " Relative change in 'tracking' step = "
920  << std::setw(15) << (endEnergy-startEnergy)/startEnergy << G4endl
921  << " Starting E= " << std::setw(12) << startEnergy / MeV << " MeV " << G4endl
922  << " Ending E= " << std::setw(12) << endEnergy / MeV << " MeV " << G4endl;
923  G4cout << " Energy has been corrected -- however, review"
924  << " field propagation parameters for accuracy." << G4endl;
925  if( (fVerboseLevel > 2 ) || (no_warnings<4) || (no_large_ediff == warnModulo * moduloFactor) )
926  {
927  G4cout << " These include EpsilonStepMax(/Min) in G4FieldManager "
928  << " which determine fractional error per step for integrated quantities. " << G4endl
929  << " Note also the influence of the permitted number of integration steps."
930  << G4endl;
931  }
932  G4cerr << "ERROR - G4CoupledTransportation::AlongStepGetPIL()" << G4endl
933  << " Bad 'endpoint'. Energy change detected"
934  << " and corrected. "
935  << " Has occurred already "
936  << no_large_ediff << " times." << G4endl;
937  if( no_large_ediff == warnModulo * moduloFactor )
938  {
939  warnModulo *= moduloFactor;
940  }
941  }
942  }
943 }