ATLAS Offline Software
Loading...
Searching...
No Matches
QuirkTransportation Class Reference

#include <QuirkTransportation.h>

Inheritance diagram for QuirkTransportation:
Collaboration diagram for QuirkTransportation:

Public Member Functions

 QuirkTransportation (G4int verbosityLevel=1)
 ~QuirkTransportation ()
G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
G4VParticleChange * AlongStepDoIt (const G4Track &track, const G4Step &stepData)
G4VParticleChange * PostStepDoIt (const G4Track &track, const G4Step &stepData)
G4double PostStepGetPhysicalInteractionLength (const G4Track &, G4double previousStepSize, G4ForceCondition *pForceCond)
G4PropagatorInField * GetPropagatorInField ()
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)
G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
G4VParticleChange * AtRestDoIt (const G4Track &, const G4Step &)
void StartTracking (G4Track *aTrack)

Protected Member Functions

G4bool DoesGlobalFieldExist ()

Private Attributes

G4Navigator * m_linearNavigator
G4PropagatorInField * m_fieldPropagator
G4ThreeVector m_transportEndPosition
G4ThreeVector m_transportEndMomentumDir
G4double m_transportEndKineticEnergy
G4ThreeVector m_transportEndSpin
G4bool m_momentumChanged
G4double m_candidateEndGlobalTime
G4bool m_particleIsLooping
G4TouchableHandle m_currentTouchableHandle
G4bool m_geometryLimitedStep
G4ThreeVector m_previousSftOrigin
G4double m_previousSafety
G4ParticleChangeForTransport m_particleChange
G4double m_endpointDistance
G4double m_threshold_Warning_Energy
G4double m_threshold_Important_Energy
G4int m_thresholdTrials
G4int m_noLooperTrials
G4double m_sumEnergyKilled
G4double m_maxEnergyKilled
G4bool m_shortStepOptimisation
G4SafetyHelper * m_safetyHelper
G4int m_verboseLevel

Detailed Description

Definition at line 59 of file QuirkTransportation.h.

Constructor & Destructor Documentation

◆ QuirkTransportation()

QuirkTransportation::QuirkTransportation ( G4int verbosityLevel = 1)

Definition at line 79 of file QuirkTransportation.cxx.

80 : G4VProcess( G4String("QuirkTransportation"), fTransportation ),
82 m_momentumChanged( false ),
83 m_particleIsLooping( false ),
84 m_currentTouchableHandle(), // Points to (G4VTouchable*) 0
85 m_geometryLimitedStep( false ),
86 m_previousSftOrigin (0.,0.,0.),
87 m_previousSafety ( 0.0 ),
88 m_endpointDistance( 0.0 ),
89 m_threshold_Warning_Energy( 100 * CLHEP::MeV ),
90 m_threshold_Important_Energy( 250 * CLHEP::MeV ),
91 m_thresholdTrials( 10 ),
94 m_shortStepOptimisation(false), // Old default: true (=fast short steps)
95 m_verboseLevel( verboseLevel )
96{
97 G4TransportationManager* transportMgr ;
98
99 transportMgr = G4TransportationManager::GetTransportationManager() ;
100
101 m_linearNavigator = transportMgr->GetNavigatorForTracking() ;
102
103 // fGlobalFieldMgr = transportMgr->GetFieldManager() ;
104
105 m_fieldPropagator = transportMgr->GetPropagatorInField() ;
106
107 m_safetyHelper = transportMgr->GetSafetyHelper(); // New
108
109 // Cannot determine whether a field exists here, as it would
110 // depend on the relative order of creating the detector's
111 // field and this process. That order is not guaranted.
112 // Instead later the method DoesGlobalFieldExist() is called
113
115}
G4ThreeVector m_previousSftOrigin
G4PropagatorInField * m_fieldPropagator
G4Navigator * m_linearNavigator
G4SafetyHelper * m_safetyHelper
G4TouchableHandle m_currentTouchableHandle

◆ ~QuirkTransportation()

QuirkTransportation::~QuirkTransportation ( )

Definition at line 119 of file QuirkTransportation.cxx.

