Geant4 11.1.1
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)
 
G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)
 
G4bool IsFirstStepInAnyVolume () const
 
G4bool IsLastStepInAnyVolume () const
 
G4bool IsFirstStepInMassVolume () const
 
G4bool IsLastStepInMassVolume () const
 
void StartTracking (G4Track *aTrack)
 
void EndTracking ()
 
- Public Member Functions inherited from G4Transportation
 G4Transportation (G4int verbosityLevel=1, const G4String &aName="Transportation")
 
 ~G4Transportation ()
 
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)
 
G4bool FieldExertedForce ()
 
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)
 
void EnableShortStepOptimisation (G4bool optimise=true)
 
G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
void StartTracking (G4Track *aTrack)
 
virtual void ProcessDescription (std::ostream &outFile) const
 
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 const G4VProcessGetCreatorProcess () const
 
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 void SetSignifyStepsInAnyVolume (G4bool anyVol)
 
static G4bool GetSignifyStepsInAnyVolume ()
 
static G4bool EnableUseMagneticMoment (G4bool useMoment=true)
 
- Static Public Member Functions inherited from G4Transportation
static G4bool EnableMagneticMoment (G4bool useMoment=true)
 
static G4bool EnableGravity (G4bool useGravity=true)
 
static void SetSilenceLooperWarnings (G4bool val)
 
static G4bool GetSilenceLooperWarnings ()
 
static G4bool EnableUseMagneticMoment (G4bool useMoment=true)
 
- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 

Protected Member Functions

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

Additional Inherited Members

- Protected Attributes inherited from G4Transportation
G4NavigatorfLinearNavigator
 
G4PropagatorInFieldfFieldPropagator
 
G4ThreeVector fTransportEndPosition = G4ThreeVector( 0.0, 0.0, 0.0 )
 
G4ThreeVector fTransportEndMomentumDir = G4ThreeVector( 0.0, 0.0, 0.0 )
 
G4double fTransportEndKineticEnergy = 0.0
 
G4ThreeVector fTransportEndSpin = G4ThreeVector( 0.0, 0.0, 0.0 )
 
G4bool fMomentumChanged = true
 
G4bool fEndGlobalTimeComputed = false
 
G4double fCandidateEndGlobalTime = 0.0
 
G4bool fAnyFieldExists = false
 
G4bool fParticleIsLooping = false
 
G4bool fNewTrack = true
 
G4bool fFirstStepInVolume = true
 
G4bool fLastStepInVolume = false
 
G4bool fGeometryLimitedStep = true
 
G4bool fFieldExertedForce = false
 
G4TouchableHandle fCurrentTouchableHandle
 
G4ThreeVector fPreviousSftOrigin
 
G4double fPreviousSafety
 
G4ParticleChangeForTransport fParticleChange
 
G4double fEndPointDistance
 
G4double fThreshold_Warning_Energy = 1.0 * CLHEP::keV
 
G4double fThreshold_Important_Energy = 1.0 * CLHEP::MeV
 
G4int fThresholdTrials = 10
 
G4int fAbandonUnstableTrials = 0
 
G4int fNoLooperTrials = 0
 
G4double fSumEnergyKilled = 0.0
 
G4double fSumEnerSqKilled = 0.0
 
G4double fMaxEnergyKilled = -1.0
 
G4int fMaxEnergyKilledPDG = 0
 
unsigned long fNumLoopersKilled = 0
 
G4double fSumEnergyKilled_NonElectron = 0.0
 
G4double fSumEnerSqKilled_NonElectron = 0.0
 
G4double fMaxEnergyKilled_NonElectron = -1.0
 
G4int fMaxEnergyKilled_NonElecPDG = 0
 
unsigned long fNumLoopersKilled_NonElectron = 0
 
G4double fSumEnergySaved = 0.0
 
G4double fMaxEnergySaved = -1.0
 
G4double fSumEnergyUnstableSaved = 0.0
 
G4bool fShortStepOptimisation
 
