78#define PrepareState() G4ITTransportationState* __state = this->GetState<G4ITTransportationState>();
82#define State(theXInfo) (__state->theXInfo)
97 fThreshold_Warning_Energy(100 * MeV),
98 fThreshold_Important_Energy(250 * MeV),
100 fUnimportant_Energy(1 * MeV),
101 fSumEnergyKilled(0.0),
102 fMaxEnergyKilled(0.0),
103 fShortStepOptimisation(false),
104 fVerboseLevel(verbose)
126 fInstantiateProcessState =
true;
172 fInstantiateProcessState = right.fInstantiateProcessState;
219 G4cout <<
" G4ITTransportation: Statistics for looping particles "
221 G4cout <<
" Sum of energy of loopers killed: "
223 G4cout <<
" Max energy of loopers killed: "
259 G4double geometryStepLength(-1.0), newSafety(-1.0);
261 State(fParticleIsLooping) =
false;
262 State(fEndGlobalTimeComputed) =
false;
263 State(fGeometryLimitedStep) =
false;
311 G4bool fieldExertsForce =
false;
312 if ((particleCharge != 0.0))
333 if (!fieldExertsForce)
340 geometryStepLength = currentMinimumStep;
341 State(fGeometryLimitedStep) =
false;
367 State(fTransportEndPosition));
372 currentSafety = newSafety;
374 State(fGeometryLimitedStep) = (linearStepLength <= currentMinimumStep);
375 if (
State(fGeometryLimitedStep))
378 geometryStepLength = linearStepLength;
383 geometryStepLength = currentMinimumStep;
386 State(fEndPointDistance) = geometryStepLength;
390 State(fTransportEndPosition) = startPosition
391 + geometryStepLength * startMomentumDir;
395 State(fTransportEndMomentumDir) = startMomentumDir;
398 State(fParticleIsLooping) =
false;
399 State(fMomentumChanged) =
false;
400 State(fEndGlobalTimeComputed) =
true;
401 State(theInteractionTimeLeft) =
State(fEndPointDistance)
403 State(fCandidateEndGlobalTime) =
State(theInteractionTimeLeft)
421 <<
"ITTransportation does not support external fields.";
423 <<
" If you are dealing with a tradiational MC simulation, ";
424 exceptionDescription <<
"please use G4Transportation.";
426 G4Exception(
"G4ITTransportation::AlongStepGetPhysicalInteractionLength",
582 if (currentMinimumStep == 0.0)
584 if (currentSafety == 0.0)
586 State(fGeometryLimitedStep) =
true;
596 if (currentSafety <
State(fEndPointDistance))
601 if (particleCharge != 0.0)
605 State(fTransportEndPosition));
606 currentSafety = endSafety;
619 State(fTransportEndPosition));
625 currentSafety +=
State(fEndPointDistance);
627#ifdef G4DEBUG_TRANSPORT
629 G4cout <<
"***G4Transportation::AlongStepGPIL ** " <<
G4endl;
630 G4cout <<
" Called Navigator->ComputeSafety at "
631 <<
State(fTransportEndPosition)
632 <<
" and it returned safety= " << endSafety <<
G4endl;
633 G4cout <<
" Adding endpoint distance " <<
State(fEndPointDistance)
634 <<
" to obtain pseudo-safety= " << currentSafety <<
G4endl;
644 return geometryStepLength;
649 const double timeStep,
650 double& oPhysicalStep)
660 State(fGeometryLimitedStep) =
false;
666 State(fEndGlobalTimeComputed) =
true;
670 if (!
State(fMomentumChanged))
674 oPhysicalStep = initialVelocity * timeStep;
678 State(fTransportEndPosition) = startPosition
679 + oPhysicalStep * startMomentumDir;
692#if defined (DEBUG_MEM)
693 MemStat mem_first, mem_second, mem_diff;
696#if defined (DEBUG_MEM)
697 mem_first = MemoryUsage();
708 if (!pdefOpticalPhoton) pdefOpticalPhoton =
735 if (
State(fEndGlobalTimeComputed) ==
false)
748 deltaTime = stepLength / finalVelocity;
750 else if (initialVelocity > 0.0)
752 deltaTime = stepLength / initialVelocity;
755 State(fCandidateEndGlobalTime) = startTime + deltaTime;
759 deltaTime =
State(fCandidateEndGlobalTime) - startTime;
782 if (
State(fParticleIsLooping))
805 <<
" G4ITTransportation is killing track that is looping or stuck "
807 <<
" MeV energy." <<
G4endl;
808 G4cout <<
" Number of trials = " <<
State(fNoLooperTrials)
809 <<
" No of calls to AlongStepDoIt = " << noCalls
813 State(fNoLooperTrials) = 0;
817 State(fNoLooperTrials)++;
821 G4cout <<
" G4ITTransportation::AlongStepDoIt(): Particle looping - "
822 <<
" Number of trials = " <<
State(fNoLooperTrials)
823 <<
" No of calls to = " << noCalls <<
G4endl;
830 State(fNoLooperTrials)=0;
841#if defined (DEBUG_MEM)
842 mem_second = MemoryUsage();
843 mem_diff = mem_second-mem_first;
844 G4cout <<
"\t || MEM || End of G4ITTransportation::AlongStepDoIt, diff is: "
877 G4bool isLastStep =
false;
888 if (
State(fGeometryLimitedStep))
893 G4cout <<
"Step is limited by geometry "
901 if (
State(fCurrentTouchableHandle)->GetVolume() == 0)
904 exceptionDescription <<
"No current touchable found ";
905 G4Exception(
" G4ITTransportation::PostStepDoIt",
"G4ITTransportation001",
912 State(fCurrentTouchableHandle),
true);
916 if (
State(fCurrentTouchableHandle)->GetVolume() == 0)
924 G4cout <<
"G4ITTransportation will killed the track because "
925 "State(fCurrentTouchableHandle)->GetVolume() == 0"<<
G4endl;
931 retCurrentTouchable =
State(fCurrentTouchableHandle);
947#ifdef G4DEBUG_TRANSPORT
952 if( ! (exiting || entering) )
954 G4cout <<
" Transport> : Proposed isLastStep= " << isLastStep
957 <<
" Track position : " << track.
GetPosition() /nanometer <<
" [nm]"
980#ifdef G4DEBUG_TRANSPORT
982 G4cout <<
" Transport> Proposed isLastStep= " << isLastStep
983 <<
" Geometry did not limit step. Position : "
1010 pNewMaterialCutsCouple =
1014 if (pNewVol != 0 && pNewMaterialCutsCouple != 0
1015 && pNewMaterialCutsCouple->
GetMaterial() != pNewMaterial)
1043 if (fInstantiateProcessState)
1052 GetIT(track)->GetTrackingInfo()->GetTrackStateManager());
1083 fieldMgrStore->ClearAllChordFindersState();
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
#define fPreviousSftOrigin
G4IT * GetIT(const G4Track *track)
CLHEP::Hep3Vector G4ThreeVector
G4ReferenceCountedHandle< G4VTouchable > G4TouchableHandle
G4GLOB_DLL std::ostream G4cout
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
G4double fMaxEnergyKilled
G4ParticleChangeForTransport fParticleChange
G4double fUnimportant_Energy
G4double fThreshold_Important_Energy
G4double fSumEnergyKilled
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
G4bool DoesGlobalFieldExist()
virtual G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData)
virtual ~G4ITTransportation()
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double, G4double currentMinimumStep, G4double ¤tSafety, G4GPILSelection *selection)
G4ITSafetyHelper * fpSafetyHelper
G4ITTransportation(const G4String &aName="ITTransportation", G4int verbosityLevel=0)
G4bool fShortStepOptimisation
G4double fThreshold_Warning_Energy
G4ITNavigator * fLinearNavigator
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &)
virtual void StartTracking(G4Track *aTrack)
G4TrackingInformation * GetTrackingInfo()
G4VSensitiveDetector * GetSensitiveDetector() const
G4Material * GetMaterial() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
const G4Material * GetMaterial() const
G4ProductionCuts * GetProductionCuts() const
void SetPointerToVectorOfAuxiliaryPoints(std::vector< G4ThreeVector > *vec)
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 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)
void ClearPropagatorState()
std::vector< G4ThreeVector > * GimmeTrajectoryVectorAndForgetIt() const
G4double GetVelocity() const
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
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
static G4TransportationManager * GetTransportationManager()
G4PropagatorInField * GetPropagatorInField() const
G4FieldManager * GetFieldManager() const
virtual void StartTracking(G4Track *)
void SetInstantiateProcessState(G4bool flag)
G4shared_ptr< G4ProcessState > fpState
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLastStepInVolume(G4bool flag)
void ProposeTrueStepLength(G4double truePathLength)
G4LogicalVolume * GetLogicalVolume() const
G4bool enableAlongStepDoIt
void SetProcessSubType(G4int)
virtual void StartTracking(G4Track *)
G4bool enablePostStepDoIt
G4VParticleChange * pParticleChange
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
G4TouchableHandle fCurrentTouchableHandle
G4bool fParticleIsLooping
virtual ~G4ITTransportationState()
G4ThreeVector fPreviousSftOrigin
G4ThreeVector fTransportEndMomentumDir
G4ITTransportationState()
Process State.
G4double fEndPointDistance
G4double fCandidateEndGlobalTime
G4ThreeVector fTransportEndPosition
G4double fTransportEndKineticEnergy
G4bool fEndGlobalTimeComputed
G4bool fGeometryLimitedStep
G4ThreeVector fTransportEndSpin