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

#include <G4CoupledTransportation.hh>

+ Inheritance diagram for G4CoupledTransportation:

Public Member Functions

 G4CoupledTransportation (G4int verbosityLevel=0)
 
 ~G4CoupledTransportation ()
 
G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
 
G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)
 
G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &, G4double previousStepSize, G4ForceCondition *pForceCond)
 
G4PropagatorInFieldGetPropagatorInField ()
 
void SetPropagatorInField (G4PropagatorInField *pFieldPropagator)
 
G4double GetThresholdWarningEnergy () const
 
G4double GetThresholdImportantEnergy () const
 
G4int GetThresholdTrials () const
 
void SetThresholdWarningEnergy (G4double newEnWarn)
 
void SetThresholdImportantEnergy (G4double newEnImp)
 
void SetThresholdTrials (G4int newMaxTrials)
 
void SetHighLooperThresholds ()
 
void SetLowLooperThresholds ()
 
void PushThresholdsToLogger ()
 
void ReportLooperThresholds ()
 
G4double GetMaxEnergyKilled () const
 
G4double GetSumEnergyKilled () const
 
void ResetKilledStatistics (G4int report=1)
 
G4bool IsFirstStepInAnyVolume () const
 
G4bool IsLastStepInAnyVolume () const
 
G4bool IsFirstStepInMassVolume () const
 
G4bool IsLastStepInMassVolume () const
 
void StartTracking (G4Track *aTrack)
 
void EndTracking ()
 
G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
void PrintStatistics (std::ostream &outStr) const
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4VProcessoperator= (const G4VProcess &)=delete
 
G4bool operator== (const G4VProcess &right) const
 
G4bool 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
 
virtual void ProcessDescription (std::ostream &outfile) const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Static Public Member Functions

static G4bool EnableMagneticMoment (G4bool useMoment=true)
 
static G4bool EnableGravity (G4bool useGravity)
 
static void SetSilenceLooperWarnings (G4bool val)
 
static G4bool GetSilenceLooperWarnings ()
 
static void SetSignifyStepsInAnyVolume (G4bool anyVol)
 
static G4bool GetSignifyStepsInAnyVolume ()
 
static G4bool EnableUseMagneticMoment (G4bool useMoment=true)
 
- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 

Protected Member Functions

G4bool DoesAnyFieldExist ()
 
void ReportInexactEnergy (G4double startEnergy, G4double endEnergy)
 
void ReportMove (G4ThreeVector OldVector, G4ThreeVector NewVector, const G4String &Quantity)
 
void ReportMissingLogger (const char *methodName)
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Friends

class G4Transportation
 

Additional Inherited Members

- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft = -1.0
 
G4double currentInteractionLength = -1.0
 
G4double theInitialNumberOfInteractionLength = -1.0
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType = fNotDefined
 
G4int theProcessSubType = -1
 
G4double thePILfactor = 1.0
 
G4int verboseLevel = 0
 
G4bool enableAtRestDoIt = true
 
G4bool enableAlongStepDoIt = true
 
G4bool enablePostStepDoIt = true
 

Detailed Description

Definition at line 65 of file G4CoupledTransportation.hh.

Constructor & Destructor Documentation

◆ G4CoupledTransportation()

G4CoupledTransportation::G4CoupledTransportation ( G4int  verbosityLevel = 0)

Definition at line 62 of file G4CoupledTransportation.cc.

63 : G4VProcess( G4String("CoupledTransportation"), fTransportation ),
64 fTransportEndPosition(0.0, 0.0, 0.0),
65 fTransportEndMomentumDir(0.0, 0.0, 0.0),
66 fTransportEndKineticEnergy(0.0),
67 fTransportEndSpin(0.0, 0.0, 0.0),
68 fMomentumChanged(false),
69 fEndGlobalTimeComputed(false),
70 fCandidateEndGlobalTime(0.0),
71 fParticleIsLooping( false ),
72 fNewTrack( true ),
73 fPreviousSftOrigin( 0.,0.,0. ),
74 fPreviousMassSafety( 0.0 ),
75 fPreviousFullSafety( 0.0 ),
76 fMassGeometryLimitedStep( false ),
77 fAnyGeometryLimitedStep( false ),
78 fEndpointDistance( -1.0 ),
79 fSumEnergyKilled( 0.0 ), fMaxEnergyKilled( 0.0 ),
80 fFirstStepInMassVolume( true ),
81 fFirstStepInAnyVolume( true )
82{
84 SetVerboseLevel(verbosity);
85
86 G4TransportationManager* transportMgr ;
87
89
90 fMassNavigator = transportMgr->GetNavigatorForTracking() ;
91 fFieldPropagator = transportMgr->GetPropagatorInField() ;
92 // fGlobalFieldMgr = transportMgr->GetFieldManager() ;
93 fNavigatorId= transportMgr->ActivateNavigator( fMassNavigator );
94 if( verboseLevel > 0 )
95 {
96 G4cout << " G4CoupledTransportation constructor: ----- " << G4endl;
97 G4cout << " Verbose level is " << verboseLevel << G4endl;
98 G4cout << " Navigator Id obtained in G4CoupledTransportation constructor "
99 << fNavigatorId << G4endl;
100 G4cout << " Reports First/Last in "
101 << (fSignifyStepInAnyVolume ? " any " : " mass " )
102 << " geometry " << G4endl;
103 }
104 fPathFinder= G4PathFinder::GetInstance();
105 fpSafetyHelper = transportMgr->GetSafetyHelper(); // New
106
107 fpLogger = new G4TransportationLogger("G4Transportation", verbosity);
108
110 // Use the old defaults: Warning = 100 MeV, Important = 250 MeV, No Trials = 10;
111
113 // Should be done by Set methods in SetHighLooperThresholds -- making sure
114
115 static G4ThreadLocal G4TouchableHandle* pNullTouchableHandle = 0;
116 if ( !pNullTouchableHandle) { pNullTouchableHandle = new G4TouchableHandle; }
117 fCurrentTouchableHandle = *pNullTouchableHandle;
118 // Points to (G4VTouchable*) 0
119
120 fAnyFieldExists= DoesAnyFieldExist();
121}
@ fTransportation
G4ReferenceCountedHandle< G4VTouchable > G4TouchableHandle
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static G4PathFinder * GetInstance()
Definition: G4PathFinder.cc:52
static G4TransportationManager * GetTransportationManager()
G4PropagatorInField * GetPropagatorInField() const
G4SafetyHelper * GetSafetyHelper() const
G4Navigator * GetNavigatorForTracking() const
G4int ActivateNavigator(G4Navigator *aNavigator)
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:412
G4int verboseLevel
Definition: G4VProcess.hh:356
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:406
#define G4ThreadLocal
Definition: tls.hh:77

◆ ~G4CoupledTransportation()

G4CoupledTransportation::~G4CoupledTransportation ( )

Definition at line 125 of file G4CoupledTransportation.cc.

126{
127 if( fSumEnergyKilled > 0.0 )
128 {
130 }
131 delete fpLogger;
132}
void PrintStatistics(std::ostream &outStr) const