G4SafetyHelperfpSafetyHelper
 
G4TransportationLoggerfpLogger
 
- 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
 
- Static Protected Attributes inherited from G4Transportation
static G4bool fUseMagneticMoment =false
 
static G4bool fUseGravity = false
 
static G4bool fSilenceLooperWarnings = false
 

Detailed Description

Definition at line 56 of file G4CoupledTransportation.hh.

Constructor & Destructor Documentation

◆ G4CoupledTransportation()

G4CoupledTransportation::G4CoupledTransportation ( G4int  verbosityLevel = 0)

Definition at line 61 of file G4CoupledTransportation.cc.

62 : G4Transportation( verbosity, "CoupledTransportation" ),
63 fPreviousMassSafety( 0.0 ),
64 fPreviousFullSafety( 0.0 ),
65 fMassGeometryLimitedStep( false ),
66 fFirstStepInMassVolume( true )
67{
69 // SetVerboseLevel is called in the constructor of G4Transportation
70
71 if( verboseLevel > 0 )
72 {
73 G4cout << " G4CoupledTransportation constructor: ----- " << G4endl;
74 G4cout << " Verbose level is " << verboseLevel << G4endl;
75 G4cout << " Reports First/Last in "
76 << (fSignifyStepInAnyVolume ? " any " : " mass " )
77 << " geometry " << G4endl;
78 }
79 fPathFinder= G4PathFinder::GetInstance();
80}
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static G4PathFinder * GetInstance()
Definition: G4PathFinder.cc:52
G4int verboseLevel
Definition: G4VProcess.hh:360
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:410

◆ ~G4CoupledTransportation()

G4CoupledTransportation::~G4CoupledTransportation ( )

Definition at line 84 of file G4CoupledTransportation.cc.

85{
86}

Member Function Documentation

◆ AlongStepGetPhysicalInteractionLength()

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

Implements G4VProcess.

Definition at line 95 of file G4CoupledTransportation.cc.