120{
121 if( (m_verboseLevel > 0) && (m_sumEnergyKilled > 0.0 ) ){
122 G4cout << " QuirkTransportation: Statistics for looping particles " << G4endl;
123 G4cout << " Sum of energy of loopers killed: " << m_sumEnergyKilled << G4endl;
124 G4cout << " Max energy of loopers killed: " << m_maxEnergyKilled << G4endl;
125 }
126}

Member Function Documentation

◆ AlongStepDoIt()

G4VParticleChange * QuirkTransportation::AlongStepDoIt ( const G4Track & track,
const G4Step & stepData )

Definition at line 353 of file QuirkTransportation.cxx.

355{
356#ifdef G4VERBOSE
357 static std::atomic<G4int> noCalls=0;
358 noCalls++;
359#endif
360
361 m_particleChange.Initialize(track) ;
362
363 // Code for specific process
364 //
366 m_particleChange.ProposeMomentumDirection(m_transportEndMomentumDir) ;
368 m_particleChange.SetMomentumChanged(m_momentumChanged) ;
369
370 m_particleChange.ProposePolarization(m_transportEndSpin);
371
372 // Calculate Lab Time of Flight (ONLY if field Equations used it!)
373 // G4double endTime = m_candidateEndGlobalTime;
374 // G4double delta_time = endTime - startTime;
375
376 G4double startTime = track.GetGlobalTime() ;
377 G4double deltaTime = m_candidateEndGlobalTime - startTime ;
378 m_particleChange.ProposeGlobalTime( m_candidateEndGlobalTime ) ;
379
380 // Now Correct by Lorentz factor to get "proper" deltaTime
381
382 G4double restMass = track.GetDynamicParticle()->GetMass() ;
383 G4double deltaProperTime = deltaTime*( restMass/track.GetTotalEnergy() ) ;
384
385 m_particleChange.ProposeProperTime(track.GetProperTime() + deltaProperTime) ;
386 //m_particleChange. ProposeTrueStepLength( track.GetStepLength() ) ;
387
388 // If the particle is caught looping or is stuck (in very difficult
389 // boundaries) in a magnetic field (doing many steps)
390 // THEN this kills it ...
391 //
393 {
394 G4double endEnergy= m_transportEndKineticEnergy;
395
396 if( (endEnergy < m_threshold_Important_Energy)
398 // Kill the looping particle
399 //
400 m_particleChange.ProposeTrackStatus( fStopAndKill ) ;
401
402 // 'Bare' statistics
403 m_sumEnergyKilled += endEnergy;
404 if( endEnergy > m_maxEnergyKilled) { m_maxEnergyKilled= endEnergy; }
405
406#ifdef G4VERBOSE
407 if( (m_verboseLevel > 1) ||
408 ( endEnergy > m_threshold_Warning_Energy ) ) {
409 G4cout << " QuirkTransportation is killing track that is looping or stuck "
410 << G4endl
411 << " This track has " << track.GetKineticEnergy() / CLHEP::MeV
412 << " MeV energy." << G4endl;
413 G4cout << " Number of trials = " << m_noLooperTrials
414 << " No of calls to AlongStepDoIt = " << noCalls
415 << G4endl;
416 }
417#endif
419 }
420 else{
422#ifdef G4VERBOSE
423 if( (m_verboseLevel > 2) ){
424 G4cout << " QuirkTransportation::AlongStepDoIt(): Particle looping - "
425 << " Number of trials = " << m_noLooperTrials
426 << " No of calls to = " << noCalls
427 << G4endl;
428 }
429#endif
430 }
431 }else{
433 }
434
435 // Another (sometimes better way) is to use a user-limit maximum Step size
436 // to alleviate this problem ..
437
438 // Introduce smooth curved trajectories to particle-change
439 //
440 m_particleChange.SetPointerToVectorOfAuxiliaryPoints
441 (m_fieldPropagator->GimmeTrajectoryVectorAndForgetIt() );
442
443 return &m_particleChange ;
444}
G4ThreeVector m_transportEndPosition
G4ThreeVector m_transportEndSpin
G4ParticleChangeForTransport m_particleChange
G4ThreeVector m_transportEndMomentumDir