Member Function Documentation

◆ AlongStepDoIt()

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

Implements G4VProcess.

Definition at line 563 of file G4CoupledTransportation.cc.

565{
566 static G4ThreadLocal G4long noCallsCT_ASDI=0;
567 const char *methodName= "AlongStepDoIt";
568
569 noCallsCT_ASDI++;
570
571 fParticleChange.Initialize(track) ;
572 // sets all its members to the value of corresponding members in G4Track
573
574 // Code specific for Transport
575 //
576 fParticleChange.ProposePosition(fTransportEndPosition) ;
577 fParticleChange.ProposeMomentumDirection(fTransportEndMomentumDir) ;
578 fParticleChange.ProposeEnergy(fTransportEndKineticEnergy) ;
579 fParticleChange.SetMomentumChanged(fMomentumChanged) ;
580
581 fParticleChange.ProposePolarization(fTransportEndSpin);
582
583 G4double deltaTime = 0.0 ;
584
585 // Calculate Lab Time of Flight (ONLY if field Equations used it!)
586 // G4double endTime = fCandidateEndGlobalTime;
587 // G4double delta_time = endTime - startTime;
588
589 G4double startTime = track.GetGlobalTime() ;
590
591 if (!fEndGlobalTimeComputed)
592 {
593 G4double finalInverseVel= DBL_MAX, initialInverseVel=DBL_MAX;
594
595 // The time was not integrated .. make the best estimate possible
596 //
597 G4double finalVelocity = track.GetVelocity() ;
598 if( finalVelocity > 0.0 ) { finalInverseVel= 1.0 / finalVelocity; }
599 G4double initialVelocity = stepData.GetPreStepPoint()->GetVelocity() ;
600 if( initialVelocity > 0.0 ) { initialInverseVel= 1.0 / initialVelocity; }
601 G4double stepLength = track.GetStepLength() ;
602
603 if (finalVelocity > 0.0)
604 {
605 // deltaTime = stepLength/finalVelocity ;
606 G4double meanInverseVelocity = 0.5
607 * (initialInverseVel + finalInverseVel);
608 deltaTime = stepLength * meanInverseVelocity ;
609 }
610 else
611 {
612 deltaTime = stepLength * initialInverseVel ;
613 } // Could do with better estimate for final step (finalVelocity = 0) ?
614
615 fCandidateEndGlobalTime = startTime + deltaTime ;
616 fParticleChange.ProposeLocalTime( track.GetLocalTime() + deltaTime) ;
617 }
618 else
619 {
620 deltaTime = fCandidateEndGlobalTime - startTime ;
621 fParticleChange.ProposeGlobalTime( fCandidateEndGlobalTime ) ;
622 }
623
624 // Now Correct by Lorentz factor to get "proper" deltaTime
625
626 G4double restMass = track.GetDynamicParticle()->GetMass() ;
627 G4double deltaProperTime = deltaTime*( restMass/track.GetTotalEnergy() ) ;
628
629 fParticleChange.ProposeProperTime(track.GetProperTime() + deltaProperTime) ;
630 // fParticleChange. ProposeTrueStepLength( track.GetStepLength() ) ;
631
632 // If the particle is caught looping or is stuck (in very difficult
633 // boundaries) in a magnetic field (doing many steps) THEN this kills it ...
634 //
635 if ( fParticleIsLooping )
636 {
637 G4double endEnergy= fTransportEndKineticEnergy;
638
639 const G4ParticleDefinition* particleType=
640 track.GetDynamicParticle() -> GetParticleDefinition();
641 G4bool stable = particleType->GetPDGStable();
642
643 G4bool candidateForEnd = (endEnergy < fThreshold_Important_Energy)
644 || (fNoLooperTrials >= fThresholdTrials) ;
645
646 if( candidateForEnd && stable )
647 {
648 const G4int electronPDG= 11; // G4Electron::G4Electron()->GetPDGEncoding();
649 G4int particlePDG= particleType->GetPDGEncoding();
650
651 // Kill the looping particle
652 //
653 fParticleChange.ProposeTrackStatus( fStopAndKill ) ;
654
655 // Simple statistics
656 fSumEnergyKilled += endEnergy;
657 fSumEnerSqKilled = endEnergy * endEnergy;
658 fNumLoopersKilled++;
659
660 if( endEnergy > fMaxEnergyKilled ) {
661 fMaxEnergyKilled = endEnergy;
662 fMaxEnergyKilledPDG = particlePDG;
663 }
664
665 if( particleType->GetPDGEncoding() != electronPDG )
666 {
667 fSumEnergyKilled_NonElectron += endEnergy;
668 fSumEnerSqKilled_NonElectron += endEnergy * endEnergy;
669 fNumLoopersKilled_NonElectron++;
670
671 if( endEnergy > fMaxEnergyKilled_NonElectron )
672 {
673 fMaxEnergyKilled_NonElectron = endEnergy;
674 fMaxEnergyKilled_NonElecPDG = particlePDG;
675 }
676 }
677
678 if( endEnergy > fThreshold_Warning_Energy && ! fSilenceLooperWarnings )
679 {
680 fpLogger->ReportLoopingTrack( track, stepData, fNoLooperTrials,
681 noCallsCT_ASDI, methodName );
682 }
683
684 fNoLooperTrials=0;
685 }
686 else
687 {
688 fNoLooperTrials ++;
689
690 fMaxEnergySaved = std::max( endEnergy, fMaxEnergySaved);
691 if( fNoLooperTrials == 1 ) {
692 fSumEnergySaved += endEnergy;
693 if ( !stable )
694 fSumEnergyUnstableSaved += endEnergy;
695 }
696#ifdef G4VERBOSE
697 if( verboseLevel > 2 && ! fSilenceLooperWarnings )
698 {
699 G4cout << " ** G4CoupledTransportation::AlongStepDoIt():"
700 << " Particle is looping but is saved ..." << G4endl
701 << " Number of trials (this track) = " << fNoLooperTrials
702 << G4endl
703 << " Steps by this track: " << track.GetCurrentStepNumber()
704 << G4endl
705 << " Total no of calls to this method (all tracks) = "
706 << noCallsCT_ASDI << G4endl;
707 }
708#endif
709 }
710 }
711 else
712 {
713 fNoLooperTrials=0;
714 }
715
716 // Another (sometimes better way) is to use a user-limit maximum Step size
717 // to alleviate this problem ..
718 // Add smooth curved trajectories to particle-change
719 //
720 // fParticleChange.SetPointerToVectorOfAuxiliaryPoints
721 // (fFieldPropagator->GimmeTrajectoryVectorAndForgetIt() );
722
723 return &fParticleChange ;
724}
@ fStopAndKill
double G4double
Definition: G4Types.hh:83
long G4long
Definition: G4Types.hh:87
bool G4bool
Definition: G4Types.hh:86
G4double GetMass() const
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 ProposeProperTime(G4double finalProperTime)
void ProposeEnergy(G4double finalEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
void ProposeGlobalTime(G4double t)
G4bool GetPDGStable() const
G4double GetVelocity() const
G4StepPoint * GetPreStepPoint() const
G4double GetVelocity() const
G4double GetGlobalTime() const
G4double GetProperTime() const
G4int GetCurrentStepNumber() const
G4double GetLocalTime() const
const G4DynamicParticle * GetDynamicParticle() const
G4double GetStepLength() const
G4double GetTotalEnergy() const
void ReportLoopingTrack(const G4Track &track, const G4Step &stepInfo, G4int numTrials, G4long noCalls, const char *methodName) const
void ProposeTrackStatus(G4TrackStatus status)
#define DBL_MAX
Definition: templates.hh:62

◆ AlongStepGetPhysicalInteractionLength()

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

Implements G4VProcess.

Definition at line 164 of file G4CoupledTransportation.cc.

170{
171 G4double geometryStepLength;
172 G4double startMassSafety= 0.0; // estimated safety for start point (mass geometry)
173 G4double startFullSafety= 0.0; // estimated safety for start point (all geometries)
174 G4double safetyProposal= -1.0; // local copy of proposal
175
176 G4ThreeVector EndUnitMomentum ;
177 G4double lengthAlongCurve = 0.0 ;
178
179 fParticleIsLooping = false ;
180
181 // Initial actions moved to StartTrack()
182 // --------------------------------------
183 // Note: in case another process changes touchable handle
184 // it will be necessary to add here (for all steps)
185 // fCurrentTouchableHandle = aTrack->GetTouchableHandle();
186
187 // GPILSelection is set to defaule value of CandidateForSelection
188 // It is a return value
189 //
190 *selection = CandidateForSelection ;
191
192 fFirstStepInMassVolume = fNewTrack | fMassGeometryLimitedStep ;
193 fFirstStepInAnyVolume = fNewTrack | fAnyGeometryLimitedStep ;
194
195#ifdef G4DEBUG_TRANSPORT
196 G4cout << " CoupledTransport::AlongStep GPIL: "
197 << " 1st-step: any= " <<fFirstStepInAnyVolume << " ( geom= "
198 << fAnyGeometryLimitedStep << " ) "
199 << " mass= " << fFirstStepInMassVolume << " ( geom= "
200 << fMassGeometryLimitedStep << " ) "
201 << " newTrack= " << fNewTrack << G4endl;
202#endif
203
204 // fLastStepInVolume= false;
205 fNewTrack = false;
206
207 // Get initial Energy/Momentum of the track
208 //
209 const G4DynamicParticle* pParticle = track.GetDynamicParticle() ;
210 const G4ParticleDefinition* pParticleDef = pParticle->GetDefinition() ;
211 G4ThreeVector startMomentumDir = pParticle->GetMomentumDirection() ;
212 G4ThreeVector startPosition = track.GetPosition() ;
213 G4VPhysicalVolume* currentVolume= track.GetVolume();
214
215#ifdef G4DEBUG_TRANSPORT
216 if( verboseLevel > 1 )
217 {
218 G4cout << "G4CoupledTransportation::AlongStepGPIL> called in volume "
219 << currentVolume->GetName() << G4endl;
220 }
221#endif
222 // G4double theTime = track.GetGlobalTime() ;
223
224 // The Step Point safety can be limited by other geometries and/or the
225 // assumptions of any process - it's not always the geometrical safety.
226 // We calculate the starting point's isotropic safety here.
227 //
228 G4ThreeVector OriginShift = startPosition - fPreviousSftOrigin ;
229 G4double MagSqShift = OriginShift.mag2() ;
230 startMassSafety = 0.0;
231 startFullSafety= 0.0;
232
233 // Recall that FullSafety <= MassSafety
234 // Original: if( MagSqShift < sqr(fPreviousMassSafety) ) {
235 if( MagSqShift < sqr(fPreviousFullSafety) )
236 {
237 G4double mag_shift= std::sqrt(MagSqShift);
238 startMassSafety = std::max( (fPreviousMassSafety - mag_shift), 0.0);
239 startFullSafety = std::max( (fPreviousFullSafety - mag_shift), 0.0);
240 // Need to be consistent between full safety with Mass safety
241 // in order reproduce results in simple case
242 // --> use same calculation method
243
244 // Only compute full safety if massSafety > 0. Else it remains 0
245 // startFullSafety = fPathFinder->ComputeSafety( startPosition );
246 }
247
248 // Is the particle charged or has it a magnetic moment?
249 //
250 G4double particleCharge = pParticle->GetCharge() ;
251 G4double magneticMoment = pParticle->GetMagneticMoment() ;
252 G4double restMass = pParticle->GetMass() ;
253
254 fMassGeometryLimitedStep = false ; // Set default - alex
255 fAnyGeometryLimitedStep = false;
256
257 // There is no need to locate the current volume. It is Done elsewhere:
258 // On track construction
259 // By the tracking, after all AlongStepDoIts, in "Relocation"
260
261 // Check if the particle has a force, EM or gravitational, exerted on it
262 //
263 G4FieldManager* fieldMgr= nullptr;
264 G4bool fieldExertsForce = false ;
265
266 const G4Field* ptrField= nullptr;
267
268 fieldMgr = fFieldPropagator->FindAndSetFieldManager( track.GetVolume() );
269 G4bool eligibleEM = (particleCharge != 0.0)
270 || ( fUseMagneticMoment && (magneticMoment != 0.0) );
271 G4bool eligibleGrav = fUseGravity && (restMass != 0.0) ;
272
273 if( (fieldMgr!=nullptr) && (eligibleEM||eligibleGrav) )
274 {
275 // Message the field Manager, to configure it for this track
276 //
277 fieldMgr->ConfigureForTrack( &track );
278
279 // The above call can transition from a null field-ptr oto a finite field.
280 // If the field manager has no field ptr, the field is zero
281 // by definition ( = there is no field ! )
282 //
283 ptrField= fieldMgr->GetDetectorField();
284
285 if( ptrField != nullptr)
286 {
287 fieldExertsForce = eligibleEM
288 || ( eligibleGrav && ptrField->IsGravityActive() );
289 }
290 }
291 G4double momentumMagnitude = pParticle->GetTotalMomentum() ;
292
293 if( fieldExertsForce )
294 {
295 auto equationOfMotion= fFieldPropagator->GetCurrentEquationOfMotion();
296
297 G4ChargeState chargeState(particleCharge, // Charge can change (dynamic)
298 magneticMoment,
299 pParticleDef->GetPDGSpin() );
300 if( equationOfMotion )
301 {
302 equationOfMotion->SetChargeMomentumMass( chargeState,
303 momentumMagnitude,
304 restMass );
305 }
306#ifdef G4DEBUG_TRANSPORT
307 else
308 {
309 G4cerr << " ERROR in G4CoupledTransportation> "
310 << "Cannot find valid Equation of motion: " << G4endl;
311 << " Unable to pass Charge, Momentum and Mass " << G4endl;
312 }
313#endif
314 }
315
316 G4ThreeVector polarizationVec = track.GetPolarization() ;
317 G4FieldTrack aFieldTrack = G4FieldTrack(startPosition,
318 track.GetGlobalTime(), // Lab.
319 track.GetMomentumDirection(),
320 track.GetKineticEnergy(),
321 restMass,
322 particleCharge,
323 polarizationVec,
324 pParticleDef->GetPDGMagneticMoment(),
325 0.0, // Length along track
326 pParticleDef->GetPDGSpin() ) ;
327 G4int stepNo= track.GetCurrentStepNumber();
328
329 ELimited limitedStep;
330 G4FieldTrack endTrackState('a'); // Default values
331
332 fMassGeometryLimitedStep = false ; // default
333 fAnyGeometryLimitedStep = false ;
334 if( currentMinimumStep > 0 )
335 {
336 G4double newMassSafety= 0.0; // temp. for recalculation
337
338 // Do the Transport in the field (non recti-linear)
339 //
340 lengthAlongCurve = fPathFinder->ComputeStep( aFieldTrack,
341 currentMinimumStep,
342 fNavigatorId,
343 stepNo,
344 newMassSafety,
345 limitedStep,
346 endTrackState,
347 currentVolume ) ;
348
349 G4double newFullSafety= fPathFinder->GetCurrentSafety();
350 // this was estimated already in step above
351
352 if( limitedStep == kUnique || limitedStep == kSharedTransport )
353 {
354 fMassGeometryLimitedStep = true ;
355 }
356
357 fAnyGeometryLimitedStep = (fPathFinder->GetNumberGeometriesLimitingStep() != 0);
358
359#ifdef G4DEBUG_TRANSPORT
360 if( fMassGeometryLimitedStep && !fAnyGeometryLimitedStep )
361 {
362 std::ostringstream message;
363 message << " ERROR in determining geometries limiting the step" << G4endl;
364 message << " Limiting: mass=" << fMassGeometryLimitedStep
365 << " any= " << fAnyGeometryLimitedStep << G4endl;
366 message << "Incompatible conditions - by which geometries was it limited ?";
367 G4Exception("G4CoupledTransportation::AlongStepGetPhysicalInteractionLength()",
368 "PathFinderConfused", FatalException, message);
369 }
370#endif
371
372 geometryStepLength = std::min( lengthAlongCurve, currentMinimumStep);
373
374 // Momentum: Magnitude and direction can be changed too now ...
375 //
376 fMomentumChanged = true ;
377 fTransportEndMomentumDir = endTrackState.GetMomentumDir() ;
378
379 // Remember last safety origin & value.
380 fPreviousSftOrigin = startPosition ;
381 fPreviousMassSafety = newMassSafety ;
382 fPreviousFullSafety = newFullSafety ;
383 // fpSafetyHelper->SetCurrentSafety( newFullSafety, startPosition);
384
385#ifdef G4DEBUG_TRANSPORT
386 if( verboseLevel > 1 )
387 {
388 G4cout << "G4Transport:CompStep> "
389 << " called the pathfinder for a new step at " << startPosition
390 << " and obtained step = " << lengthAlongCurve << G4endl;
391 G4cout << " New safety (preStep) = " << newMassSafety
392 << " versus precalculated = " << startMassSafety << G4endl;
393 }
394#endif
395
396 // Store as best estimate value
397 startMassSafety = newMassSafety ;
398 startFullSafety = newFullSafety ;
399
400 // Get the End-Position and End-Momentum (Dir-ection)
401 fTransportEndPosition = endTrackState.GetPosition() ;
402 fTransportEndKineticEnergy = endTrackState.GetKineticEnergy() ;
403 }
404 else
405 {
406 geometryStepLength = lengthAlongCurve= 0.0 ;
407 fMomentumChanged = false ;
408 // fMassGeometryLimitedStep = false ; // --- ???
409 // fAnyGeometryLimitedStep = true;
410 fTransportEndMomentumDir = track.GetMomentumDirection();
411 fTransportEndKineticEnergy = track.GetKineticEnergy();
412
413 fTransportEndPosition = startPosition;
414
415 endTrackState= aFieldTrack; // Ensures that time is updated
416
417 // If the step length requested is 0, and we are on a boundary
418 // then a boundary will also limit the step.
419 if( startMassSafety == 0.0 )
420 {
421 fMassGeometryLimitedStep = true ;
422 fAnyGeometryLimitedStep = true;
423 }
424 // TODO: Add explicit logical status for being at a boundary
425 }
426 // G4FieldTrack aTrackState(endTrackState);
427
428 if( !fieldExertsForce )
429 {
430 fParticleIsLooping = false ;
431 fMomentumChanged = false ;
432 fEndGlobalTimeComputed = false ;
433 }
434 else
435 {
436 fParticleIsLooping = fFieldPropagator->IsParticleLooping() ;
437
438#ifdef G4DEBUG_TRANSPORT
439 if( verboseLevel > 1 )
440 {
441 G4cout << " G4CT::CS End Position = "
442 << fTransportEndPosition << G4endl;
443 G4cout << " G4CT::CS End Direction = "
444 << fTransportEndMomentumDir << G4endl;
445 }
446#endif
447 if( fFieldPropagator->GetCurrentFieldManager()->DoesFieldChangeEnergy() )
448 {
449 // If the field can change energy, then the time must be integrated
450 // - so this should have been updated
451 //
452 fCandidateEndGlobalTime = endTrackState.GetLabTimeOfFlight();
453 fEndGlobalTimeComputed = true;
454
455 // was ( fCandidateEndGlobalTime != track.GetGlobalTime() );
456 // a cleaner way is to have FieldTrack knowing whether time
457 // is updated
458 }
459 else
460 {
461 // The energy should be unchanged by field transport,
462 // - so the time changed will be calculated elsewhere
463 //
464 fEndGlobalTimeComputed = false;
465
466 // Check that the integration preserved the energy
467 // - and if not correct this!
468 G4double startEnergy= track.GetKineticEnergy();
469 G4double endEnergy= fTransportEndKineticEnergy;
470
471 static G4ThreadLocal G4int no_inexact_steps=0; // , no_large_ediff;
472 G4double absEdiff = std::fabs(startEnergy- endEnergy);
473 if( absEdiff > perMillion * endEnergy )
474 {
475 no_inexact_steps++;
476 // Possible statistics keeping here ...
477 }
478#ifdef G4VERBOSE
479 if( (verboseLevel > 1) && ( absEdiff > perThousand * endEnergy) )
480 {
481 ReportInexactEnergy(startEnergy, endEnergy);
482 } // end of if (verboseLevel)
483#endif
484 // Correct the energy for fields that conserve it
485 // This - hides the integration error
486 // - but gives a better physical answer
487 fTransportEndKineticEnergy= track.GetKineticEnergy();
488 }
489 }
490
491 fEndpointDistance = (fTransportEndPosition - startPosition).mag() ;
492 fTransportEndSpin = endTrackState.GetSpin();
493
494 // Calculate the safety
495
496 safetyProposal= startFullSafety; // used to be startMassSafety
497 // Changed to accomodate processes that cannot update the safety
498
499 // Update safety for the end-point, if becomes negative at the end-point.
500
501 if( (startFullSafety < fEndpointDistance )
502 && ( particleCharge != 0.0 ) ) // Only needed to prepare for MSC
503 // && !fAnyGeometryLimitedStep ) // To-Try: No safety update if at boundary
504 {
505 G4double endFullSafety =
506 fPathFinder->ComputeSafety( fTransportEndPosition);
507 // Expected mission -- only mass geometry's safety
508 // fMassNavigator->ComputeSafety( fTransportEndPosition) ;
509 // Yet discrete processes only have poststep -- and this cannot
510 // currently revise the safety
511 // ==> so we use the all-geometry safety as a precaution
512
513 fpSafetyHelper->SetCurrentSafety( endFullSafety, fTransportEndPosition);
514 // Pushing safety to Helper avoids recalculation at this point
515
516 G4ThreeVector centerPt= G4ThreeVector(0.0, 0.0, 0.0); // Used for return value
517 G4double endMassSafety= fPathFinder->ObtainSafety( fNavigatorId, centerPt);
518 // Retrieves the mass value from PathFinder (it calculated it)
519
520 fPreviousMassSafety = endMassSafety ;
521 fPreviousFullSafety = endFullSafety;
522 fPreviousSftOrigin = fTransportEndPosition ;
523
524 // The convention (Stepping Manager's) is safety from the start point
525 //
526 safetyProposal = endFullSafety + fEndpointDistance;
527 // --> was endMassSafety
528 // Changed to accomodate processes that cannot update the safety
529
530#ifdef G4DEBUG_TRANSPORT
531 G4int prec= G4cout.precision(12) ;
532 G4cout << "***CoupledTransportation::AlongStepGPIL ** " << G4endl ;
533 G4cout << " Revised Safety at endpoint " << fTransportEndPosition
534 << " give safety values: Mass= " << endMassSafety
535 << " All= " << endFullSafety << G4endl ;
536 G4cout << " Adding endpoint distance " << fEndpointDistance
537 << " to obtain pseudo-safety= " << safetyProposal << G4endl ;
538 G4cout.precision(prec);
539 }
540 else
541 {
542 G4int prec= G4cout.precision(12) ;
543 G4cout << "***CoupledTransportation::AlongStepGPIL ** " << G4endl ;
544 G4cout << " Quick Safety estimate at endpoint "
545 << fTransportEndPosition
546 << " gives safety endpoint value = "
547 << startFullSafety - fEndpointDistance
548 << " using start-point value " << startFullSafety
549 << " and endpointDistance " << fEndpointDistance << G4endl;
550 G4cout.precision(prec);
551#endif
552 }
553
554 proposedSafetyForStart= safetyProposal;
555 fParticleChange.ProposeTrueStepLength(geometryStepLength) ;
556
557 return geometryStepLength ;
558}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
@ CandidateForSelection
#define fPreviousSftOrigin
ELimited
@ kUnique
@ kSharedTransport
CLHEP::Hep3Vector G4ThreeVector
G4GLOB_DLL std::ostream G4cerr
double mag2() const
void ReportInexactEnergy(G4double startEnergy, G4double endEnergy)
const G4ThreeVector & GetMomentumDirection() const
G4double GetCharge() const
G4ParticleDefinition * GetDefinition() const
G4double GetMagneticMoment() const
G4double GetTotalMomentum() const
G4bool DoesFieldChangeEnergy() const
virtual void ConfigureForTrack(const G4Track *)
const G4Field * GetDetectorField() const
G4bool IsGravityActive() const
Definition: G4Field.hh:101
G4double GetPDGMagneticMoment() const
G4double ComputeSafety(const G4ThreeVector &globalPoint)
unsigned int GetNumberGeometriesLimitingStep() const
G4double ObtainSafety(G4int navId, G4ThreeVector &globalCenterPoint)
G4double ComputeStep(const G4FieldTrack &pFieldTrack, G4double pCurrentProposedStepLength, G4int navigatorId, G4int stepNo, G4double &pNewSafety, ELimited &limitedStep, G4FieldTrack &EndState, G4VPhysicalVolume *currentVolume)
G4double GetCurrentSafety() const
G4FieldManager * FindAndSetFieldManager(G4VPhysicalVolume *pCurrentPhysVol)
G4FieldManager * GetCurrentFieldManager()
G4bool IsParticleLooping() const
G4EquationOfMotion * GetCurrentEquationOfMotion()
void SetCurrentSafety(G4double val, const G4ThreeVector &pos)
G4VPhysicalVolume * GetVolume() const
const G4ThreeVector & GetPosition() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
const G4ThreeVector & GetPolarization() const
void ProposeTrueStepLength(G4double truePathLength)
const G4String & GetName() const
T sqr(const T &x)
Definition: templates.hh:128

◆ AtRestDoIt()

G4VParticleChange * G4CoupledTransportation::AtRestDoIt ( const G4Track ,
const G4Step  
)
inlinevirtual

Implements G4VProcess.

Definition at line 167 of file G4CoupledTransportation.hh.

168 { return 0; } // No operation in AtRestDoIt

◆ AtRestGetPhysicalInteractionLength()

G4double G4CoupledTransportation::AtRestGetPhysicalInteractionLength ( const G4Track ,
G4ForceCondition  
)
inlinevirtual

Implements G4VProcess.

Definition at line 163 of file G4CoupledTransportation.hh.

165 { return -1.0; } // No operation in AtRestGPIL

◆ DoesAnyFieldExist()

G4bool G4CoupledTransportation::DoesAnyFieldExist ( )
protected

◆ EnableGravity()

G4bool G4CoupledTransportation::EnableGravity ( G4bool  useGravity)
static

Definition at line 1076 of file G4CoupledTransportation.cc.

1077{
1078 G4bool lastValue= fUseGravity;
1079 fUseGravity= useGravity;
1080 G4Transportation::fUseGravity= useGravity;
1081 return lastValue;
1082}

◆ EnableMagneticMoment()

G4bool G4CoupledTransportation::EnableMagneticMoment ( G4bool  useMoment = true)
static

Definition at line 1065 of file G4CoupledTransportation.cc.

1066{
1067 G4bool lastValue= fUseMagneticMoment;
1068 fUseMagneticMoment= useMoment;
1069 G4Transportation::fUseMagneticMoment= useMoment;
1070 return lastValue;
1071}

Referenced by EnableUseMagneticMoment().

◆ EnableUseMagneticMoment()

static G4bool G4CoupledTransportation::EnableUseMagneticMoment ( G4bool  useMoment = true)
inlinestatic

Definition at line 159 of file G4CoupledTransportation.hh.

160 { return EnableMagneticMoment(useMoment); }
static G4bool EnableMagneticMoment(G4bool useMoment=true)

◆ EndTracking()

void G4CoupledTransportation::EndTracking ( )
virtual

Reimplemented from G4VProcess.

Definition at line 1006 of file G4CoupledTransportation.cc.

1007{
1009 fPathFinder->EndTrack();
1010 // Resets TransportationManager to use ordinary Navigator
1011}

◆ GetMaxEnergyKilled()

G4double G4CoupledTransportation::GetMaxEnergyKilled ( ) const
inline

◆ GetPropagatorInField()

G4PropagatorInField * G4CoupledTransportation::GetPropagatorInField ( )

◆ GetSignifyStepsInAnyVolume()

static G4bool G4CoupledTransportation::GetSignifyStepsInAnyVolume ( )
inlinestatic

Definition at line 140 of file G4CoupledTransportation.hh.

141 { return fSignifyStepInAnyVolume; }

◆ GetSilenceLooperWarnings()

G4bool G4CoupledTransportation::GetSilenceLooperWarnings ( )
static

Definition at line 1096 of file G4CoupledTransportation.cc.

1097{
1098 return fSilenceLooperWarnings;
1099}

◆ GetSumEnergyKilled()

G4double G4CoupledTransportation::GetSumEnergyKilled ( ) const
inline

◆ GetThresholdImportantEnergy()

G4double G4CoupledTransportation::GetThresholdImportantEnergy ( ) const
inline

◆ GetThresholdTrials()

G4int G4CoupledTransportation::GetThresholdTrials ( ) const
inline

◆ GetThresholdWarningEnergy()

G4double G4CoupledTransportation::GetThresholdWarningEnergy ( ) const
inline

◆ IsFirstStepInAnyVolume()

G4bool G4CoupledTransportation::IsFirstStepInAnyVolume ( ) const
inline

Definition at line 149 of file G4CoupledTransportation.hh.

149{ return fFirstStepInAnyVolume; }

◆ IsFirstStepInMassVolume()

G4bool G4CoupledTransportation::IsFirstStepInMassVolume ( ) const
inline

Definition at line 151 of file G4CoupledTransportation.hh.

151{ return fFirstStepInMassVolume; }

◆ IsLastStepInAnyVolume()

G4bool G4CoupledTransportation::IsLastStepInAnyVolume ( ) const
inline

Definition at line 150 of file G4CoupledTransportation.hh.

150{ return fAnyGeometryLimitedStep; }

◆ IsLastStepInMassVolume()

G4bool G4CoupledTransportation::IsLastStepInMassVolume ( ) const
inline

Definition at line 152 of file G4CoupledTransportation.hh.

152{ return fMassGeometryLimitedStep; }

◆ PostStepDoIt()

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

Implements G4VProcess.

Definition at line 764 of file G4CoupledTransportation.cc.

766{
767 G4TouchableHandle retCurrentTouchable ; // The one to return
768
769 // Initialize ParticleChange (by setting all its members equal
770 // to corresponding members in G4Track)
771 // fParticleChange.Initialize(track) ; // To initialise TouchableChange
772
773 fParticleChange.ProposeTrackStatus(track.GetTrackStatus()) ;
774
775 if( fSignifyStepInAnyVolume )
776 {
777 fParticleChange.ProposeFirstStepInVolume( fFirstStepInAnyVolume );
778 }
779 else
780 {
781 fParticleChange.ProposeFirstStepInVolume( fFirstStepInMassVolume );
782 }
783
784 // Check that the end position and direction are preserved
785 // since call to AlongStepDoIt
786
787#ifdef G4DEBUG_TRANSPORT
788 if( ( verboseLevel > 0 )
789 && ((fTransportEndPosition - track.GetPosition()).mag2() >= 1.0e-16) )
790 {
791 ReportMove( track.GetPosition(), fTransportEndPosition,
792 "End of Step Position" );
793 G4cerr << " Problem in G4CoupledTransportation::PostStepDoIt " << G4endl;
794 }
795
796 // If the Step was determined by the volume boundary, relocate the particle
797 // The pathFinder will know that the geometry limited the step (!?)
798
799 if( verboseLevel > 0 )
800 {
801 G4cout << " Calling PathFinder::Locate() from "
802 << " G4CoupledTransportation::PostStepDoIt " << G4endl;
803 G4cout << " fAnyGeometryLimitedStep is " << fAnyGeometryLimitedStep
804 << G4endl;
805
806 }
807#endif
808
809 if(fAnyGeometryLimitedStep)
810 {
811 fPathFinder->Locate( track.GetPosition(),
812 track.GetMomentumDirection(),
813 true);
814
815 // fCurrentTouchable will now become the previous touchable,
816 // and what was the previous will be freed.
817 // (Needed because the preStepPoint can point to the previous touchable)
818
819 fCurrentTouchableHandle=
820 fPathFinder->CreateTouchableHandle( fNavigatorId );
821
822#ifdef G4DEBUG_TRANSPORT
823 if( verboseLevel > 0 )
824 {
825 G4cout << "G4CoupledTransportation::PostStepDoIt --- fNavigatorId = "
826 << fNavigatorId << G4endl;
827 }
828 if( verboseLevel > 1 )
829 {
830 G4VPhysicalVolume* vol= fCurrentTouchableHandle->GetVolume();
831 G4cout << "CHECK !!!!!!!!!!! fCurrentTouchableHandle->GetVolume() = "
832 << vol;
833 if( vol ) { G4cout << "Name=" << vol->GetName(); }
834 G4cout << G4endl;
835 }
836#endif
837
838 // Check whether the particle is out of the world volume
839 // If so it has exited and must be killed.
840 //
841 if( fCurrentTouchableHandle->GetVolume() == 0 )
842 {
843 fParticleChange.ProposeTrackStatus( fStopAndKill ) ;
844 }
845 retCurrentTouchable = fCurrentTouchableHandle ;
846 // fParticleChange.SetTouchableHandle( fCurrentTouchableHandle ) ;
847 }
848 else // fAnyGeometryLimitedStep is false
849 {
850#ifdef G4DEBUG_TRANSPORT
851 if( verboseLevel > 1 )
852 {
853 G4cout << "G4CoupledTransportation::PostStepDoIt -- "
854 << " fAnyGeometryLimitedStep = " << fAnyGeometryLimitedStep
855 << " must be false " << G4endl;
856 }
857#endif
858 // This serves only to move each of the Navigator's location
859 //
860 // fLinearNavigator->LocateGlobalPointWithinVolume( track.GetPosition() ) ;
861
862 fPathFinder->ReLocate( track.GetPosition() );
863 // track.GetMomentumDirection() );
864
865 // Keep the value of the track's current Touchable is retained,
866 // and use it to overwrite the (unset) one in particle change.
867 // Expect this must be fCurrentTouchable too
868 // - could it be different, eg at the start of a step ?
869 //
870 retCurrentTouchable = track.GetTouchableHandle() ;
871 // fParticleChange.SetTouchableHandle( track.GetTouchableHandle() ) ;
872 } // endif ( fAnyGeometryLimitedStep )
873
874#ifdef G4DEBUG_NAVIGATION
875 G4cout << " CoupledTransport::AlongStep GPIL: "
876 << " last-step: any= " << fAnyGeometryLimitedStep
877 << " . ..... x . "
878 << " mass= " << fMassGeometryLimitedStep
879 << G4endl;
880#endif
881
882 if( fSignifyStepInAnyVolume )
883 fParticleChange.ProposeLastStepInVolume(fAnyGeometryLimitedStep);
884 else
885 fParticleChange.ProposeLastStepInVolume(fMassGeometryLimitedStep);
886
887 const G4VPhysicalVolume* pNewVol = retCurrentTouchable->GetVolume() ;
888 const G4Material* pNewMaterial = 0 ;
889 const G4VSensitiveDetector* pNewSensitiveDetector = 0 ;
890
891 if( pNewVol != 0 )
892 {
893 pNewMaterial= pNewVol->GetLogicalVolume()->GetMaterial();
894 pNewSensitiveDetector= pNewVol->GetLogicalVolume()->GetSensitiveDetector();
895 }
896
897 // ( const_cast<G4Material *> pNewMaterial ) ;
898 // ( const_cast<G4VSensitiveDetetor *> pNewSensitiveDetector) ;
899
900 fParticleChange.SetMaterialInTouchable( (G4Material *) pNewMaterial ) ;
901 fParticleChange.SetSensitiveDetectorInTouchable( (G4VSensitiveDetector *) pNewSensitiveDetector ) ;
902 // "temporarily" until Get/Set Material of ParticleChange,
903 // and StepPoint can be made const.
904
905 const G4MaterialCutsCouple* pNewMaterialCutsCouple = 0;
906 if( pNewVol != 0 )
907 {
908 pNewMaterialCutsCouple=pNewVol->GetLogicalVolume()->GetMaterialCutsCouple();
909 if( pNewMaterialCutsCouple!=0
910 && pNewMaterialCutsCouple->GetMaterial()!=pNewMaterial )
911 {
912 // for parametrized volume
913 //
914 pNewMaterialCutsCouple = G4ProductionCutsTable::GetProductionCutsTable()
915 ->GetMaterialCutsCouple(pNewMaterial,
916 pNewMaterialCutsCouple->GetProductionCuts());
917 }
918 }
919 fParticleChange.SetMaterialCutsCoupleInTouchable( pNewMaterialCutsCouple );
920
921 // Must always set the touchable in ParticleChange, whether relocated or not
922 //
923 fParticleChange.SetTouchableHandle(retCurrentTouchable) ;
924
925 return &fParticleChange ;
926}
void ReportMove(G4ThreeVector OldVector, G4ThreeVector NewVector, const G4String &Quantity)
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)
void ReLocate(const G4ThreeVector &position)
void Locate(const G4ThreeVector &position, const G4ThreeVector &direction, G4bool relativeSearch=true)
G4TouchableHandle CreateTouchableHandle(G4int navId) const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
static G4ProductionCutsTable * GetProductionCutsTable()
G4TrackStatus GetTrackStatus() const
const G4TouchableHandle & GetTouchableHandle() const
void ProposeLastStepInVolume(G4bool flag)
void ProposeFirstStepInVolume(G4bool flag)
G4LogicalVolume * GetLogicalVolume() const
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:41