101{
102 G4double geometryStepLength;
103 G4double startFullSafety= 0.0; // estimated safety for start point (all geometries)
104 G4double safetyProposal= -1.0; // local copy of proposal
105
106 G4ThreeVector EndUnitMomentum ;
107 G4double lengthAlongCurve = 0.0 ;
108
109 fParticleIsLooping = false ;
110
111 // Initial actions moved to StartTrack()
112 // --------------------------------------
113 // Note: in case another process changes touchable handle
114 // it will be necessary to add here (for all steps)
115 // fCurrentTouchableHandle = aTrack->GetTouchableHandle();
116
117 // GPILSelection is set to defaule value of CandidateForSelection
118 // It is a return value
119 //
120 *selection = CandidateForSelection ;
121
122 fFirstStepInMassVolume = fNewTrack || fMassGeometryLimitedStep ;
124
125#ifdef G4DEBUG_TRANSPORT
126 G4cout << " CoupledTransport::AlongStep GPIL: "
127 << " 1st-step: any= " <<fFirstStepInVolume << " ( geom= "
128 << fGeometryLimitedStep << " ) "
129 << " mass= " << fFirstStepInMassVolume << " ( geom= "
130 << fMassGeometryLimitedStep << " ) "
131 << " newTrack= " << fNewTrack << G4endl;
132#endif
133
134 // fLastStepInVolume= false;
135 fNewTrack = false;
136
137 // Get initial Energy/Momentum of the track
138 //
139 const G4DynamicParticle* pParticle = track.GetDynamicParticle() ;
140 const G4ParticleDefinition* pParticleDef = pParticle->GetDefinition() ;
141 G4ThreeVector startMomentumDir = pParticle->GetMomentumDirection() ;
142 G4ThreeVector startPosition = track.GetPosition() ;
143 G4VPhysicalVolume* currentVolume= track.GetVolume();
144
145#ifdef G4DEBUG_TRANSPORT
146 if( verboseLevel > 1 )
147 {
148 G4cout << "G4CoupledTransportation::AlongStepGPIL> called in volume "
149 << currentVolume->GetName() << G4endl;
150 }
151#endif
152 // G4double theTime = track.GetGlobalTime() ;
153
154 // The Step Point safety can be limited by other geometries and/or the
155 // assumptions of any process - it's not always the geometrical safety.
156 // We calculate the starting point's isotropic safety here.
157 //
158 G4ThreeVector OriginShift = startPosition - fPreviousSftOrigin ;
159 G4double MagSqShift = OriginShift.mag2() ;
160 startFullSafety= 0.0;
161
162 // Recall that FullSafety <= MassSafety
163 // Original: if( MagSqShift < sqr(fPreviousMassSafety) ) {
164 if( MagSqShift < sqr(fPreviousFullSafety) )
165 {
166 G4double mag_shift= std::sqrt(MagSqShift);
167 startFullSafety = std::max( (fPreviousFullSafety - mag_shift), 0.0);
168 // Need to be consistent between full safety with Mass safety
169 // in order reproduce results in simple case
170 // --> use same calculation method
171
172 // Only compute full safety if massSafety > 0. Else it remains 0
173 // startFullSafety = fPathFinder->ComputeSafety( startPosition );
174 }
175
176 // Is the particle charged or has it a magnetic moment?
177 //
178 G4double particleCharge = pParticle->GetCharge() ;
179 G4double magneticMoment = pParticle->GetMagneticMoment() ;
180 G4double restMass = pParticle->GetMass() ;
181
182 fMassGeometryLimitedStep = false ; // Set default - alex
183 fGeometryLimitedStep = false;
184
185 // There is no need to locate the current volume. It is Done elsewhere:
186 // On track construction
187 // By the tracking, after all AlongStepDoIts, in "Relocation"
188
189 // Check if the particle has a force, EM or gravitational, exerted on it
190 //
191 G4FieldManager* fieldMgr= nullptr;
192 G4bool fieldExertsForce = false ;
193
194 const G4Field* ptrField= nullptr;
195
197 G4bool eligibleEM = (particleCharge != 0.0)
198 || ( fUseMagneticMoment && (magneticMoment != 0.0) );
199 G4bool eligibleGrav = fUseGravity && (restMass != 0.0) ;
200
201 if( (fieldMgr!=nullptr) && (eligibleEM||eligibleGrav) )
202 {
203 // Message the field Manager, to configure it for this track
204 //
205 fieldMgr->ConfigureForTrack( &track );
206
207 // The above call can transition from a null field-ptr oto a finite field.
208 // If the field manager has no field ptr, the field is zero
209 // by definition ( = there is no field ! )
210 //
211 ptrField= fieldMgr->GetDetectorField();
212
213 if( ptrField != nullptr)
214 {
215 fieldExertsForce = eligibleEM
216 || ( eligibleGrav && ptrField->IsGravityActive() );
217 }
218 }
219 G4double momentumMagnitude = pParticle->GetTotalMomentum() ;
220
221 if( fieldExertsForce )
222 {
223 auto equationOfMotion= fFieldPropagator->GetCurrentEquationOfMotion();
224
225 G4ChargeState chargeState(particleCharge, // Charge can change (dynamic)
226 magneticMoment,
227 pParticleDef->GetPDGSpin() );
228 if( equationOfMotion )
229 {
230 equationOfMotion->SetChargeMomentumMass( chargeState,
231 momentumMagnitude,
232 restMass );
233 }
234#ifdef G4DEBUG_TRANSPORT
235 else
236 {
237 G4cerr << " ERROR in G4CoupledTransportation> "
238 << "Cannot find valid Equation of motion: " << G4endl;
239 << " Unable to pass Charge, Momentum and Mass " << G4endl;
240 }
241#endif
242 }
243
244 G4ThreeVector polarizationVec = track.GetPolarization() ;
245 G4FieldTrack aFieldTrack = G4FieldTrack(startPosition,
246 track.GetGlobalTime(), // Lab.
247 track.GetMomentumDirection(),
248 track.GetKineticEnergy(),
249 restMass,
250 particleCharge,
251 polarizationVec,
252 pParticleDef->GetPDGMagneticMoment(),
253 0.0, // Length along track
254 pParticleDef->GetPDGSpin() ) ;
255 G4int stepNo= track.GetCurrentStepNumber();
256
257 ELimited limitedStep;
258 G4FieldTrack endTrackState('a'); // Default values
259
260 fMassGeometryLimitedStep = false ; // default
261 fGeometryLimitedStep = false;
262 if( currentMinimumStep > 0 )
263 {
264 G4double newMassSafety= 0.0; // temp. for recalculation
265
266 // Do the Transport in the field (non recti-linear)
267 //
268 lengthAlongCurve = fPathFinder->ComputeStep( aFieldTrack,
269 currentMinimumStep,
271 stepNo,
272 newMassSafety,
273 limitedStep,
274 endTrackState,
275 currentVolume ) ;
276
277 G4double newFullSafety= fPathFinder->GetCurrentSafety();
278 // this was estimated already in step above
279
280 if( limitedStep == kUnique || limitedStep == kSharedTransport )
281 {
282 fMassGeometryLimitedStep = true ;
283 }
284
286
287#ifdef G4DEBUG_TRANSPORT
288 if( fMassGeometryLimitedStep && !fGeometryLimitedStep )
289 {
290 std::ostringstream message;
291 message << " ERROR in determining geometries limiting the step" << G4endl;
292 message << " Limiting: mass=" << fMassGeometryLimitedStep
293 << " any= " << fGeometryLimitedStep << G4endl;
294 message << "Incompatible conditions - by which geometries was it limited ?";
295 G4Exception("G4CoupledTransportation::AlongStepGetPhysicalInteractionLength()",
296 "PathFinderConfused", FatalException, message);
297 }
298#endif
299
300 geometryStepLength = std::min( lengthAlongCurve, currentMinimumStep);
301
302 // Momentum: Magnitude and direction can be changed too now ...
303 //
304 fMomentumChanged = true ;
305 fTransportEndMomentumDir = endTrackState.GetMomentumDir() ;
306
307 // Remember last safety origin & value.
308 fPreviousSftOrigin = startPosition ;
309 fPreviousMassSafety = newMassSafety ;
310 fPreviousFullSafety = newFullSafety ;
311 // fpSafetyHelper->SetCurrentSafety( newFullSafety, startPosition);
312
313#ifdef G4DEBUG_TRANSPORT
314 if( verboseLevel > 1 )
315 {
316 G4cout << "G4Transport:CompStep> "
317 << " called the pathfinder for a new step at " << startPosition
318 << " and obtained step = " << lengthAlongCurve << G4endl;
319 G4cout << " New safety (preStep) = " << newMassSafety << G4endl;
320 }
321#endif
322
323 // Store as best estimate value
324 startFullSafety = newFullSafety ;
325
326 // Get the End-Position and End-Momentum (Dir-ection)
327 fTransportEndPosition = endTrackState.GetPosition() ;
328 fTransportEndKineticEnergy = endTrackState.GetKineticEnergy() ;
329 }
330 else
331 {
332 geometryStepLength = lengthAlongCurve= 0.0 ;
333 fMomentumChanged = false ;
334 // fMassGeometryLimitedStep = false ; // --- ???
335 // fGeometryLimitedStep = true;
338
339 fTransportEndPosition = startPosition;
340
341 endTrackState= aFieldTrack; // Ensures that time is updated
342 }
343 // G4FieldTrack aTrackState(endTrackState);
344
345 if( !fieldExertsForce )
346 {
347 fParticleIsLooping = false ;
348 fMomentumChanged = false ;
349 fEndGlobalTimeComputed = false ;
350 }
351 else
352 {
354
355#ifdef G4DEBUG_TRANSPORT
356 if( verboseLevel > 1 )
357 {
358 G4cout << " G4CT::CS End Position = "
360 G4cout << " G4CT::CS End Direction = "
362 }
363#endif
365 {
366 // If the field can change energy, then the time must be integrated
367 // - so this should have been updated
368 //
369 fCandidateEndGlobalTime = endTrackState.GetLabTimeOfFlight();
371
372 // was ( fCandidateEndGlobalTime != track.GetGlobalTime() );
373 // a cleaner way is to have FieldTrack knowing whether time
374 // is updated
375 }
376 else
377 {
378 // The energy should be unchanged by field transport,
379 // - so the time changed will be calculated elsewhere
380 //
382
383#ifdef G4VERBOSE
384 // Check that the integration preserved the energy
385 // - and if not correct this!
386 G4double startEnergy= track.GetKineticEnergy();
388
389 G4double absEdiff = std::fabs(startEnergy- endEnergy);
390 if( (verboseLevel > 1) && ( absEdiff > perThousand * endEnergy) )
391 {
392 ReportInexactEnergy(startEnergy, endEnergy);
393 } // end of if (verboseLevel)
394#endif
395 // Correct the energy for fields that conserve it
396 // This - hides the integration error
397 // - but gives a better physical answer
399 }
400 }
401
402 fEndPointDistance = (fTransportEndPosition - startPosition).mag() ;
403 fTransportEndSpin = endTrackState.GetSpin();
404
405 // Calculate the safety
406
407 safetyProposal= startFullSafety; // used to be startMassSafety
408 // Changed to accomodate processes that cannot update the safety
409
410 // Update safety for the end-point, if becomes negative at the end-point.
411
412 if( (startFullSafety < fEndPointDistance )
413 && ( particleCharge != 0.0 ) ) // Only needed to prepare for MSC
414 // && !fGeometryLimitedStep ) // To-Try: No safety update if at boundary
415 {
416 G4double endFullSafety =
417 fPathFinder->ComputeSafety( fTransportEndPosition);
418 // Expected mission -- only mass geometry's safety
419 // fLinearNavigator->ComputeSafety( fTransportEndPosition) ;
420 // Yet discrete processes only have poststep -- and this cannot
421 // currently revise the safety
422 // ==> so we use the all-geometry safety as a precaution
423
425 // Pushing safety to Helper avoids recalculation at this point
426
427 G4ThreeVector centerPt= G4ThreeVector(0.0, 0.0, 0.0); // Used for return value
428 G4double endMassSafety= fPathFinder->ObtainSafety( G4TransportationManager::kMassNavigatorId, centerPt);
429 // Retrieves the mass value from PathFinder (it calculated it)
430
431 fPreviousMassSafety = endMassSafety ;
432 fPreviousFullSafety = endFullSafety;
434
435 // The convention (Stepping Manager's) is safety from the start point
436 //
437 safetyProposal = endFullSafety + fEndPointDistance;
438 // --> was endMassSafety
439 // Changed to accomodate processes that cannot update the safety
440
441#ifdef G4DEBUG_TRANSPORT
442 G4int prec= G4cout.precision(12) ;
443 G4cout << "***CoupledTransportation::AlongStepGPIL ** " << G4endl ;
444 G4cout << " Revised Safety at endpoint " << fTransportEndPosition
445 << " give safety values: Mass= " << endMassSafety
446 << " All= " << endFullSafety << G4endl ;
447 G4cout << " Adding endpoint distance " << fEndPointDistance
448 << " to obtain pseudo-safety= " << safetyProposal << G4endl ;
449 G4cout.precision(prec);
450 }
451 else
452 {
453 G4int prec= G4cout.precision(12) ;
454 G4cout << "***CoupledTransportation::AlongStepGPIL ** " << G4endl ;
455 G4cout << " Quick Safety estimate at endpoint "
457 << " gives safety endpoint value = "
458 << startFullSafety - fEndPointDistance
459 << " using start-point value " << startFullSafety
460 << " and endpointDistance " << fEndPointDistance << G4endl;
461 G4cout.precision(prec);
462#endif
463 }
464
465 proposedSafetyForStart= safetyProposal;
466 fParticleChange.ProposeTrueStepLength(geometryStepLength) ;
467
468 return geometryStepLength ;
469}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
@ CandidateForSelection
ELimited
@ kUnique
@ kSharedTransport
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
G4GLOB_DLL std::ostream G4cerr
double mag2() const
void ReportInexactEnergy(G4double startEnergy, G4double endEnergy)
G4double GetMass() const
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
G4double GetGlobalTime() const
G4int GetCurrentStepNumber() const
const G4DynamicParticle * GetDynamicParticle() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
const G4ThreeVector & GetPolarization() const
static constexpr G4int kMassNavigatorId
G4ThreeVector fPreviousSftOrigin
G4ThreeVector fTransportEndSpin
G4double fTransportEndKineticEnergy
G4PropagatorInField * fFieldPropagator
G4ThreeVector fTransportEndPosition
G4ParticleChangeForTransport fParticleChange
static G4bool fUseGravity
static G4bool fUseMagneticMoment
G4ThreeVector fTransportEndMomentumDir
G4SafetyHelper * fpSafetyHelper
G4double fCandidateEndGlobalTime
void ProposeTrueStepLength(G4double truePathLength)
const G4String & GetName() const
T sqr(const T &x)
Definition: templates.hh:128