◆ AlongStepGetPhysicalInteractionLength()

G4double QuirkTransportation::AlongStepGetPhysicalInteractionLength ( const G4Track & track,
G4double previousStepSize,
G4double currentMinimumStep,
G4double & currentSafety,
G4GPILSelection * selection )

Definition at line 135 of file QuirkTransportation.cxx.

141{
142 G4double geometryStepLength;
143 m_particleIsLooping = false ;
144
145 // Initial actions moved to StartTrack()
146 // --------------------------------------
147 // Note: in case another process changes touchable handle
148 // it will be necessary to add here (for all steps)
149 // m_currentTouchableHandle = aTrack->GetTouchableHandle();
150
151 // GPILSelection is set to defaule value of CandidateForSelection
152 // It is a return value
153 //
154 *selection = CandidateForSelection ;
155
156 // Get initial Energy/Momentum of the track
157 //
158 const G4DynamicParticle* pParticle = track.GetDynamicParticle() ;
159 G4ParticleDefinition* pParticleDef = pParticle->GetDefinition() ;
160 G4ThreeVector startPosition = track.GetPosition() ;
161
162 // G4double theTime = track.GetGlobalTime() ;
163
164 // The Step Point safety can be limited by other geometries and/or the
165 // assumptions of any process - it's not always the geometrical safety.
166 // We calculate the starting point's isotropic safety here.
167 //
168 G4ThreeVector OriginShift = startPosition - m_previousSftOrigin ;
169 G4double MagSqShift = OriginShift.mag2() ;
170 if( MagSqShift >= sqr(m_previousSafety) )
171 {
172 currentSafety = 0.0 ;
173 }
174 else
175 {
176 currentSafety = m_previousSafety - std::sqrt(MagSqShift) ;
177 }
178
179 // Is the particle charged ?
180 //
181 G4double particleCharge = pParticle->GetCharge() ;
182
183 // There is no need to locate the current volume. It is Done elsewhere:
184 // On track construction
185 // By the tracking, after all AlongStepDoIts, in "Relocation"
186
187 // Set up field manager
188 G4FieldManager* fieldMgr = m_fieldPropagator->FindAndSetFieldManager( track.GetVolume() );
189 if (fieldMgr != 0) {
190 // Message the field Manager, to configure it for this track
191 fieldMgr->ConfigureForTrack( &track );
192 } else {
193 G4Exception("QuirkTransportation::AlongStepGetPhysicalInteractionLength", "QuirkNoFieldMgr", RunMustBeAborted, "no field manager");
194 }
195
196 // Set up stepper to calculate quirk trajectory
197 Quirk* quirkDef = dynamic_cast<Quirk*>(pParticleDef);
198 if (quirkDef == 0) {
199 G4Exception("QuirkTransportation::AlongStepGetPhysicalInteractionLength", "NonQuirk", FatalErrorInArgument, "QuirkTransportation run on non-quirk particle");
200 std::abort();
201 }
202 HyperbolaStepper quirkStepper(
203 quirkDef->GetStringIn(),
204 track,
205 fieldMgr->DoesFieldExist() ? fieldMgr->GetDetectorField() : 0
206 );
207
208 // Switch out chord finder with quirk chord finder
209 G4ChordFinder* oldChordFinder = fieldMgr->GetChordFinder();
210 G4ChordFinder quirkChordFinder(
211 new G4MagInt_Driver(0.0, &quirkStepper, quirkStepper.GetNumberOfVariables())
212 );
213 fieldMgr->SetChordFinder(&quirkChordFinder);
214
215 G4double momentumMagnitude = pParticle->GetTotalMomentum() ;
216 G4ThreeVector EndUnitMomentum ;
217 G4double restMass = pParticleDef->GetPDGMass() ;
218
219#if G4VERSION_NUMBER > 1009
220 //FIXME untested!!!
221 G4ChargeState chargeState(particleCharge, // in e+ units
222 quirkDef->GetPDGSpin(),
223 0,
224 0,
225 0); //no magnetic charge?
226 G4EquationOfMotion* equationOfMotion = (m_fieldPropagator->GetChordFinder()->GetIntegrationDriver()->GetStepper())->GetEquationOfMotion();
227 equationOfMotion->SetChargeMomentumMass( chargeState,
228 momentumMagnitude, // in Mev/c
229 restMass ) ;
230#else
231 m_fieldPropagator->SetChargeMomentumMass( particleCharge, // in e+ units
232 momentumMagnitude, // in Mev/c
233 restMass ) ;
234#endif
235
236 G4ThreeVector spin = track.GetPolarization() ;
237 G4FieldTrack aFieldTrack = G4FieldTrack( startPosition,
238 track.GetMomentumDirection(),
239 0.0,
240 track.GetKineticEnergy(),
241 restMass,
242 track.GetVelocity(),
243 track.GetGlobalTime(), // Lab.
244 0.0, // substituting step length for proper time
245 &spin ) ;
246
247 //if (currentMinimumStep <= 0) {
248 // G4cout << "QuirkTransportation: currentMinimumStep = " << currentMinimumStep << G4endl;
249 // G4Exception(
250 // "QuirkTransportation::AlongStepGetPhysicalInteractionLength",
251 // "BadMinimumStep",
252 // EventMustBeAborted,
253 // "QuirkTransportation: currentMinimumStep <= 0"
254 // );
255 //}
256 G4bool dbg = false; //(track.GetCurrentStepNumber() % 1000000 == 0);
257 quirkStepper.SetDebug(dbg);
258 if (dbg) G4cout << "QuirkTransportation: start = " << aFieldTrack.GetPosition() << G4endl;
259 if (dbg) G4cout << "QuirkTransportation: currentMinimumStep = " << currentMinimumStep << G4endl;
260 if (dbg) G4cout << "QuirkTransportation: maxlength = " << quirkStepper.GetMaxLength() << G4endl;
261 currentMinimumStep = std::min(currentMinimumStep, quirkStepper.GetMaxLength());
262 if (dbg) G4cout << "QuirkTransportation: currentMinimumStep = " << currentMinimumStep << G4endl;
263 if( currentMinimumStep > 0 )
264 {
265 // Call field propagator to handle boundary crossings
266 G4double lengthAlongCurve = m_fieldPropagator->ComputeStep( aFieldTrack,
267 currentMinimumStep,
268 currentSafety,
269 track.GetVolume() ) ;
270 if (dbg) G4cout << "QuirkTransportation: moved " << lengthAlongCurve << G4endl;
271 m_geometryLimitedStep = lengthAlongCurve < currentMinimumStep;
273 geometryStepLength = lengthAlongCurve ;
274 } else {
275 geometryStepLength = currentMinimumStep ;
276 }
277 // Update stepper, string vectors with step length from field propagator
278 quirkStepper.Update(aFieldTrack, aFieldTrack.GetCurveLength() == 0 && !m_geometryLimitedStep);
279 }
280 else
281 {
282 geometryStepLength = 0.0 ;
283 m_geometryLimitedStep = false ;
284 }
285 if (dbg) G4cout << "QuirkTransportation: moved " << aFieldTrack.GetCurveLength() << G4endl;
286 if (dbg) G4cout << "QuirkTransportation: end = " << aFieldTrack.GetPosition() << " [" << aFieldTrack.GetProperTimeOfFlight() << "]" << G4endl;
287
288 // Remember last safety origin & value.
289 //
290 m_previousSftOrigin = startPosition ;
291 m_previousSafety = currentSafety ;
292 // m_safetyHelper->SetCurrentSafety( newSafety, startPosition);
293
294 // Get end-of-step quantities
295 geometryStepLength = aFieldTrack.GetCurveLength();
296 m_transportEndPosition = aFieldTrack.GetPosition() ;
297 m_momentumChanged = true ;
298 m_transportEndMomentumDir = aFieldTrack.GetMomentumDir() ;
299 m_transportEndKineticEnergy = aFieldTrack.GetKineticEnergy() ;
300 m_candidateEndGlobalTime = aFieldTrack.GetLabTimeOfFlight();
301 m_transportEndSpin = aFieldTrack.GetSpin();
302 m_particleIsLooping = m_fieldPropagator->IsParticleLooping() ;
303 m_endpointDistance = (m_transportEndPosition - startPosition).mag() ;
304
305 // If we are asked to go a step length of 0, and we are on a boundary
306 // then a boundary will also limit the step -> we must flag this.
307 //
308 if( currentMinimumStep == 0.0 )
309 {
310 if( currentSafety == 0.0 ) m_geometryLimitedStep = true ;
311 }
312
313 // Update the safety starting from the end-point,
314 // if it will become negative at the end-point.
315 //
316 if( currentSafety < m_endpointDistance )
317 {
318 G4double endSafety =
320 currentSafety = endSafety ;
322 m_previousSafety = currentSafety ;
323 m_safetyHelper->SetCurrentSafety( currentSafety, m_transportEndPosition);
324
325 // Because the Stepping Manager assumes it is from the start point,
326 // add the StepLength
327 //
328 currentSafety += m_endpointDistance ;
329
330#ifdef G4DEBUG_TRANSPORT
331 G4cout.precision(12) ;
332 G4cout << "***QuirkTransportation::AlongStepGPIL ** " << G4endl ;
333 G4cout << " Called Navigator->ComputeSafety at " << m_transportEndPosition
334 << " and it returned safety= " << endSafety << G4endl ;
335 G4cout << " Adding endpoint distance " << m_endpointDistance
336 << " to obtain pseudo-safety= " << currentSafety << G4endl ;
337#endif
338 }
339
340 m_particleChange.ProposeTrueStepLength(geometryStepLength) ;
341
342 // Restore original stepper
343 fieldMgr->SetChordFinder(oldChordFinder);
344
345 return geometryStepLength ;
346}
Scalar mag() const
mag method
double spin(const T &p)
Definition AtlasPID.h:1147
#define sqr(t)
const InfracolorForce & GetStringIn() const
Definition Quirk.h:29
const std::string selection
bool dbg
Definition listroot.cxx:36