◆ PostStepGetPhysicalInteractionLength()

G4double G4CoupledTransportation::PostStepGetPhysicalInteractionLength ( const G4Track ,
G4double  previousStepSize,
G4ForceCondition pForceCond 
)
virtual

Implements G4VProcess.

Definition at line 732 of file G4CoupledTransportation.cc.

736{
737 // Must act as PostStep action -- to relocate particle
738 *pForceCond = Forced ;
739 return DBL_MAX ;
740}
@ Forced

◆ PrintStatistics()

void G4CoupledTransportation::PrintStatistics ( std::ostream &  outStr) const

Definition at line 137 of file G4CoupledTransportation.cc.

138{
139 if( fSumEnergyKilled > 0.0 )
140 {
141 outStr << " G4CoupledTransportation: Statistics for looping particles "
142 << G4endl;
143 outStr << " Sum of energy of loopers killed: "
144 << fSumEnergyKilled / MeV << " MeV " << G4endl;
145 outStr << " Max energy of loopers killed: "
146 << fMaxEnergyKilled / MeV << " MeV " << G4endl;
147
148
149 outStr << " Max energy of loopers 'saved': " << fMaxEnergySaved << G4endl;
150 outStr << " Sum of energy of loopers 'saved': "
151 << fSumEnergySaved << G4endl;
152 outStr << " Sum of energy of unstable loopers 'saved': "
153 << fSumEnergyUnstableSaved << G4endl;
154 }
155}

