Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ITTransportation Class Reference

#include <G4ITTransportation.hh>

+ Inheritance diagram for G4ITTransportation:

Classes

struct  G4ITTransportationState
 

Public Member Functions

 G4ITTransportation (const G4String &aName="ITTransportation", G4int verbosityLevel=1)
 
virtual ~G4ITTransportation ()
 
 G4ITTransportation (const G4ITTransportation &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void ComputeStep (const G4Track &, const G4Step &, const double timeStep, double &spaceStep)
 
virtual void StartTracking (G4Track *aTrack)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &, G4double, G4ForceCondition *pForceCond)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &)
 
G4PropagatorInFieldGetPropagatorInField ()
 
void SetPropagatorInField (G4PropagatorInField *pFieldPropagator)
 
void SetVerboseLevel (G4int verboseLevel)
 
G4int GetVerboseLevel () const
 
G4double GetThresholdWarningEnergy () const
 
G4double GetThresholdImportantEnergy () const
 
G4int GetThresholdTrials () const
 
void SetThresholdWarningEnergy (G4double newEnWarn)
 
void SetThresholdImportantEnergy (G4double newEnImp)
 
void SetThresholdTrials (G4int newMaxTrials)
 
G4double GetMaxEnergyKilled () const
 
G4double GetSumEnergyKilled () const
 
void ResetKilledStatistics (G4int report=1)
 
void EnableShortStepOptimisation (G4bool optimise=true)
 
- Public Member Functions inherited from G4VITProcess
 G4VITProcess (const G4String &name, G4ProcessType type=fNotDefined)
 
virtual ~G4VITProcess ()
 
 G4VITProcess (const G4VITProcess &other)
 
G4VITProcessoperator= (const G4VITProcess &other)
 
G4int operator== (const G4VITProcess &right) const
 
G4int operator!= (const G4VITProcess &right) const
 
size_t GetProcessID () const
 
G4ProcessState_LockGetProcessState ()
 
void SetProcessState (G4ProcessState_Lock *aProcInfo)
 
virtual void StartTracking (G4Track *)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
G4double GetInteractionTimeLeft ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4bool ProposesTimeStep () const
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4int operator== (const G4VProcess &right) const
 
G4int operator!= (const G4VProcess &right) const
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual void StartTracking (G4Track *)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 

Protected Member Functions

G4bool DoesGlobalFieldExist ()
 
void SetInstantiateProcessState (G4bool flag)
 
G4bool InstantiateProcessState ()
 
- Protected Member Functions inherited from G4VITProcess
void RetrieveProcessInfo ()
 
void CreateInfo ()
 
virtual void ClearInteractionTimeLeft ()
 
virtual void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
virtual void ClearNumberOfInteractionLengthLeft ()
 
void SetInstantiateProcessState (G4bool flag)
 
G4bool InstantiateProcessState ()
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Protected Attributes

G4ITTransportationState *constfTransportationState
 
G4ITNavigatorfLinearNavigator
 
G4PropagatorInFieldfFieldPropagator
 
G4ParticleChangeForTransport fParticleChange
 
G4double fThreshold_Warning_Energy
 
G4double fThreshold_Important_Energy
 
G4int fThresholdTrials
 
G4double fUnimportant_Energy
 
G4double fSumEnergyKilled
 
G4double fMaxEnergyKilled
 
G4bool fShortStepOptimisation
 
G4SafetyHelperfpSafetyHelper
 
G4int fVerboseLevel
 
- Protected Attributes inherited from G4VITProcess
G4ProcessStatefpState
 
G4bool fProposesTimeStep
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChangepParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VITProcess
static const size_t & GetMaxProcessIndex ()
 
- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 

Detailed Description

Definition at line 57 of file G4ITTransportation.hh.

Constructor & Destructor Documentation

◆ G4ITTransportation() [1/2]

G4ITTransportation::G4ITTransportation ( const G4String aName = "ITTransportation",
G4int  verbosityLevel = 1 
)

Definition at line 81 of file G4ITTransportation.cc.

81 :
84 fThreshold_Warning_Energy( 100 * MeV ),
85 fThreshold_Important_Energy( 250 * MeV ),
86 fThresholdTrials( 10 ),
87 fUnimportant_Energy( 1 * MeV ), // Not used
89 fShortStepOptimisation(false), // Old default: true (=fast short steps)
90 fVerboseLevel( verbose )
91{
93 G4TransportationManager* transportMgr ;
95 G4ITTransportationManager* ITtransportMgr ;
97 fLinearNavigator = ITtransportMgr->GetNavigatorForTracking() ;
98 fFieldPropagator = transportMgr->GetPropagatorInField() ;
99 fpSafetyHelper = 0; // transportMgr->GetSafetyHelper(); // New
100
101 // Cannot determine whether a field exists here, as it would
102 // depend on the relative order of creating the detector's
103 // field and this process. That order is not guaranted.
104 // Instead later the method DoesGlobalFieldExist() is called
105
106 enableAtRestDoIt = false;
107 enableAlongStepDoIt = true;
108 enablePostStepDoIt = true;
112 fInstantiateProcessState = true;
113}
@ fTransportation
#define InitProcessState(destination, source)
Definition: G4VITProcess.hh:53
static G4ITTransportationManager * GetTransportationManager()
G4ParticleChangeForTransport fParticleChange
G4double fThreshold_Important_Energy
G4SafetyHelper * fpSafetyHelper
void SetInstantiateProcessState(G4bool flag)
G4PropagatorInField * fFieldPropagator
G4ITTransportationState *const & fTransportationState
G4ITNavigator * fLinearNavigator
static G4TransportationManager * GetTransportationManager()
G4PropagatorInField * GetPropagatorInField() const
void SetInstantiateProcessState(G4bool flag)
G4ProcessState * fpState
G4bool enableAtRestDoIt
Definition: G4VProcess.hh:350
G4bool enableAlongStepDoIt
Definition: G4VProcess.hh:351
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:403
G4bool enablePostStepDoIt
Definition: G4VProcess.hh:352
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283