◆ AtRestDoIt()

G4VParticleChange * QuirkTransportation::AtRestDoIt ( const G4Track & ,
const G4Step &  )
inline

Definition at line 133 of file QuirkTransportation.h.

136 {return 0;};

◆ AtRestGetPhysicalInteractionLength()

G4double QuirkTransportation::AtRestGetPhysicalInteractionLength ( const G4Track & ,
G4ForceCondition *  )
inline

Definition at line 127 of file QuirkTransportation.h.

130 { return -1.0; };

◆ DoesGlobalFieldExist()

G4bool QuirkTransportation::DoesGlobalFieldExist ( )
protected

◆ EnableShortStepOptimisation()

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

◆ GetMaxEnergyKilled()

G4double QuirkTransportation::GetMaxEnergyKilled ( ) const
inline

◆ GetPropagatorInField()

G4PropagatorInField * QuirkTransportation::GetPropagatorInField ( )

◆ GetSumEnergyKilled()

G4double QuirkTransportation::GetSumEnergyKilled ( ) const
inline

◆ GetThresholdImportantEnergy()

G4double QuirkTransportation::GetThresholdImportantEnergy ( ) const
inline

◆ GetThresholdTrials()

G4int QuirkTransportation::GetThresholdTrials ( ) const
inline

◆ GetThresholdWarningEnergy()