◆ EnableUseMagneticMoment()

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

Definition at line 104 of file G4CoupledTransportation.hh.

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

◆ EndTracking()

void G4CoupledTransportation::EndTracking ( )
virtual

Reimplemented from G4VProcess.

Definition at line 639 of file G4CoupledTransportation.cc.

640{
642 fPathFinder->EndTrack();
643 // Resets TransportationManager to use ordinary Navigator
644}
static G4TransportationManager * GetTransportationManager()

◆ GetSignifyStepsInAnyVolume()

static G4bool G4CoupledTransportation::GetSignifyStepsInAnyVolume ( )
inlinestatic

Definition at line 85 of file G4CoupledTransportation.hh.

86 { return fSignifyStepInAnyVolume; }

◆ IsFirstStepInAnyVolume()

G4bool G4CoupledTransportation::IsFirstStepInAnyVolume ( ) const
inline

Definition at line 94 of file G4CoupledTransportation.hh.

94{ return fFirstStepInVolume; }

◆ IsFirstStepInMassVolume()

G4bool G4CoupledTransportation::IsFirstStepInMassVolume ( ) const
inline

Definition at line 96 of file G4CoupledTransportation.hh.

96{ return fFirstStepInMassVolume; }

◆ IsLastStepInAnyVolume()