◆ ~G4ITTransportation()

G4ITTransportation::~G4ITTransportation ( )
virtual

Definition at line 189 of file G4ITTransportation.cc.

190{
191#ifdef G4VERBOSE
192 if( (fVerboseLevel > 0) && (fSumEnergyKilled > 0.0 ) )
193 {
194 G4cout << " G4ITTransportation: Statistics for looping particles " << G4endl;
195 G4cout << " Sum of energy of loopers killed: " << fSumEnergyKilled << G4endl;
196 G4cout << " Max energy of loopers killed: " << fMaxEnergyKilled << G4endl;
197 }
198#endif
199}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout

◆ G4ITTransportation() [2/2]

G4ITTransportation::G4ITTransportation ( const G4ITTransportation right)

Definition at line 116 of file G4ITTransportation.cc.

116 :
117 G4VITProcess(right),
119{
120 // Copy attributes
129
130 // Setup Navigators
131 G4TransportationManager* transportMgr ;
133 G4ITTransportationManager* ITtransportMgr ;
135 fLinearNavigator = ITtransportMgr->GetNavigatorForTracking() ;
136 fFieldPropagator = transportMgr->GetPropagatorInField() ;
137 fpSafetyHelper = 0; //transportMgr->GetSafetyHelper(); // New
138
139 // Cannot determine whether a field exists here, as it would
140 // depend on the relative order of creating the detector's
141 // field and this process. That order is not guaranted.
142 // Instead later the method DoesGlobalFieldExist() is called
143
144 enableAtRestDoIt = false;
145 enableAlongStepDoIt = true;
146 enablePostStepDoIt = true;
147
151 fInstantiateProcessState = right.fInstantiateProcessState;
152}

Member Function Documentation

◆ AlongStepDoIt()

G4VParticleChange * G4ITTransportation::AlongStepDoIt ( const G4Track track,
const G4Step stepData 
)
virtual

!!!!!!! A REVOIR !!!!


Implements G4VProcess.

Reimplemented in G4DNABrownianTransportation.

Definition at line 610 of file G4ITTransportation.cc.

612{
613
614 // G4cout << "G4ITTransportation::AlongStepDoIt" << G4endl;
615 // set pdefOpticalPhoton
616 static G4ParticleDefinition* pdefOpticalPhoton =
618
619 static G4int noCalls=0;
620 noCalls++;
621
623
624 // Code for specific process
625 //
626 fParticleChange.ProposePosition(State(fTransportEndPosition)) ;
627 fParticleChange.ProposeMomentumDirection(State(fTransportEndMomentumDir)) ;
628 fParticleChange.ProposeEnergy(State(fTransportEndKineticEnergy)) ;
629 fParticleChange.SetMomentumChanged(State(fMomentumChanged)) ;
630
631 fParticleChange.ProposePolarization(State(fTransportEndSpin));
632
633 G4double deltaTime = 0.0 ;
634
635 // Calculate Lab Time of Flight (ONLY if field Equations used it!)
636 // G4double endTime = State(fCandidateEndGlobalTime);
637 // G4double delta_time = endTime - startTime;
638
639 G4double startTime = track.GetGlobalTime() ;
640 ///___________________________________________________________________________
641 /// !!!!!!!
642 /// A REVOIR !!!!
643 if (State(fEndGlobalTimeComputed) == false)
644 {
645 // The time was not integrated .. make the best estimate possible
646 //
647 G4double initialVelocity = stepData.GetPreStepPoint()->GetVelocity() ;
648 G4double stepLength = track.GetStepLength() ;
649
650 deltaTime= 0.0; // in case initialVelocity = 0
651 if (track.GetParticleDefinition() == pdefOpticalPhoton)
652 {
653 // For only Optical Photon, final velocity is used
654 double finalVelocity = track.CalculateVelocityForOpticalPhoton();
655 fParticleChange.ProposeVelocity(finalVelocity);
656 deltaTime = stepLength/finalVelocity ;
657 }
658 else if( initialVelocity > 0.0 )
659 {
660 deltaTime = stepLength/initialVelocity ;
661 }
662
663 State(fCandidateEndGlobalTime) = startTime + deltaTime ;
664 }
665 else
666 {
667 deltaTime = State(fCandidateEndGlobalTime) - startTime ;
668 }
669
670 fParticleChange.ProposeGlobalTime( State(fCandidateEndGlobalTime) ) ;
671 fParticleChange.ProposeLocalTime( track.GetLocalTime() + deltaTime) ;
672 /*
673 // Now Correct by Lorentz factor to get delta "proper" Time
674
675 G4double restMass = track.GetDynamicParticle()->GetMass() ;
676 G4double deltaProperTime = deltaTime*( restMass/track.GetTotalEnergy() ) ;
677
678 fParticleChange.ProposeProperTime(track.GetProperTime() + deltaProperTime) ;
679*/
680
681 fParticleChange. ProposeTrueStepLength( track.GetStepLength() ) ;
682
683 ///___________________________________________________________________________
684 ///
685
686 // If the particle is caught looping or is stuck (in very difficult
687 // boundaries) in a magnetic field (doing many steps)
688 // THEN this kills it ...
689 //
690 if ( State(fParticleIsLooping) )
691 {
692 G4double endEnergy= State(fTransportEndKineticEnergy);
693
694 if( (endEnergy < fThreshold_Important_Energy)
695 || (State(fNoLooperTrials) >= fThresholdTrials ) )
696 {
697 // Kill the looping particle
698 //
699 // G4cout << "G4ITTransportation will killed the molecule"<< G4endl;
701
702 // 'Bare' statistics
703 fSumEnergyKilled += endEnergy;
704 if( endEnergy > fMaxEnergyKilled)
705 {
706 fMaxEnergyKilled= endEnergy;
707 }
708
709#ifdef G4VERBOSE
710 if( (fVerboseLevel > 1) ||
711 ( endEnergy > fThreshold_Warning_Energy ) )
712 {
713 G4cout << " G4ITTransportation is killing track that is looping or stuck "
714 << G4endl
715 << " This track has " << track.GetKineticEnergy() / MeV
716 << " MeV energy." << G4endl;
717 G4cout << " Number of trials = " << State(fNoLooperTrials)
718 << " No of calls to AlongStepDoIt = " << noCalls
719 << G4endl;
720 }
721#endif
722 State(fNoLooperTrials)=0;
723 }
724 else
725 {
726 State(fNoLooperTrials) ++;
727#ifdef G4VERBOSE
728 if( (fVerboseLevel > 2) )
729 {
730 G4cout << " G4ITTransportation::AlongStepDoIt(): Particle looping - "
731 << " Number of trials = " << State(fNoLooperTrials)
732 << " No of calls to = " << noCalls
733 << G4endl;
734 }
735#endif
736 }
737 }
738 else
739 {
740 State(fNoLooperTrials)=0;
741 }
742
743 // Another (sometimes better way) is to use a user-limit maximum Step size
744 // to alleviate this problem ..
745
746 // Introduce smooth curved trajectories to particle-change
747 //
750
751 return &fParticleChange ;
752}
#define State(theXInfo)
@ fStopAndKill
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
void SetPointerToVectorOfAuxiliaryPoints(std::vector< G4ThreeVector > *theNewVectorPointer)
virtual void Initialize(const G4Track &)
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()
std::vector< G4ThreeVector > * GimmeTrajectoryVectorAndForgetIt() const
G4double GetVelocity() const
G4StepPoint * GetPreStepPoint() const
G4double CalculateVelocityForOpticalPhoton() const
Definition: G4Track.cc:243
const G4ParticleDefinition * GetParticleDefinition() const
G4double GetGlobalTime() const
G4double GetLocalTime() const
G4double GetKineticEnergy() const
G4double GetStepLength() const
void ProposeTrackStatus(G4TrackStatus status)

