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