G4double QuirkTransportation::GetThresholdWarningEnergy ( ) const
inline

◆ GetVerboseLevel()

G4int QuirkTransportation::GetVerboseLevel ( ) const
inline

◆ PostStepDoIt()

G4VParticleChange * QuirkTransportation::PostStepDoIt ( const G4Track & track,
const G4Step & stepData )

Definition at line 464 of file QuirkTransportation.cxx.

466{
467 G4TouchableHandle retCurrentTouchable ; // The one to return
468 G4bool isLastStep= false;
469
470 // Initialize ParticleChange (by setting all its members equal
471 // to corresponding members in G4Track)
472 // m_particleChange.Initialize(track) ; // To initialise TouchableChange
473
474 m_particleChange.ProposeTrackStatus(track.GetTrackStatus()) ;
475
476 // If the Step was determined by the volume boundary,
477 // logically relocate the particle
478
480 {
481 // fCurrentTouchable will now become the previous touchable,
482 // and what was the previous will be freed.
483 // (Needed because the preStepPoint can point to the previous touchable)
484
485 m_linearNavigator->SetGeometricallyLimitedStep() ;
487 LocateGlobalPointAndUpdateTouchableHandle( track.GetPosition(),
488 track.GetMomentumDirection(),
490 true ) ;
491 // Check whether the particle is out of the world volume
492 // If so it has exited and must be killed.
493 //
494 if( m_currentTouchableHandle->GetVolume() == 0 )
495 {
496 m_particleChange.ProposeTrackStatus( fStopAndKill ) ;
497 G4cout << "QuirkTransportation: number of steps = " << track.GetCurrentStepNumber() << G4endl;
498 }
499 retCurrentTouchable = m_currentTouchableHandle ;
500 m_particleChange.SetTouchableHandle( m_currentTouchableHandle ) ;
501
502 // Update the Step flag which identifies the Last Step in a volume
503 isLastStep = m_linearNavigator->ExitedMotherVolume()
504 || m_linearNavigator->EnteredDaughterVolume() ;
505
506#ifdef G4DEBUG_TRANSPORT
507 // Checking first implementation of flagging Last Step in Volume
508 G4bool exiting = m_linearNavigator->ExitedMotherVolume();
509 G4bool entering = m_linearNavigator->EnteredDaughterVolume();
510 if( ! (exiting || entering) ) {
511 G4cout << " Transport> : Proposed isLastStep= " << isLastStep
512 << " Exiting " << m_linearNavigator->ExitedMotherVolume()
513 << " Entering " << m_linearNavigator->EnteredDaughterVolume()
514 << G4endl;
515 }
516#endif
517 }
518 else // m_geometryLimitedStep is false
519 {
520 // This serves only to move the Navigator's location
521 //
522 m_linearNavigator->LocateGlobalPointWithinVolume( track.GetPosition() ) ;
523
524 // The value of the track's current Touchable is retained.
525 // (and it must be correct because we must use it below to
526 // overwrite the (unset) one in particle change)
527 // It must be fCurrentTouchable too ??
528 //
529 m_particleChange.SetTouchableHandle( track.GetTouchableHandle() ) ;
530 retCurrentTouchable = track.GetTouchableHandle() ;
531
532 isLastStep= false;
533#ifdef G4DEBUG_TRANSPORT
534 // Checking first implementation of flagging Last Step in Volume
535 G4cout << " Transport> Proposed isLastStep= " << isLastStep
536 << " Geometry did not limit step. " << G4endl;
537#endif
538 } // endif ( m_geometryLimitedStep )
539
540 m_particleChange.ProposeLastStepInVolume(isLastStep);
541
542 const G4VPhysicalVolume* pNewVol = retCurrentTouchable->GetVolume() ;
543 G4Material* pNewMaterial = 0 ;
544 G4VSensitiveDetector* pNewSensitiveDetector = 0 ;
545
546 if( pNewVol != 0 )
547 {
548 pNewMaterial= pNewVol->GetLogicalVolume()->GetMaterial();
549 pNewSensitiveDetector= pNewVol->GetLogicalVolume()->GetSensitiveDetector();
550 }
551
552 // ( <const_cast> pNewMaterial ) ;
553 // ( <const_cast> pNewSensitiveDetector) ;
554
555 m_particleChange.SetMaterialInTouchable( pNewMaterial ) ;
556 m_particleChange.SetSensitiveDetectorInTouchable( pNewSensitiveDetector ) ;
557
558 const G4MaterialCutsCouple* pNewMaterialCutsCouple = 0;
559 if( pNewVol != 0 )
560 {
561 pNewMaterialCutsCouple=pNewVol->GetLogicalVolume()->GetMaterialCutsCouple();
562 }
563
564 if( pNewVol!=0 && pNewMaterialCutsCouple!=0 && pNewMaterialCutsCouple->GetMaterial()!=pNewMaterial )
565 {
566 // for parametrized volume
567 //
568 pNewMaterialCutsCouple =
569 G4ProductionCutsTable::GetProductionCutsTable()
570 ->GetMaterialCutsCouple(pNewMaterial,
571 pNewMaterialCutsCouple->GetProductionCuts());
572 }
573 m_particleChange.SetMaterialCutsCoupleInTouchable( pNewMaterialCutsCouple );
574
575 // temporarily until Get/Set Material of ParticleChange,
576 // and StepPoint can be made const.
577 // Set the touchable in ParticleChange
578 // this must always be done because the particle change always
579 // uses this value to overwrite the current touchable pointer.
580 //
581 m_particleChange.SetTouchableHandle(retCurrentTouchable) ;
582
583 return &m_particleChange ;
584}