Referenced by ~G4CoupledTransportation().

◆ PushThresholdsToLogger()

void G4CoupledTransportation::PushThresholdsToLogger ( )

◆ ReportInexactEnergy()

void G4CoupledTransportation::ReportInexactEnergy ( G4double  startEnergy,
G4double  endEnergy 
)
protected

Definition at line 1016 of file G4CoupledTransportation.cc.

1018{
1019 static G4ThreadLocal G4int no_warnings= 0, warnModulo=1,
1020 moduloFactor= 10, no_large_ediff= 0;
1021
1022 if( std::fabs(startEnergy- endEnergy) > perThousand * endEnergy )
1023 {
1024 no_large_ediff ++;
1025 if( (no_large_ediff% warnModulo) == 0 )
1026 {
1027 no_warnings++;
1028 std::ostringstream message;
1029 message << "Energy change in Step is above 1^-3 relative value. "
1030 << G4endl
1031 << " Relative change in 'tracking' step = "
1032 << std::setw(15) << (endEnergy-startEnergy)/startEnergy
1033 << G4endl
1034 << " Starting E= " << std::setw(12) << startEnergy / MeV
1035 << " MeV " << G4endl
1036 << " Ending E= " << std::setw(12) << endEnergy / MeV
1037 << " MeV " << G4endl
1038 << "Energy has been corrected -- however, review"
1039 << " field propagation parameters for accuracy." << G4endl;
1040 if ( (verboseLevel > 2 ) || (no_warnings<4)
1041 || (no_large_ediff == warnModulo * moduloFactor) )
1042 {
1043 message << "These include EpsilonStepMax(/Min) in G4FieldManager,"
1044 << G4endl
1045 << "which determine fractional error per step for integrated quantities."
1046 << G4endl
1047 << "Note also the influence of the permitted number of integration steps."
1048 << G4endl;
1049 }
1050 message << "Bad 'endpoint'. Energy change detected and corrected."
1051 << G4endl
1052 << "Has occurred already " << no_large_ediff << " times.";
1053 G4Exception("G4CoupledTransportation::AlongStepGetPIL()",
1054 "EnergyChange", JustWarning, message);
1055 if( no_large_ediff == warnModulo * moduloFactor )
1056 {
1057 warnModulo *= moduloFactor;
1058 }
1059 }
1060 }
1061}
@ JustWarning