Referenced by G4DNABrownianTransportation::AlongStepDoIt().

◆ AlongStepGetPhysicalInteractionLength()

G4double G4ITTransportation::AlongStepGetPhysicalInteractionLength ( const G4Track track,
G4double  ,
G4double  currentMinimumStep,
G4double currentSafety,
G4GPILSelection selection 
)
virtual

Implements G4VProcess.

Reimplemented in G4DNABrownianTransportation.

Definition at line 218 of file G4ITTransportation.cc.

224{
225 G4double geometryStepLength(-1.0), newSafety(-1.0) ;
226
227 State(fParticleIsLooping) = false ;
228 State(fEndGlobalTimeComputed) = false ;
229 State(fGeometryLimitedStep) = false ;
230
231 // Initial actions moved to StartTrack()
232 // --------------------------------------
233 // Note: in case another process changes touchable handle
234 // it will be necessary to add here (for all steps)
235 State(fCurrentTouchableHandle) = track.GetTouchableHandle();
236
237 // GPILSelection is set to defaule value of CandidateForSelection
238 // It is a return value
239 //
240 *selection = CandidateForSelection ;
241
242 // Get initial Energy/Momentum of the track
243 //
244 const G4DynamicParticle* pParticle = track.GetDynamicParticle() ;
245// const G4ParticleDefinition* pParticleDef = pParticle->GetDefinition() ;
246 G4ThreeVector startMomentumDir = pParticle->GetMomentumDirection() ;
247 G4ThreeVector startPosition = track.GetPosition() ;
248
249 // G4double theTime = track.GetGlobalTime() ;
250
251 // The Step Point safety can be limited by other geometries and/or the
252 // assumptions of any process - it's not always the geometrical safety.
253 // We calculate the starting point's isotropic safety here.
254 //
255 G4ThreeVector OriginShift = startPosition - State(fPreviousSftOrigin) ;
256 G4double MagSqShift = OriginShift.mag2() ;
257 if( MagSqShift >= sqr(State(fPreviousSafety)) )
258 {
259 currentSafety = 0.0 ;
260 }
261 else
262 {
263 currentSafety = State(fPreviousSafety) - std::sqrt(MagSqShift) ;
264 }
265
266 // Is the particle charged ?
267 //
268 G4double particleCharge = pParticle->GetCharge() ;
269
270
271 // There is no need to locate the current volume. It is Done elsewhere:
272 // On track construction
273 // By the tracking, after all AlongStepDoIts, in "Relocation"
274
275 // Check whether the particle have an (EM) field force exerting upon it
276 //
277 G4FieldManager* fieldMgr=0;
278 G4bool fieldExertsForce = false ;
279 if( (particleCharge != 0.0) )
280 {
282 if (fieldMgr != 0)
283 {
284 // Message the field Manager, to configure it for this track
285 fieldMgr->ConfigureForTrack( &track );
286 // Moved here, in order to allow a transition
287 // from a zero-field status (with fieldMgr->(field)0
288 // to a finite field status
289
290 // If the field manager has no field, there is no field !
291 fieldExertsForce = (fieldMgr->GetDetectorField() != 0);
292 }
293 }
294
295 // G4cout << " G4Transport: field exerts force= " << fieldExertsForce
296 // << " fieldMgr= " << fieldMgr << G4endl;
297
298 // Choose the calculation of the transportation: Field or not
299 //
300 if( !fieldExertsForce )
301 {
302 G4double linearStepLength ;
303 if( fShortStepOptimisation && (currentMinimumStep <= currentSafety) )
304 {
305 // The Step is guaranteed to be taken
306 //
307 geometryStepLength = currentMinimumStep ;
308 State(fGeometryLimitedStep) = false ;
309 }
310 else
311 {
312 // Find whether the straight path intersects a volume
313 //
314 linearStepLength = fLinearNavigator->ComputeStep( startPosition,
315 startMomentumDir,
316 currentMinimumStep,
317 newSafety) ;
318 // Remember last safety origin & value.
319 //
320 State(fPreviousSftOrigin) = startPosition ;
321 State(fPreviousSafety) = newSafety ;
322 // fpSafetyHelper->SetCurrentSafety( newSafety, startPosition);
323
324 // The safety at the initial point has been re-calculated:
325 //
326 currentSafety = newSafety ;
327
328 State(fGeometryLimitedStep)= (linearStepLength <= currentMinimumStep);
329 if( State(fGeometryLimitedStep) )
330 {
331 // The geometry limits the Step size (an intersection was found.)
332 geometryStepLength = linearStepLength ;
333 }
334 else
335 {
336 // The full Step is taken.
337 geometryStepLength = currentMinimumStep ;
338 }
339 }
340 State(endpointDistance) = geometryStepLength ;
341
342 // Calculate final position
343 //
344 State(fTransportEndPosition) = startPosition+geometryStepLength*startMomentumDir ;
345
346 // Momentum direction, energy and polarisation are unchanged by transport
347 //
348 State(fTransportEndMomentumDir) = startMomentumDir ;
349 State(fTransportEndKineticEnergy) = track.GetKineticEnergy() ;
350 State(fTransportEndSpin) = track.GetPolarization();
351 State(fParticleIsLooping) = false ;
352 State(fMomentumChanged) = false ;
353 State(fEndGlobalTimeComputed) = true ;
354 State(theInteractionTimeLeft) = State(endpointDistance)/track.GetVelocity();
355 State(fCandidateEndGlobalTime) = State(theInteractionTimeLeft)+track.GetGlobalTime();
356
357 // G4cout << "track.GetVelocity() : " << track.GetVelocity() << G4endl;
358 // G4cout << "State(endpointDistance) : " << G4BestUnit(State(endpointDistance),"Length") << G4endl;
359 // G4cout << "State(theInteractionTimeLeft) : " << G4BestUnit(State(theInteractionTimeLeft),"Time") << G4endl;
360 // G4cout << "track.GetGlobalTime() : " << G4BestUnit(track.GetGlobalTime(),"Time") << G4endl;
361
362 }
363 else // A field exerts force
364 {
365
366 G4ExceptionDescription exceptionDescription;
367 exceptionDescription << "ITTransportation does not support external fields.";
368 exceptionDescription << " If you are dealing with a tradiational MC simulation, ";
369 exceptionDescription << "please use G4Transportation.";
370
371 G4Exception("G4ITTransportation::AlongStepGetPhysicalInteractionLength","NoExternalFieldSupport",
372 FatalException,exceptionDescription);
373 /*
374 G4double momentumMagnitude = pParticle->GetTotalMomentum() ;
375 // G4ThreeVector EndUnitMomentum ;
376 G4double lengthAlongCurve ;
377 G4double restMass = pParticleDef->GetPDGMass() ;
378
379 fFieldPropagator->SetChargeMomentumMass( particleCharge, // in e+ units
380 momentumMagnitude, // in Mev/c
381 restMass ) ;
382
383 G4ThreeVector spin = track.GetPolarization() ;
384 G4FieldTrack aFieldTrack = G4FieldTrack( startPosition,
385 track.GetMomentumDirection(),
386 0.0,
387 track.GetKineticEnergy(),
388 restMass,
389 track.GetVelocity(),
390 track.GetGlobalTime(), // Lab.
391 track.GetProperTime(), // Part.
392 &spin ) ;
393 if( currentMinimumStep > 0 )
394 {
395 // Do the Transport in the field (non recti-linear)
396 //
397 lengthAlongCurve = fFieldPropagator->ComputeStep( aFieldTrack,
398 currentMinimumStep,
399 currentSafety,
400 track.GetVolume() ) ;
401 State(fGeometryLimitedStep)= lengthAlongCurve < currentMinimumStep;
402 if( State(fGeometryLimitedStep) )
403 {
404 geometryStepLength = lengthAlongCurve ;
405 }
406 else
407 {
408 geometryStepLength = currentMinimumStep ;
409 }
410
411 // Remember last safety origin & value.
412 //
413 State(fPreviousSftOrigin) = startPosition ;
414 State(fPreviousSafety) = currentSafety ;
415 fpSafetyHelper->SetCurrentSafety( newSafety, startPosition);
416 }
417 else
418 {
419 geometryStepLength = lengthAlongCurve= 0.0 ;
420 State(fGeometryLimitedStep) = false ;
421 }
422
423 // Get the End-Position and End-Momentum (Dir-ection)
424 //
425 State(fTransportEndPosition) = aFieldTrack.GetPosition() ;
426
427 // Momentum: Magnitude and direction can be changed too now ...
428 //
429 State(fMomentumChanged) = true ;
430 State(fTransportEndMomentumDir) = aFieldTrack.GetMomentumDir() ;
431
432 State(fTransportEndKineticEnergy) = aFieldTrack.GetKineticEnergy() ;
433
434 if( fFieldPropagator->GetCurrentFieldManager()->DoesFieldChangeEnergy() )
435 {
436 // If the field can change energy, then the time must be integrated
437 // - so this should have been updated
438 //
439 State(fCandidateEndGlobalTime) = aFieldTrack.GetLabTimeOfFlight();
440 State(fEndGlobalTimeComputed) = true;
441
442 State(theInteractionTimeLeft) = State(fCandidateEndGlobalTime) - track.GetGlobalTime() ;
443
444 // was ( State(fCandidateEndGlobalTime) != track.GetGlobalTime() );
445 // a cleaner way is to have FieldTrack knowing whether time is updated.
446 }
447 else
448 {
449 // The energy should be unchanged by field transport,
450 // - so the time changed will be calculated elsewhere
451 //
452 State(fEndGlobalTimeComputed) = false;
453
454 // Check that the integration preserved the energy
455 // - and if not correct this!
456 G4double startEnergy= track.GetKineticEnergy();
457 G4double endEnergy= State(fTransportEndKineticEnergy);
458
459 static G4int no_inexact_steps=0, no_large_ediff;
460 G4double absEdiff = std::fabs(startEnergy- endEnergy);
461 if( absEdiff > perMillion * endEnergy )
462 {
463 no_inexact_steps++;
464 // Possible statistics keeping here ...
465 }
466#ifdef G4VERBOSE
467 if( fVerboseLevel > 1 )
468 {
469 if( std::fabs(startEnergy- endEnergy) > perThousand * endEnergy )
470 {
471 static G4int no_warnings= 0, warnModulo=1, moduloFactor= 10;
472 no_large_ediff ++;
473 if( (no_large_ediff% warnModulo) == 0 )
474 {
475 no_warnings++;
476 G4cout << "WARNING - G4Transportation::AlongStepGetPIL() "
477 << " Energy change in Step is above 1^-3 relative value. " << G4endl
478 << " Relative change in 'tracking' step = "
479 << std::setw(15) << (endEnergy-startEnergy)/startEnergy << G4endl
480 << " Starting E= " << std::setw(12) << startEnergy / MeV << " MeV " << G4endl
481 << " Ending E= " << std::setw(12) << endEnergy / MeV << " MeV " << G4endl;
482 G4cout << " Energy has been corrected -- however, review"
483 << " field propagation parameters for accuracy." << G4endl;
484 if( (fVerboseLevel > 2 ) || (no_warnings<4) || (no_large_ediff == warnModulo * moduloFactor) )
485 {
486 G4cout << " These include EpsilonStepMax(/Min) in G4FieldManager "
487 << " which determine fractional error per step for integrated quantities. " << G4endl
488 << " Note also the influence of the permitted number of integration steps."
489 << G4endl;
490 }
491 G4cerr << "ERROR - G4Transportation::AlongStepGetPIL()" << G4endl
492 << " Bad 'endpoint'. Energy change detected"
493 << " and corrected. "
494 << " Has occurred already "
495 << no_large_ediff << " times." << G4endl;
496 if( no_large_ediff == warnModulo * moduloFactor )
497 {
498 warnModulo *= moduloFactor;
499 }
500 }
501 }
502 } // end of if (fVerboseLevel)
503#endif
504 // Correct the energy for fields that conserve it
505 // This - hides the integration error
506 // - but gives a better physical answer
507 State(fTransportEndKineticEnergy)= track.GetKineticEnergy();
508 }
509
510 State(fTransportEndSpin) = aFieldTrack.GetSpin();
511 State(fParticleIsLooping) = fFieldPropagator->IsParticleLooping() ;
512 State(endpointDistance) = (State(fTransportEndPosition) - startPosition).mag() ;
513 // State(theInteractionTimeLeft) = track.GetVelocity()/State(endpointDistance) ;
514*/
515 }
516
517 // If we are asked to go a step length of 0, and we are on a boundary
518 // then a boundary will also limit the step -> we must flag this.
519 //
520 if( currentMinimumStep == 0.0 )
521 {
522 if( currentSafety == 0.0 )
523 {
524 State(fGeometryLimitedStep) = true ;
525// G4cout << "!!!! Safety is NULL, on the Boundary !!!!!" << G4endl;
526// G4cout << " Track position : " << track.GetPosition() /nanometer << G4endl;
527 }
528 }
529
530 // Update the safety starting from the end-point,
531 // if it will become negative at the end-point.
532 //
533 if( currentSafety < State(endpointDistance) )
534 {
535 // if( particleCharge == 0.0 )
536 // G4cout << " Avoiding call to ComputeSafety : charge = 0.0 " << G4endl;
537
538 if( particleCharge != 0.0 )
539 {
540
541 G4double endSafety =
542 fLinearNavigator->ComputeSafety( State(fTransportEndPosition)) ;
543 currentSafety = endSafety ;
544 State(fPreviousSftOrigin) = State(fTransportEndPosition) ;
545 State(fPreviousSafety) = currentSafety ;
546 fpSafetyHelper->SetCurrentSafety( currentSafety, State(fTransportEndPosition));
547
548 // Because the Stepping Manager assumes it is from the start point,
549 // add the StepLength
550 //
551 currentSafety += State(endpointDistance) ;
552
553#ifdef G4DEBUG_TRANSPORT
554 G4cout.precision(12) ;
555 G4cout << "***G4Transportation::AlongStepGPIL ** " << G4endl ;
556 G4cout << " Called Navigator->ComputeSafety at " << State(fTransportEndPosition)
557 << " and it returned safety= " << endSafety << G4endl ;
558 G4cout << " Adding endpoint distance " << State(endpointDistance)
559 << " to obtain pseudo-safety= " << currentSafety << G4endl ;
560#endif
561 }
562 }
563
564 // fParticleChange.ProposeTrueStepLength(geometryStepLength) ;
565
566 return geometryStepLength ;
567}
@ FatalException
@ CandidateForSelection
bool G4bool
Definition: G4Types.hh:67
double mag2() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetCharge() const
virtual void ConfigureForTrack(const G4Track *)
const G4Field * GetDetectorField() const
virtual G4double ComputeStep(const G4ThreeVector &pGlobalPoint, const G4ThreeVector &pDirection, const G4double pCurrentProposedStepLength, G4double &pNewSafety)
virtual G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4double pProposedMaxLength=DBL_MAX, const G4bool keepState=false)
G4FieldManager * FindAndSetFieldManager(G4VPhysicalVolume *pCurrentPhysVol)
void SetCurrentSafety(G4double val, const G4ThreeVector &pos)
G4double GetVelocity() const
G4VPhysicalVolume * GetVolume() const
const G4ThreeVector & GetPosition() const
const G4DynamicParticle * GetDynamicParticle() const
const G4TouchableHandle & GetTouchableHandle() const
const G4ThreeVector & GetPolarization() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
T sqr(const T &x)
Definition: templates.hh:145