◆ PostStepGetPhysicalInteractionLength()

G4double QuirkTransportation::PostStepGetPhysicalInteractionLength ( const G4Track & ,
G4double previousStepSize,
G4ForceCondition * pForceCond )

Definition at line 452 of file QuirkTransportation.cxx.

456{
457 *pForceCond = Forced ;
458 return DBL_MAX ; // was kInfinity ; but convention now is DBL_MAX
459}

◆ ResetKilledStatistics()

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

◆ SetPropagatorInField()

void QuirkTransportation::SetPropagatorInField ( G4PropagatorInField * pFieldPropagator)

◆ SetThresholdImportantEnergy()

void QuirkTransportation::SetThresholdImportantEnergy ( G4double newEnImp)
inline

◆ SetThresholdTrials()

void QuirkTransportation::SetThresholdTrials ( G4int newMaxTrials)
inline

◆ SetThresholdWarningEnergy()

void QuirkTransportation::SetThresholdWarningEnergy ( G4double newEnWarn)
inline

◆ SetVerboseLevel()

void QuirkTransportation::SetVerboseLevel ( G4int verboseLevel)
inline

◆ StartTracking()

void QuirkTransportation::StartTracking ( G4Track * aTrack)

Definition at line 590 of file QuirkTransportation.cxx.