Referenced by AlongStepGetPhysicalInteractionLength().

◆ ReportLooperThresholds()

void G4CoupledTransportation::ReportLooperThresholds ( )

Definition at line 1144 of file G4CoupledTransportation.cc.

1145{
1146 PushThresholdsToLogger(); // To be absolutely certain they are in sync
1147 fpLogger->ReportLooperThresholds("G4CoupledTransportation");
1148}
void ReportLooperThresholds(const char *className)

Referenced by ReportMissingLogger(), SetHighLooperThresholds(), and SetLowLooperThresholds().

◆ ReportMissingLogger()

void G4CoupledTransportation::ReportMissingLogger ( const char *  methodName)
protected

Definition at line 1132 of file G4CoupledTransportation.cc.

1133{
1134 const char* message= "Logger object missing from G4CoupledTransportation";
1135 G4String classAndMethod= G4String("G4CoupledTransportation") + G4String( methodName );
1136 G4Exception(classAndMethod, "Missing Logger", JustWarning, message);
1137
1139}

◆ ReportMove()

void G4CoupledTransportation::ReportMove ( G4ThreeVector  OldVector,
G4ThreeVector  NewVector,
const G4String Quantity 
)
protected

Definition at line 744 of file G4CoupledTransportation.cc.

747{
748 G4ThreeVector moveVec = ( NewVector - OldVector );
749
750 G4cerr << G4endl
751 << "**************************************************************"
752 << G4endl;
753 G4cerr << "Endpoint has moved between value expected from TransportEndPosition "
754 << " and value from Track in PostStepDoIt. " << G4endl
755 << "Change of " << Quantity << " is " << moveVec.mag() / mm
756 << " mm long, "
757 << " and its vector is " << (1.0/mm) * moveVec << " mm " << G4endl
758 << "Endpoint of ComputeStep was " << OldVector
759 << " and current position to locate is " << NewVector << G4endl;
760}
double mag() const