Referenced by G4DNABrownianTransportation::AlongStepGetPhysicalInteractionLength().

◆ AtRestDoIt()

virtual G4VParticleChange * G4ITTransportation::AtRestDoIt ( const G4Track ,
const G4Step  
)
inlinevirtual

Implements G4VProcess.

Definition at line 91 of file G4ITTransportation.hh.

95 {
96 return 0;
97 }

◆ AtRestGetPhysicalInteractionLength()

virtual G4double G4ITTransportation::AtRestGetPhysicalInteractionLength ( const G4Track ,
G4ForceCondition  
)
inlinevirtual

Implements G4VProcess.

Definition at line 82 of file G4ITTransportation.hh.

86 {
87 return -1.0;
88 }

◆ BuildPhysicsTable()

virtual void G4ITTransportation::BuildPhysicsTable ( const G4ParticleDefinition )
inlinevirtual

Reimplemented from G4VITProcess.

Reimplemented in G4DNABrownianTransportation.

Definition at line 69 of file G4ITTransportation.hh.

69{;}

◆ ComputeStep()

void G4ITTransportation::ComputeStep ( const G4Track track,
const G4Step ,
const double  timeStep,
double &  spaceStep 
)
virtual

Reimplemented in G4DNABrownianTransportation.

Definition at line 570 of file G4ITTransportation.cc.

