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