G4bool G4CoupledTransportation::IsLastStepInAnyVolume ( ) const
inline

Definition at line 95 of file G4CoupledTransportation.hh.

95{ return fGeometryLimitedStep; }

◆ IsLastStepInMassVolume()

G4bool G4CoupledTransportation::IsLastStepInMassVolume ( ) const
inline

Definition at line 97 of file G4CoupledTransportation.hh.

97{ return fMassGeometryLimitedStep; }

◆ PostStepDoIt()

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

Implements G4VProcess.

Definition at line 493 of file G4CoupledTransportation.cc.

495{
496 G4TouchableHandle retCurrentTouchable ; // The one to return
497
498 // Initialize ParticleChange (by setting all its members equal
499 // to corresponding members in G4Track)
500 // fParticleChange.Initialize(track) ; // To initialise TouchableChange
501
503
504 if( fSignifyStepInAnyVolume )
505 {
507 }
508 else
509 {
510 fParticleChange.ProposeFirstStepInVolume( fFirstStepInMassVolume );
511 }
512
513 // Check that the end position and direction are preserved
514 // since call to AlongStepDoIt
515
516#ifdef G4DEBUG_TRANSPORT
517 if( ( verboseLevel > 0 )
518 && ((fTransportEndPosition - track.GetPosition()).mag2() >= 1.0e-16) )
519 {
521 "End of Step Position" );
522 G4cerr << " Problem in G4CoupledTransportation::PostStepDoIt " << G4endl;
523 }
524
525 // If the Step was determined by the volume boundary, relocate the particle
526 // The pathFinder will know that the geometry limited the step (!?)
527
528 if( verboseLevel > 0 )
529 {
530 G4cout << " Calling PathFinder::Locate() from "
531 << " G4CoupledTransportation::PostStepDoIt " << G4endl;
532 G4cout << " fGeometryLimitedStep is " << fGeometryLimitedStep << G4endl;
533 }
534#endif
535
537 {
538 fPathFinder->Locate( track.GetPosition(),
539 track.GetMomentumDirection(),
540 true);
541
542 // fCurrentTouchable will now become the previous touchable,
543 // and what was the previous will be freed.
544 // (Needed because the preStepPoint can point to the previous touchable)
545
548
549#ifdef G4DEBUG_TRANSPORT
550 if( verboseLevel > 1 )
551 {
553 G4cout << "CHECK !!!!!!!!!!! fCurrentTouchableHandle->GetVolume() = "
554 << vol;
555 if( vol ) { G4cout << "Name=" << vol->GetName(); }
556 G4cout << G4endl;
557 }
558#endif
559
560 // Check whether the particle is out of the world volume
561 // If so it has exited and must be killed.
562 //
564 {
566 }
567 retCurrentTouchable = fCurrentTouchableHandle ;
568 // fParticleChange.SetTouchableHandle( fCurrentTouchableHandle ) ;
569 }
570 else // fGeometryLimitedStep is false
571 {
572#ifdef G4DEBUG_TRANSPORT
573 if( verboseLevel > 1 )
574 {
575 G4cout << "G4CoupledTransportation::PostStepDoIt -- "
576 << " fGeometryLimitedStep = " << fGeometryLimitedStep
577 << " must be false " << G4endl;
578 }
579#endif
580 // This serves only to move each of the Navigator's location
581 //
582 // fLinearNavigator->LocateGlobalPointWithinVolume( track.GetPosition() ) ;
583
584 fPathFinder->ReLocate( track.GetPosition() );
585 // track.GetMomentumDirection() );
586
587 // Keep the value of the track's current Touchable is retained,
588 // and use it to overwrite the (unset) one in particle change.
589 // Expect this must be fCurrentTouchable too
590 // - could it be different, eg at the start of a step ?
591 //
592 retCurrentTouchable = track.GetTouchableHandle() ;
593 // fParticleChange.SetTouchableHandle( track.GetTouchableHandle() ) ;
594 } // endif ( fGeometryLimitedStep )
595
596#ifdef G4DEBUG_NAVIGATION
597 G4cout << " CoupledTransport::AlongStep GPIL: "
598 << " last-step: any= " << fGeometryLimitedStep << " . ..... x . "
599 << " mass= " << fMassGeometryLimitedStep << G4endl;
600#endif
601
602 if( fSignifyStepInAnyVolume )
604 else
605 fParticleChange.ProposeLastStepInVolume(fMassGeometryLimitedStep);
606
607 SetTouchableInformation(retCurrentTouchable);
608
609 return &fParticleChange ;
610}
@ fStopAndKill
void ReportMove(G4ThreeVector OldVector, G4ThreeVector NewVector, const G4String &Quantity)
void ReLocate(const G4ThreeVector &position)
void Locate(const G4ThreeVector &position, const G4ThreeVector &direction, G4bool relativeSearch=true)
G4TouchableHandle CreateTouchableHandle(G4int navId) const
G4TrackStatus GetTrackStatus() const
const G4TouchableHandle & GetTouchableHandle() const
void SetTouchableInformation(const G4TouchableHandle &touchable)
G4TouchableHandle fCurrentTouchableHandle
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLastStepInVolume(G4bool flag)
void ProposeFirstStepInVolume(G4bool flag)
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:34