574{
575 const G4DynamicParticle* pParticle = track.GetDynamicParticle() ;
576 G4ThreeVector startMomentumDir = pParticle->GetMomentumDirection() ;
577 G4ThreeVector startPosition = track.GetPosition() ;
578
579 track.CalculateVelocity();
580 G4double initialVelocity = track.GetVelocity() ;
581
582 State(fGeometryLimitedStep) = false;
583
584 /////////////////////////
585 // !!! CASE NO FIELD !!!
586 /////////////////////////
587 State(fCandidateEndGlobalTime) = timeStep + track.GetGlobalTime();
588 State(fEndGlobalTimeComputed) = true ;
589
590 // Choose the calculation of the transportation: Field or not
591 //
592 if( !State(fMomentumChanged) )
593 {
594 // G4cout << "Momentum has not changed" << G4endl;
595 fParticleChange.ProposeVelocity(initialVelocity);
596 oPhysicalStep = initialVelocity*timeStep ;
597
598 // Calculate final position
599 //
600 State(fTransportEndPosition) = startPosition + oPhysicalStep*startMomentumDir ;
601 }
602}
G4double CalculateVelocity() const
Definition: G4Track.cc:210

Referenced by G4ITStepProcessor::DoDefinePhysicalStepLength(), and G4ITStepProcessor::FindTransportationStep().