Referenced by PostStepDoIt().

◆ ResetKilledStatistics()

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

◆ SetHighLooperThresholds()

void G4CoupledTransportation::SetHighLooperThresholds ( )

Definition at line 1103 of file G4CoupledTransportation.cc.

1104{
1105 // Setting old 'high' values for thresholds - old values, potentially appropriate
1106 // for the energy frontier HEP experiments
1107 SetThresholdWarningEnergy( 100.0 * CLHEP::MeV ); // Warn above this energy
1108 SetThresholdImportantEnergy( 250.0 * CLHEP::MeV ); // Give a few trial above this E);
1109
1110 G4int maxTrials = 10;
1111 SetThresholdTrials( maxTrials );
1112
1114}
void SetThresholdTrials(G4int newMaxTrials)
void SetThresholdImportantEnergy(G4double newEnImp)
void SetThresholdWarningEnergy(G4double newEnWarn)

Referenced by G4CoupledTransportation().

◆ SetLowLooperThresholds()

void G4CoupledTransportation::SetLowLooperThresholds ( )

Definition at line 1117 of file G4CoupledTransportation.cc.

1118{
1119 // These values were the default in Geant4 10.5 - beta
1120 SetThresholdWarningEnergy( 1.0 * CLHEP::keV ); // Warn above this En
1121 SetThresholdImportantEnergy( 1.0 * CLHEP::MeV ); // Extra trials above it
1122
1123 G4int maxTrials = 30; // A new value - was 10
1124 SetThresholdTrials( maxTrials );
1125
1127}