◆ ReportInexactEnergy()

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

Definition at line 649 of file G4CoupledTransportation.cc.

651{
652 static G4ThreadLocal G4int no_warnings= 0, warnModulo=1,
653 moduloFactor= 10, no_large_ediff= 0;
654
655 if( std::fabs(startEnergy- endEnergy) > perThousand * endEnergy )
656 {
657 no_large_ediff ++;
658 if( (no_large_ediff% warnModulo) == 0 )
659 {
660 no_warnings++;
661 std::ostringstream message;
662 message << "Energy change in Step is above 1^-3 relative value. "
663 << G4endl
664 << " Relative change in 'tracking' step = "
665 << std::setw(15) << (endEnergy-startEnergy)/startEnergy
666 << G4endl
667 << " Starting E= " << std::setw(12) << startEnergy / MeV
668 << " MeV " << G4endl
669 << " Ending E= " << std::setw(12) << endEnergy / MeV
670 << " MeV " << G4endl
671 << "Energy has been corrected -- however, review"
672 << " field propagation parameters for accuracy." << G4endl;
673 if ( (verboseLevel > 2 ) || (no_warnings<4)
674 || (no_large_ediff == warnModulo * moduloFactor) )
675 {
676 message << "These include EpsilonStepMax(/Min) in G4FieldManager,"
677 << G4endl
678 << "which determine fractional error per step for integrated quantities."
679 << G4endl
680 << "Note also the influence of the permitted number of integration steps."
681 << G4endl;
682 }
683 message << "Bad 'endpoint'. Energy change detected and corrected."
684 << G4endl
685 << "Has occurred already " << no_large_ediff << " times.";
686 G4Exception("G4CoupledTransportation::AlongStepGetPIL()",
687 "EnergyChange", JustWarning, message);
688 if( no_large_ediff == warnModulo * moduloFactor )
689 {
690 warnModulo *= moduloFactor;
691 }
692 }
693 }
694}
@ JustWarning
#define G4ThreadLocal
Definition: tls.hh:77