◆ DoesGlobalFieldExist()

G4bool G4ITTransportation::DoesGlobalFieldExist ( )
protected

Definition at line 201 of file G4ITTransportation.cc.

202{
203 G4TransportationManager* transportMgr;
205
206 // fFieldExists= transportMgr->GetFieldManager()->DoesFieldExist();
207 // return fFieldExists;
208 return transportMgr->GetFieldManager()->DoesFieldExist();
209}
G4bool DoesFieldExist() const
G4FieldManager * GetFieldManager() const

Referenced by StartTracking().

◆ EnableShortStepOptimisation()

void G4ITTransportation::EnableShortStepOptimisation ( G4bool  optimise = true)
inline

◆ GetMaxEnergyKilled()

G4double G4ITTransportation::GetMaxEnergyKilled ( ) const
inline

◆ GetPropagatorInField()

G4PropagatorInField * G4ITTransportation::GetPropagatorInField ( )

◆ GetSumEnergyKilled()

G4double G4ITTransportation::GetSumEnergyKilled ( ) const
inline

◆ GetThresholdImportantEnergy()

G4double G4ITTransportation::GetThresholdImportantEnergy ( ) const
inline

◆ GetThresholdTrials()

G4int G4ITTransportation::GetThresholdTrials ( ) const
inline

◆ GetThresholdWarningEnergy()

G4double G4ITTransportation::GetThresholdWarningEnergy ( ) const
inline

◆ GetVerboseLevel()

G4int G4ITTransportation::GetVerboseLevel ( ) const
inline

◆ InstantiateProcessState()

G4bool G4ITTransportation::InstantiateProcessState ( )
inlineprotected

Definition at line 239 of file G4ITTransportation.hh.

239{ return fInstantiateProcessState; }

◆ PostStepDoIt()

G4VParticleChange * G4ITTransportation::PostStepDoIt ( const G4Track track,
const G4Step  
)
virtual

Implements G4VProcess.

Reimplemented in G4DNABrownianTransportation.

Definition at line 774 of file G4ITTransportation.cc.