◆ SetPropagatorInField()

void G4CoupledTransportation::SetPropagatorInField ( G4PropagatorInField pFieldPropagator)

◆ SetSignifyStepsInAnyVolume()

static void G4CoupledTransportation::SetSignifyStepsInAnyVolume ( G4bool  anyVol)
inlinestatic

Definition at line 138 of file G4CoupledTransportation.hh.

139 { fSignifyStepInAnyVolume = anyVol; }

◆ SetSilenceLooperWarnings()

void G4CoupledTransportation::SetSilenceLooperWarnings ( G4bool  val)
static

Definition at line 1088 of file G4CoupledTransportation.cc.

1089{
1090 fSilenceLooperWarnings= val; // Flag to *Supress* all 'looper' warnings
1091 // G4CoupledTransportation::fSilenceLooperWarnings= val;
1092}

◆ SetThresholdImportantEnergy()

void G4CoupledTransportation::SetThresholdImportantEnergy ( G4double  newEnImp)
inline

◆ SetThresholdTrials()

void G4CoupledTransportation::SetThresholdTrials ( G4int  newMaxTrials)
inline

◆ SetThresholdWarningEnergy()

void G4CoupledTransportation::SetThresholdWarningEnergy ( G4double  newEnWarn)
inline

◆ StartTracking()