591{
592 G4VProcess::StartTracking(aTrack);
593
594// The actions here are those that were taken in AlongStepGPIL
595// when track.GetCurrentStepNumber()==1
596
597 // reset safety value and center
598 //
599 m_previousSafety = 0.0 ;
600 m_previousSftOrigin = G4ThreeVector(0.,0.,0.) ;
601
602 // reset looping counter -- for motion in field
604 // Must clear this state .. else it depends on last track's value
605 // --> a better solution would set this from state of suspended track TODO ?
606 // Was if( aTrack->GetCurrentStepNumber()==1 ) { .. }
607
608 // ChordFinder reset internal state
609 //
610 m_fieldPropagator->ClearPropagatorState();
611 // Resets all state of field propagator class (ONLY)
612 // including safety values (in case of overlaps and to wipe for first track).
613
614 // G4ChordFinder* chordF= m_fieldPropagator->GetChordFinder();
615 // if( chordF ) chordF->ResetStepEstimate();
616
617 // Make sure to clear the chord finders of all fields (ie managers)
618 G4FieldManagerStore::GetInstance()->ClearAllChordFindersState();
619
620 // Update the current touchable handle (from the track's)
621 //
622 m_currentTouchableHandle = aTrack->GetTouchableHandle();
623
624 // Initialize infracolor string
625 auto part_nc ATLAS_THREAD_SAFE = // aTrack is non-const but GetParticleDefinition only returns const
626 const_cast<G4ParticleDefinition*>(aTrack->GetParticleDefinition());
627 Quirk* quirkDef = dynamic_cast<Quirk*>(part_nc);
628 if (quirkDef == 0) {
629 G4Exception("QuirkTransportation::StartTracking", "NonQuirk", FatalErrorInArgument, "QuirkTransportation run on non-quirk particle");
630 std::abort();
631 }
632 quirkDef->GetStringIn().StartTracking(aTrack);
633}
#define ATLAS_THREAD_SAFE
void StartTracking(const G4Track *dest)

Member Data Documentation

◆ m_candidateEndGlobalTime

G4double QuirkTransportation::m_candidateEndGlobalTime
private

Definition at line 161 of file QuirkTransportation.h.

◆ m_currentTouchableHandle

G4TouchableHandle QuirkTransportation::m_currentTouchableHandle
private

Definition at line 166 of file QuirkTransportation.h.

◆ m_endpointDistance

G4double QuirkTransportation::m_endpointDistance
private

Definition at line 183 of file QuirkTransportation.h.

◆ m_fieldPropagator

G4PropagatorInField* QuirkTransportation::m_fieldPropagator
private

Definition at line 150 of file QuirkTransportation.h.

◆ m_geometryLimitedStep

G4bool QuirkTransportation::m_geometryLimitedStep
private

Definition at line 173 of file QuirkTransportation.h.

◆ m_linearNavigator

G4Navigator* QuirkTransportation::m_linearNavigator
private

Definition at line 149 of file QuirkTransportation.h.

◆ m_maxEnergyKilled

G4double QuirkTransportation::m_maxEnergyKilled
private

Definition at line 196 of file QuirkTransportation.h.

◆ m_momentumChanged

G4bool QuirkTransportation::m_momentumChanged
private

Definition at line 160 of file QuirkTransportation.h.

◆ m_noLooperTrials

G4int QuirkTransportation::m_noLooperTrials
private

Definition at line 193 of file QuirkTransportation.h.

◆ m_particleChange

G4ParticleChangeForTransport QuirkTransportation::m_particleChange
private

Definition at line 180 of file QuirkTransportation.h.

◆ m_particleIsLooping

G4bool QuirkTransportation::m_particleIsLooping
private

Definition at line 164 of file QuirkTransportation.h.

◆ m_previousSafety

G4double QuirkTransportation::m_previousSafety
private

Definition at line 177 of file QuirkTransportation.h.

◆ m_previousSftOrigin

G4ThreeVector QuirkTransportation::m_previousSftOrigin
private

Definition at line 176 of file QuirkTransportation.h.

◆ m_safetyHelper

G4SafetyHelper* QuirkTransportation::m_safetyHelper
private

Definition at line 202 of file QuirkTransportation.h.

◆ m_shortStepOptimisation

G4bool QuirkTransportation::m_shortStepOptimisation
private

Definition at line 200 of file QuirkTransportation.h.

◆ m_sumEnergyKilled

G4double QuirkTransportation::m_sumEnergyKilled
private

Definition at line 195 of file QuirkTransportation.h.

◆ m_threshold_Important_Energy

G4double QuirkTransportation::m_threshold_Important_Energy
private

Definition at line 188 of file QuirkTransportation.h.

◆ m_threshold_Warning_Energy

G4double QuirkTransportation::m_threshold_Warning_Energy
private

Definition at line 187 of file QuirkTransportation.h.

◆ m_thresholdTrials

G4int QuirkTransportation::m_thresholdTrials
private

Definition at line 189 of file QuirkTransportation.h.

◆ m_transportEndKineticEnergy

G4double QuirkTransportation::m_transportEndKineticEnergy
private

Definition at line 158 of file QuirkTransportation.h.

◆ m_transportEndMomentumDir

G4ThreeVector QuirkTransportation::m_transportEndMomentumDir
private

Definition at line 157 of file QuirkTransportation.h.

◆ m_transportEndPosition

G4ThreeVector QuirkTransportation::m_transportEndPosition
private

Definition at line 156 of file QuirkTransportation.h.

◆ m_transportEndSpin

G4ThreeVector QuirkTransportation::m_transportEndSpin
private

Definition at line 159 of file QuirkTransportation.h.

◆ m_verboseLevel

G4int QuirkTransportation::m_verboseLevel
private

Definition at line 205 of file QuirkTransportation.h.


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