776{
777 // G4cout << "G4ITTransportation::PostStepDoIt" << G4endl;
778 G4TouchableHandle retCurrentTouchable ; // The one to return
779 G4bool isLastStep= false;
780
781 // Initialize ParticleChange (by setting all its members equal
782 // to corresponding members in G4Track)
783 fParticleChange.Initialize(track) ; // To initialise TouchableChange
784
786
787 // If the Step was determined by the volume boundary,
788 // logically relocate the particle
789
790 if(State(fGeometryLimitedStep))
791 {
792
793 // G4cout << "Step is limited by geometry " << "track ID : " << track.GetTrackID() << G4endl;
794
795 // fCurrentTouchable will now become the previous touchable,
796 // and what was the previous will be freed.
797 // (Needed because the preStepPoint can point to the previous touchable)
798
799 if( State(fCurrentTouchableHandle)->GetVolume() == 0 )
800 {
801 G4ExceptionDescription exceptionDescription ;
802 exceptionDescription << "No current touchable found " ;
803 G4Exception(" G4ITTransportation::PostStepDoIt","G4ITTransportation001",
804 FatalErrorInArgument,exceptionDescription);
805 }
806
809 LocateGlobalPointAndUpdateTouchableHandle( track.GetPosition(),
810 track.GetMomentumDirection(),
811 State(fCurrentTouchableHandle),
812 true ) ;
813 // Check whether the particle is out of the world volume
814 // If so it has exited and must be killed.
815 //
816 if( State(fCurrentTouchableHandle)->GetVolume() == 0 )
817 {
818#ifdef G4VERBOSE
819 if(fVerboseLevel > 0)
820 {
821 G4cout << "Track position : " << track.GetPosition() / nanometer << " [nm]"
822 << " Track ID : " << track.GetTrackID()<< G4endl;
823 G4cout << "G4ITTransportation will killed the track because State(fCurrentTouchableHandle)->GetVolume() == 0"<< G4endl;
824 }
825#endif
827 }
828
829 retCurrentTouchable = State(fCurrentTouchableHandle) ;
830
831// G4cout << "Current volume : " << track.GetVolume()->GetName()
832// << " Next volume : "
833// << (State(fCurrentTouchableHandle)->GetVolume() ? State(fCurrentTouchableHandle)->GetVolume()->GetName():"OutWorld")
834// << " Position : " << track.GetPosition() / nanometer
835// << " track ID : " << track.GetTrackID()
836// << G4endl;
837
838 fParticleChange.SetTouchableHandle( State(fCurrentTouchableHandle) ) ;
839
840 // Update the Step flag which identifies the Last Step in a volume
843
844#ifdef G4DEBUG_TRANSPORT
845 // Checking first implementation of flagging Last Step in Volume
848
849 if( ! (exiting || entering) )
850 {
851 G4cout << " Transport> : Proposed isLastStep= " << isLastStep
852 << " Exiting " << fLinearNavigator->ExitedMotherVolume()
853 << " Entering " << fLinearNavigator->EnteredDaughterVolume()
854 << " Track position : " << track.GetPosition() /nanometer << " [nm]"
855 << G4endl;
856 G4cout << " Track position : " << track.GetPosition() /nanometer << G4endl;
857 }
858#endif
859 }
860 else // fGeometryLimitedStep is false
861 {
862 // This serves only to move the Navigator's location
863 //
865
866 // The value of the track's current Touchable is retained.
867 // (and it must be correct because we must use it below to
868 // overwrite the (unset) one in particle change)
869 // It must be fCurrentTouchable too ??
870 //
872 retCurrentTouchable = track.GetTouchableHandle() ;
873
874 isLastStep= false;
875#ifdef G4DEBUG_TRANSPORT
876 // Checking first implementation of flagging Last Step in Volume
877 G4cout << " Transport> Proposed isLastStep= " << isLastStep
878 << " Geometry did not limit step. Position : "
879 << track.GetPosition()/ nanometer << G4endl;
880#endif
881 } // endif ( fGeometryLimitedStep )
882
884
885 const G4VPhysicalVolume* pNewVol = retCurrentTouchable->GetVolume() ;
886 const G4Material* pNewMaterial = 0 ;
887 const G4VSensitiveDetector* pNewSensitiveDetector = 0 ;
888
889 if( pNewVol != 0 )
890 {
891 pNewMaterial= pNewVol->GetLogicalVolume()->GetMaterial();
892 pNewSensitiveDetector= pNewVol->GetLogicalVolume()->GetSensitiveDetector();
893 }
894
895 // ( <const_cast> pNewMaterial ) ;
896 // ( <const_cast> pNewSensitiveDetector) ;
897
900
901 const G4MaterialCutsCouple* pNewMaterialCutsCouple = 0;
902 if( pNewVol != 0 )
903 {
904 pNewMaterialCutsCouple=pNewVol->GetLogicalVolume()->GetMaterialCutsCouple();
905 }
906
907 if( pNewVol!=0 && pNewMaterialCutsCouple!=0 && pNewMaterialCutsCouple->GetMaterial()!=pNewMaterial )
908 {
909 // for parametrized volume
910 //
911 pNewMaterialCutsCouple =
913 ->GetMaterialCutsCouple(pNewMaterial,
914 pNewMaterialCutsCouple->GetProductionCuts());
915 }
916 fParticleChange.SetMaterialCutsCoupleInTouchable( pNewMaterialCutsCouple );
917
918 // temporarily until Get/Set Material of ParticleChange,
919 // and StepPoint can be made const.
920 // Set the touchable in ParticleChange
921 // this must always be done because the particle change always
922 // uses this value to overwrite the current touchable pointer.
923 //
924 fParticleChange.SetTouchableHandle(retCurrentTouchable) ;
925
926 return &fParticleChange ;
927}
@ FatalErrorInArgument
G4bool ExitedMotherVolume() const
G4bool EnteredDaughterVolume() const
virtual void LocateGlobalPointWithinVolume(const G4ThreeVector &position)
void SetGeometricallyLimitedStep()
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)
void SetSensitiveDetectorInTouchable(G4VSensitiveDetector *fSensitiveDetector)
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
static G4ProductionCutsTable * GetProductionCutsTable()
G4TrackStatus GetTrackStatus() const
G4int GetTrackID() const
const G4ThreeVector & GetMomentumDirection() const
void ProposeLastStepInVolume(G4bool flag)
G4LogicalVolume * GetLogicalVolume() const
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:44