void G4CoupledTransportation::StartTracking ( G4Track aTrack)
virtual

Reimplemented from G4VProcess.

Definition at line 935 of file G4CoupledTransportation.cc.

936{
937
938 G4TransportationManager* transportMgr =
940
941 // G4VProcess::StartTracking(aTrack);
942 fNewTrack= true;
943
944 // The 'initialising' actions
945 // once taken in AlongStepGPIL -- if ( track.GetCurrentStepNumber()==1 )
946
947 // fStartedNewTrack= true;
948
949 fMassNavigator = transportMgr->GetNavigatorForTracking() ;
950 fNavigatorId= transportMgr->ActivateNavigator( fMassNavigator );
951
953 G4ThreeVector direction = aTrack->GetMomentumDirection();
954
955 fPathFinder->PrepareNewTrack( position, direction);
956 // This implies a call to fPathFinder->Locate( position, direction );
957
958 // Whether field exists should be determined at run level -- TODO
959 fAnyFieldExists= DoesAnyFieldExist();
960
961 // reset safety value and center
962 //
963 fPreviousMassSafety = 0.0 ;
964 fPreviousFullSafety = 0.0 ;
965 fPreviousSftOrigin = G4ThreeVector(0.,0.,0.) ;
966
967 // reset looping counter -- for motion in field
968 fNoLooperTrials= 0;
969
970 // Must clear this state .. else it depends on last track's value
971 // --> a better solution would set this from state of suspended track TODO ?
972 // Was if( aTrack->GetCurrentStepNumber()==1 ) { .. }
973
974 // ChordFinder reset internal state
975 //
976 if( fFieldPropagator && fAnyFieldExists )
977 {
978 fFieldPropagator->ClearPropagatorState();
979 // Resets safety values, in case of overlaps.
980
981 G4ChordFinder* chordF= fFieldPropagator->GetChordFinder();
982 if( chordF ) { chordF->ResetStepEstimate(); }
983 }
984
985 // Clear the chord finders of all fields (ie managers) derived objects
986 //
988 fieldMgrStore->ClearAllChordFindersState();
989
990#ifdef G4DEBUG_TRANSPORT
991 if( verboseLevel > 1 )
992 {
993 G4cout << " Returning touchable handle " << fCurrentTouchableHandle
994 << G4endl;
995 }
996#endif
997
998 // Update the current touchable handle (from the track's)
999 //
1000 fCurrentTouchableHandle = aTrack->GetTouchableHandle();
1001}
void ResetStepEstimate()
static G4FieldManagerStore * GetInstance()
void PrepareNewTrack(const G4ThreeVector &position, const G4ThreeVector &direction, G4VPhysicalVolume *massStartVol=nullptr)
G4ChordFinder * GetChordFinder()

Friends And Related Function Documentation

◆ G4Transportation

friend class G4Transportation
friend

Definition at line 270 of file G4CoupledTransportation.hh.


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