Referenced by AlongStepGetPhysicalInteractionLength().

◆ ReportMove()

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

Definition at line 473 of file G4CoupledTransportation.cc.

476{
477 G4ThreeVector moveVec = ( NewVector - OldVector );
478
479 G4cerr << G4endl
480 << "**************************************************************"
481 << G4endl;
482 G4cerr << "Endpoint has moved between value expected from TransportEndPosition "
483 << " and value from Track in PostStepDoIt. " << G4endl
484 << "Change of " << Quantity << " is " << moveVec.mag() / mm
485 << " mm long, "
486 << " and its vector is " << (1.0/mm) * moveVec << " mm " << G4endl
487 << "Endpoint of ComputeStep was " << OldVector
488 << " and current position to locate is " << NewVector << G4endl;
489}
double mag() const

Referenced by PostStepDoIt().

◆ SetSignifyStepsInAnyVolume()

static void G4CoupledTransportation::SetSignifyStepsInAnyVolume ( G4bool  anyVol)
inlinestatic

Definition at line 83 of file G4CoupledTransportation.hh.

84 { fSignifyStepInAnyVolume = anyVol; }

◆ StartTracking()

void G4CoupledTransportation::StartTracking ( G4Track aTrack)
virtual

Reimplemented from G4VProcess.

Definition at line 619 of file G4CoupledTransportation.cc.

620{
622
624 G4ThreeVector direction = aTrack->GetMomentumDirection();
625
626 fPathFinder->PrepareNewTrack( position, direction);
627 // This implies a call to fPathFinder->Locate( position, direction );
628
629 // reset safety value and center
630 //
631 fPreviousMassSafety = 0.0 ;
632 fPreviousFullSafety = 0.0 ;
634}
void PrepareNewTrack(const G4ThreeVector &position, const G4ThreeVector &direction, G4VPhysicalVolume *massStartVol=nullptr)
void StartTracking(G4Track *aTrack)

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