Referenced by G4DNABrownianTransportation::PostStepDoIt().

◆ PostStepGetPhysicalInteractionLength()

G4double G4ITTransportation::PostStepGetPhysicalInteractionLength ( const G4Track ,
G4double  ,
G4ForceCondition pForceCond 
)
virtual

Implements G4VProcess.

Definition at line 760 of file G4ITTransportation.cc.

766{
767 *pForceCond = Forced ;
768 return DBL_MAX ; // was kInfinity ; but convention now is DBL_MAX
769}
@ Forced
#define DBL_MAX
Definition: templates.hh:83

◆ ResetKilledStatistics()

void G4ITTransportation::ResetKilledStatistics ( G4int  report = 1)
inline

◆ SetInstantiateProcessState()

void G4ITTransportation::SetInstantiateProcessState ( G4bool  flag)
inlineprotected

Definition at line 236 of file G4ITTransportation.hh.

237 { fInstantiateProcessState = flag; }

Referenced by G4ITTransportation(), and G4DNABrownianTransportation::StartTracking().

◆ SetPropagatorInField()

void G4ITTransportation::SetPropagatorInField ( G4PropagatorInField pFieldPropagator)

◆ SetThresholdImportantEnergy()

void G4ITTransportation::SetThresholdImportantEnergy ( G4double  newEnImp)
inline

◆ SetThresholdTrials()

void G4ITTransportation::SetThresholdTrials ( G4int  newMaxTrials)
inline

◆ SetThresholdWarningEnergy()

void G4ITTransportation::SetThresholdWarningEnergy ( G4double  newEnWarn)
inline

◆ SetVerboseLevel()

void G4ITTransportation::SetVerboseLevel ( G4int  verboseLevel)
inline

◆ StartTracking()

void G4ITTransportation::StartTracking ( G4Track aTrack)
virtual

Reimplemented from G4VITProcess.

Reimplemented in G4DNABrownianTransportation.

Definition at line 933 of file G4ITTransportation.cc.

934{
936 if(fInstantiateProcessState)
937 G4VITProcess::fpState = new G4ITTransportationState();
938 // Will set in the same time fTransportationState
939
940 // The actions here are those that were taken in AlongStepGPIL
941 // when track.GetCurrentStepNumber()==1
942
943 // reset safety value and center
944 //
945 // State(fPreviousSafety) = 0.0 ;
946 // State(fPreviousSftOrigin) = G4ThreeVector(0.,0.,0.) ;
947
948 // reset looping counter -- for motion in field
949 // State(fNoLooperTrials)= 0;
950 // Must clear this state .. else it depends on last track's value
951 // --> a better solution would set this from state of suspended track TODO ?
952 // Was if( aTrack->GetCurrentStepNumber()==1 ) { .. }
953
954 // ChordFinder reset internal state
955 //
957 {
959 // Resets all state of field propagator class (ONLY)
960 // including safety values (in case of overlaps and to wipe for first track).
961
962 // G4ChordFinder* chordF= fFieldPropagator->GetChordFinder();
963 // if( chordF ) chordF->ResetStepEstimate();
964 }
965
966 // Make sure to clear the chord finders of all fields (ie managers)
968 fieldMgrStore->ClearAllChordFindersState();
969
970 // Update the current touchable handle (from the track's)
971 //
972 State(fCurrentTouchableHandle) = track->GetTouchableHandle();
973
975}
static G4FieldManagerStore * GetInstance()
virtual void StartTracking(G4Track *)
Definition: G4VITProcess.cc:81
virtual void StartTracking(G4Track *)
Definition: G4VProcess.cc:125

Referenced by G4DNABrownianTransportation::StartTracking().

Member Data Documentation

◆ fFieldPropagator

G4PropagatorInField* G4ITTransportation::fFieldPropagator
protected

◆ fLinearNavigator

G4ITNavigator* G4ITTransportation::fLinearNavigator
protected

◆ fMaxEnergyKilled

G4double G4ITTransportation::fMaxEnergyKilled
protected

Definition at line 223 of file G4ITTransportation.hh.

Referenced by AlongStepDoIt(), G4ITTransportation(), and ~G4ITTransportation().

◆ fParticleChange

◆ fpSafetyHelper

G4SafetyHelper* G4ITTransportation::fpSafetyHelper
protected

◆ fShortStepOptimisation

G4bool G4ITTransportation::fShortStepOptimisation
protected

◆ fSumEnergyKilled

G4double G4ITTransportation::fSumEnergyKilled
protected

Definition at line 222 of file G4ITTransportation.hh.

Referenced by AlongStepDoIt(), G4ITTransportation(), and ~G4ITTransportation().

◆ fThreshold_Important_Energy

G4double G4ITTransportation::fThreshold_Important_Energy
protected

Definition at line 213 of file G4ITTransportation.hh.

Referenced by AlongStepDoIt(), and G4ITTransportation().

◆ fThreshold_Warning_Energy

G4double G4ITTransportation::fThreshold_Warning_Energy
protected

Definition at line 212 of file G4ITTransportation.hh.

Referenced by AlongStepDoIt(), and G4ITTransportation().

◆ fThresholdTrials

G4int G4ITTransportation::fThresholdTrials
protected

Definition at line 214 of file G4ITTransportation.hh.

Referenced by AlongStepDoIt(), and G4ITTransportation().

◆ fTransportationState

G4ITTransportationState* const& G4ITTransportation::fTransportationState
protected

Definition at line 196 of file G4ITTransportation.hh.

◆ fUnimportant_Energy

G4double G4ITTransportation::fUnimportant_Energy
protected

Definition at line 217 of file G4ITTransportation.hh.

Referenced by G4ITTransportation().

◆ fVerboseLevel


The documentation for this class was generated from the following files: