Belle II Software development
G4MonopoleTransportation.cc
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8
9// modified from GEANT4 exoticphysics/monopole/*
10
11#include <simulation/monopoles/G4MonopoleTransportation.h>
12#include <simulation/monopoles/G4Monopole.h>
13
14#include <G4ProductionCutsTable.hh>
15#include <G4ParticleTable.hh>
16#include <G4ChordFinder.hh>
17#include <G4SafetyHelper.hh>
18#include <G4FieldManagerStore.hh>
19#include <G4TransportationProcessType.hh>
20#include <CLHEP/Units/SystemOfUnits.h>
21#include <framework/logging/Logger.h>
22
23class G4VSensitiveDetector;
24
25using namespace std;
26using namespace Belle2;
27using namespace Belle2::Monopoles;
28using namespace CLHEP;
29
31 G4int verb)
32 : G4VProcess(G4String("MonopoleTransportation"), fTransportation),
33 fParticleDef(mpl),
34 fMagSetup(0),
35 fLinearNavigator(0),
36 fFieldPropagator(0),
37 fParticleIsLooping(false),
38 fPreviousSftOrigin(0., 0., 0.),
39 fPreviousSafety(0.0),
40 fThreshold_Trap_Energy(10 * MeV),
41 fThreshold_Warning_Energy(100 * MeV),
42 fThreshold_Important_Energy(250 * MeV),
43 fThresholdTrials(10),
44 fNoLooperTrials(0),
45 fSumEnergyKilled(0.0), fMaxEnergyKilled(0.0),
46 fShortStepOptimisation(false), // Old default: true (=fast short steps)
47 fpSafetyHelper(0)
48{
49 verboseLevel = verb;
50
51 // set Process Sub Type
52 SetProcessSubType(TRANSPORTATION);
53
55
56 G4TransportationManager* transportMgr ;
57
58 transportMgr = G4TransportationManager::GetTransportationManager() ;
59
60 fLinearNavigator = transportMgr->GetNavigatorForTracking() ;
61
62 // fGlobalFieldMgr = transportMgr->GetFieldManager() ;
63
64 fFieldPropagator = transportMgr->GetPropagatorInField() ;
65
66 fpSafetyHelper = transportMgr->GetSafetyHelper(); // New
67
68 // Cannot determine whether a field exists here,
69 // because it would only work if the field manager has informed
70 // about the detector's field before this transportation process
71 // is constructed.
72 // Instead later the method DoesGlobalFieldExist() is called
73
74 static G4TouchableHandle nullTouchableHandle; // Points to (G4VTouchable*) 0
75 fCurrentTouchableHandle = nullTouchableHandle;
76
79}
80
82{
83 if ((verboseLevel > 0) && (fSumEnergyKilled > 0.0)) {
84 B2DEBUG(25, " G4MonopoleTransportation: Statistics for looping particles ");
85 B2DEBUG(25, " Sum of energy of loopers killed: " << fSumEnergyKilled);
86 B2DEBUG(25, " Max energy of loopers killed: " << fMaxEnergyKilled);
87 }
88}
89
90// Responsibilities:
91// Find whether the geometry limits the Step, and to what length
92// Calculate the new value of the safety and return it.
93// Store the final time, position and momentum.
94
96AlongStepGetPhysicalInteractionLength(const G4Track& track,
97 G4double, // previousStepSize
98 G4double currentMinimumStep,
99 G4double& currentSafety,
100 G4GPILSelection* selection)
101{
103 // change to monopole equation
104
105 G4double geometryStepLength, newSafety ;
106 fParticleIsLooping = false ;
107
108 // Initial actions moved to StartTrack()
109 // --------------------------------------
110 // Note: in case another process changes touchable handle
111 // it will be necessary to add here (for all steps)
112 // fCurrentTouchableHandle = aTrack->GetTouchableHandle();
113
114 // GPILSelection is set to defaule value of CandidateForSelection
115 // It is a return value
116 //
117 *selection = CandidateForSelection ;
118
119 // Get initial Energy/Momentum of the track
120 //
121 const G4DynamicParticle* pParticle = track.GetDynamicParticle() ;
122 G4ThreeVector startMomentumDir = pParticle->GetMomentumDirection() ;
123 G4ThreeVector startPosition = track.GetPosition() ;
124
125 // G4double theTime = track.GetGlobalTime() ;
126
127 // The Step Point safety can be limited by other geometries and/or the
128 // assumptions of any process - it's not always the geometrical safety.
129 // We calculate the starting point's isotropic safety here.
130 //
131 G4ThreeVector OriginShift = startPosition - fPreviousSftOrigin ;
132 G4double MagSqShift = OriginShift.mag2() ;
133 if (MagSqShift >= sqr(fPreviousSafety)) {
134 currentSafety = 0.0 ;
135 } else {
136 currentSafety = fPreviousSafety - std::sqrt(MagSqShift) ;
137 }
138
139 // Is the monopole charged ?
140 //
141 G4double particleMagneticCharge = fParticleDef->MagneticCharge() ;
142 G4double particleElectricCharge = pParticle->GetCharge();
143
144 fGeometryLimitedStep = false ;
145 // fEndGlobalTimeComputed = false ;
146
147 // There is no need to locate the current volume. It is Done elsewhere:
148 // On track construction
149 // By the tracking, after all AlongStepDoIts, in "Relocation"
150
151 // Check whether the particle have an (EM) field force exerting upon it
152 //
153 G4FieldManager* fieldMgr = nullptr;
154 G4bool fieldExertsForce = false ;
155
156 if ((particleMagneticCharge != 0.0)) {
157 fieldMgr = fFieldPropagator->FindAndSetFieldManager(track.GetVolume());
158 if (fieldMgr != 0) {
159 // Message the field Manager, to configure it for this track
160 fieldMgr->ConfigureForTrack(&track);
161 // Moved here, in order to allow a transition
162 // from a zero-field status (with fieldMgr->(field)0
163 // to a finite field status
164
165 // If the field manager has no field, there is no field !
166 fieldExertsForce = (fieldMgr->GetDetectorField() != 0);
167 }
168 }
169
170 // Choose the calculation of the transportation: Field or not
171 //
172 if (!fieldExertsForce) {
173 G4double linearStepLength ;
174 if (fShortStepOptimisation && (currentMinimumStep <= currentSafety)) {
175 // The Step is guaranteed to be taken
176 //
177 geometryStepLength = currentMinimumStep ;
178 fGeometryLimitedStep = false ;
179 } else {
180 // Find whether the straight path intersects a volume
181 //
182 linearStepLength = fLinearNavigator->ComputeStep(startPosition,
183 startMomentumDir,
184 currentMinimumStep,
185 newSafety) ;
186 // Remember last safety origin & value.
187 //
188 fPreviousSftOrigin = startPosition ;
189 fPreviousSafety = newSafety ;
190 // fpSafetyHelper->SetCurrentSafety( newSafety, startPosition);
191
192 // The safety at the initial point has been re-calculated:
193 //
194 currentSafety = newSafety ;
195
196 fGeometryLimitedStep = (linearStepLength <= currentMinimumStep);
198 // The geometry limits the Step size (an intersection was found.)
199 geometryStepLength = linearStepLength ;
200 } else {
201 // The full Step is taken.
202 geometryStepLength = currentMinimumStep ;
203 }
204 }
205 endpointDistance = geometryStepLength ;
206
207 // Calculate final position
208 //
209 fTransportEndPosition = startPosition + geometryStepLength * startMomentumDir ;
210
211 // Momentum direction, energy and polarisation are unchanged by transport
212 //
213 fTransportEndMomentumDir = startMomentumDir ;
214 fTransportEndKineticEnergy = track.GetKineticEnergy() ;
215 fTransportEndSpin = track.GetPolarization();
216 fParticleIsLooping = false ;
217 fMomentumChanged = false ;
218 fEndGlobalTimeComputed = false ;
219 } else { // A field exerts force
220 G4double momentumMagnitude = pParticle->GetTotalMomentum() ;
221 G4ThreeVector EndUnitMomentum ;
222 G4double restMass = fParticleDef->GetPDGMass() ;
223
224 G4ChargeState chargeState(particleElectricCharge, // The charge can change (dynamic)
225 fParticleDef->GetPDGSpin(),
226 0, // Magnetic moment: pParticleDef->GetMagneticMoment(),
227 0, // Electric Dipole moment - not in Particle Definition
228 particleMagneticCharge); // in units of the positron charge
229
230 G4EquationOfMotion* equationOfMotion =
231 (fFieldPropagator->GetChordFinder()->GetIntegrationDriver()->GetStepper())
232 ->GetEquationOfMotion();
233
234 equationOfMotion
235 ->SetChargeMomentumMass(chargeState,
236 momentumMagnitude,
237 restMass) ;
238 // SetChargeMomentumMass now passes both the electric and magnetic charge - in chargeState
239
240 G4ThreeVector spin = track.GetPolarization() ;
241 G4FieldTrack aFieldTrack = G4FieldTrack(startPosition,
242 track.GetMomentumDirection(),
243 0.0,
244 track.GetKineticEnergy(),
245 restMass,
246 track.GetVelocity(),
247 track.GetGlobalTime(), // Lab.
248 track.GetProperTime(), // Part.
249 &spin) ;
250 if (currentMinimumStep > 0) {
251 // Do the Transport in the field (non recti-linear)
252 //
253 G4double lengthAlongCurve = fFieldPropagator->ComputeStep(
254 aFieldTrack,
255 currentMinimumStep,
256 currentSafety,
257 track.GetVolume()) ;
258 fGeometryLimitedStep = lengthAlongCurve < currentMinimumStep;
260 geometryStepLength = lengthAlongCurve ;
261 } else {
262 geometryStepLength = currentMinimumStep ;
263 }
264 } else {
265 geometryStepLength = 0.0 ;
266 fGeometryLimitedStep = false ;
267 }
268
269 // Remember last safety origin & value.
270 //
271 fPreviousSftOrigin = startPosition ;
272 fPreviousSafety = currentSafety ;
273 // fpSafetyHelper->SetCurrentSafety( newSafety, startPosition);
274
275 // Get the End-Position and End-Momentum (Dir-ection)
276 //
277 fTransportEndPosition = aFieldTrack.GetPosition() ;
278
279 // Momentum: Magnitude and direction can be changed too now ...
280 //
281 fMomentumChanged = true ;
282 fTransportEndMomentumDir = aFieldTrack.GetMomentumDir() ;
283
284 fTransportEndKineticEnergy = aFieldTrack.GetKineticEnergy() ;
285
286 fCandidateEndGlobalTime = aFieldTrack.GetLabTimeOfFlight();
288
289 fTransportEndSpin = aFieldTrack.GetSpin();
290 fParticleIsLooping = fFieldPropagator->IsParticleLooping() ;
291 endpointDistance = (fTransportEndPosition - startPosition).mag() ;
292 }
293
294 // If we are asked to go a step length of 0, and we are on a boundary
295 // then a boundary will also limit the step -> we must flag this.
296 //
297 if (currentMinimumStep == 0.0) {
298 if (currentSafety == 0.0) fGeometryLimitedStep = true ;
299 }
300
301 // Update the safety starting from the end-point,
302 // if it will become negative at the end-point.
303 //
304 if (currentSafety < endpointDistance) {
305 // if( particleMagneticCharge == 0.0 )
306 // G4cout << " Avoiding call to ComputeSafety : charge = 0.0 " << G4endl;
307
308 if (particleMagneticCharge != 0.0) {
309
310 G4double endSafety =
312 currentSafety = endSafety ;
314 fPreviousSafety = currentSafety ;
315 fpSafetyHelper->SetCurrentSafety(currentSafety, fTransportEndPosition);
316
317 // Because the Stepping Manager assumes it is from the start point,
318 // add the StepLength
319 //
320 currentSafety += endpointDistance ;
321
322#ifdef G4DEBUG_TRANSPORT
323 G4cout.precision(12) ;
324 G4cout << "***G4MonopoleTransportation::AlongStepGPIL ** " << G4endl ;
325 G4cout << " Called Navigator->ComputeSafety at " << fTransportEndPosition
326 << " and it returned safety= " << endSafety << G4endl ;
327 G4cout << " Adding endpoint distance " << endpointDistance
328 << " to obtain pseudo-safety= " << currentSafety << G4endl ;
329#endif
330 }
331 }
332
333 fParticleChange.ProposeTrueStepLength(geometryStepLength) ;
334
336 // change back to usual equation
337
338 return geometryStepLength ;
339}
340
341//
342// Initialize ParticleChange (by setting all its members equal
343// to corresponding members in G4Track)
344
345G4VParticleChange* G4MonopoleTransportation::AlongStepDoIt(const G4Track& track,
346 const G4Step& stepData)
347{
348#ifdef G4VERBOSE
349 static G4int noCalls = 0;
350 noCalls++;
351#endif
352
353 static const G4ParticleDefinition* fOpticalPhoton =
354 G4ParticleTable::GetParticleTable()->FindParticle("opticalphoton");
355
356 fParticleChange.Initialize(track) ;
357
358 // Code for specific process
359 //
360 fParticleChange.ProposePosition(fTransportEndPosition) ;
361 fParticleChange.ProposeMomentumDirection(fTransportEndMomentumDir) ;
363 fParticleChange.SetMomentumChanged(fMomentumChanged) ;
364
365 fParticleChange.ProposePolarization(fTransportEndSpin);
366
367 G4double deltaTime = 0.0 ;
368
369 // Calculate Lab Time of Flight (ONLY if field Equations used it!)
370 // G4double endTime = fCandidateEndGlobalTime;
371 // G4double delta_time = endTime - startTime;
372
373 G4double startTime = track.GetGlobalTime() ;
374
376 // The time was not integrated .. make the best estimate possible
377 //
378 G4double finalVelocity = track.GetVelocity() ;
379 G4double initialVelocity = stepData.GetPreStepPoint()->GetVelocity() ;
380 G4double stepLength = track.GetStepLength() ;
381
382 deltaTime = 0.0; // in case initialVelocity = 0
383 const G4DynamicParticle* fpDynamicParticle = track.GetDynamicParticle();
384 if (fpDynamicParticle->GetDefinition() == fOpticalPhoton) {
385 // A photon is in the medium of the final point
386 // during the step, so it has the final velocity.
387 deltaTime = stepLength / finalVelocity ;
388 } else if (finalVelocity > 0.0) {
389 G4double meanInverseVelocity ;
390 // deltaTime = stepLength/finalVelocity ;
391 meanInverseVelocity = 0.5
392 * (1.0 / initialVelocity + 1.0 / finalVelocity) ;
393 deltaTime = stepLength * meanInverseVelocity ;
394 } else if (initialVelocity > 0.0) {
395 deltaTime = stepLength / initialVelocity ;
396 }
397 fCandidateEndGlobalTime = startTime + deltaTime ;
398 } else {
399 deltaTime = fCandidateEndGlobalTime - startTime ;
400 }
401
402 fParticleChange.ProposeGlobalTime(fCandidateEndGlobalTime) ;
403
404 // Now Correct by Lorentz factor to get "proper" deltaTime
405
406 G4double restMass = track.GetDynamicParticle()->GetMass() ;
407 G4double deltaProperTime = deltaTime * (restMass / track.GetTotalEnergy()) ;
408
409 fParticleChange.ProposeProperTime(track.GetProperTime() + deltaProperTime) ;
410 //fParticleChange. ProposeTrueStepLength( track.GetStepLength() ) ;
411
412 // If the particle is caught looping or is stuck (in very difficult
413 // boundaries) in a magnetic field (doing many steps)
414 // THEN this kills it ...
415 //
416 if (fParticleIsLooping) {
417 G4double endEnergy = fTransportEndKineticEnergy;
418
419 if ((endEnergy < fThreshold_Important_Energy)
421 // Kill the looping particle
422 //
423 fParticleChange.ProposeTrackStatus(fStopAndKill) ;
424
425 // 'Bare' statistics
426 fSumEnergyKilled += endEnergy;
427 if (endEnergy > fMaxEnergyKilled) { fMaxEnergyKilled = endEnergy; }
428
429#ifdef G4VERBOSE
430 if ((verboseLevel > 1) ||
431 (endEnergy > fThreshold_Warning_Energy)) {
432 B2DEBUG(25, " G4MonopoleTransportation is killing track that is looping or stuck ");
433 B2DEBUG(25, " This track has " << track.GetKineticEnergy() / MeV << " MeV energy.");
434 B2DEBUG(25, " Number of trials = " << fNoLooperTrials
435 << " No of calls to AlongStepDoIt = " << noCalls);
436 }
437#endif
438 fNoLooperTrials = 0;
439 } else {
441#ifdef G4VERBOSE
442 if ((verboseLevel > 2)) {
443 B2DEBUG(25, " G4MonopoleTransportation::AlongStepDoIt(): Particle looping - "
444 << " Number of trials = " << fNoLooperTrials
445 << " No of calls to = " << noCalls);
446 }
447#endif
448 }
449 } else {
450 fNoLooperTrials = 0;
451 }
452
453 // Another (sometimes better way) is to use a user-limit maximum Step size
454 // to alleviate this problem ..
455
456 // Introduce smooth curved trajectories to particle-change
457 //
458 fParticleChange.SetPointerToVectorOfAuxiliaryPoints
459 (fFieldPropagator->GimmeTrajectoryVectorAndForgetIt());
460
461 //Monopole bounding
463 fParticleChange.ProposeTrackStatus(fStopAndKill);
464 B2DEBUG(150, "Monopole bound in " << fCurrentTouchableHandle->GetVolume()->GetName());
465 B2DEBUG(150, "with energy " << fTransportEndKineticEnergy << " MeV");
466 }
467
468 return &fParticleChange ;
469}
470
471// This ensures that the PostStep action is always called,
472// so that it can do the relocation if it is needed.
473//
474
477 G4double, // previousStepSize
478 G4ForceCondition* pForceCond)
479{
480 *pForceCond = Forced ;
481 return DBL_MAX ; // was kInfinity ; but convention now is DBL_MAX
482}
483
484G4VParticleChange* G4MonopoleTransportation::PostStepDoIt(const G4Track& track,
485 const G4Step&)
486{
487 G4TouchableHandle retCurrentTouchable ; // The one to return
488
489 // Initialize ParticleChange (by setting all its members equal
490 // to corresponding members in G4Track)
491 // fParticleChange.Initialize(track) ; // To initialise TouchableChange
492
493 fParticleChange.ProposeTrackStatus(track.GetTrackStatus()) ;
494
495 // If the Step was determined by the volume boundary,
496 // logically relocate the particle
497
499 // fCurrentTouchable will now become the previous touchable,
500 // and what was the previous will be freed.
501 // (Needed because the preStepPoint can point to the previous touchable)
502
503 fLinearNavigator->SetGeometricallyLimitedStep() ;
505 LocateGlobalPointAndUpdateTouchableHandle(track.GetPosition(),
506 track.GetMomentumDirection(),
508 true) ;
509 // Check whether the particle is out of the world volume
510 // If so it has exited and must be killed.
511 //
512 if (fCurrentTouchableHandle->GetVolume() == 0) {
513 fParticleChange.ProposeTrackStatus(fStopAndKill) ;
514 }
515 retCurrentTouchable = fCurrentTouchableHandle ;
516 fParticleChange.SetTouchableHandle(fCurrentTouchableHandle) ;
517 } else { // fGeometryLimitedStep is false
518 // This serves only to move the Navigator's location
519 //
520 fLinearNavigator->LocateGlobalPointWithinVolume(track.GetPosition()) ;
521
522 // The value of the track's current Touchable is retained.
523 // (and it must be correct because we must use it below to
524 // overwrite the (unset) one in particle change)
525 // It must be fCurrentTouchable too ??
526 //
527 fParticleChange.SetTouchableHandle(track.GetTouchableHandle()) ;
528 retCurrentTouchable = track.GetTouchableHandle() ;
529 } // endif ( fGeometryLimitedStep )
530
531 const G4VPhysicalVolume* pNewVol = retCurrentTouchable->GetVolume() ;
532 G4Material* pNewMaterial = nullptr ;
533 G4VSensitiveDetector* pNewSensitiveDetector = nullptr ;
534
535 if (pNewVol != 0) {
536 pNewMaterial = pNewVol->GetLogicalVolume()->GetMaterial();
537 pNewSensitiveDetector = pNewVol->GetLogicalVolume()->GetSensitiveDetector();
538 }
539
540 fParticleChange.SetMaterialInTouchable(
541 pNewMaterial) ;
542 fParticleChange.SetSensitiveDetectorInTouchable(
543 pNewSensitiveDetector) ;
544
545 const G4MaterialCutsCouple* pNewMaterialCutsCouple = 0;
546 if (pNewVol != 0) {
547 pNewMaterialCutsCouple = pNewVol->GetLogicalVolume()->GetMaterialCutsCouple();
548 }
549
550 if (pNewVol != 0 && pNewMaterialCutsCouple != 0 &&
551 pNewMaterialCutsCouple->GetMaterial() != pNewMaterial) {
552 // for parametrized volume
553 //
554 pNewMaterialCutsCouple =
555 G4ProductionCutsTable::GetProductionCutsTable()
556 ->GetMaterialCutsCouple(pNewMaterial,
557 pNewMaterialCutsCouple->GetProductionCuts());
558 }
559 fParticleChange.SetMaterialCutsCoupleInTouchable(pNewMaterialCutsCouple);
560
561 // temporarily until Get/Set Material of ParticleChange,
562 // and StepPoint can be made const.
563 // Set the touchable in ParticleChange
564 // this must always be done because the particle change always
565 // uses this value to overwrite the current touchable pointer.
566 //
567 fParticleChange.SetTouchableHandle(retCurrentTouchable) ;
568
569 return &fParticleChange ;
570}
571
572// New method takes over the responsibility to reset the state
573// of G4MonopoleTransportation object at the start of a new track
574// or the resumption of a suspended track.
575
576void
578{
579 G4VProcess::StartTracking(aTrack);
580
581// The actions here are those that were taken in AlongStepGPIL
582// when track.GetCurrentStepNumber()==1
583
584 // reset safety value and center
585 //
586 fPreviousSafety = 0.0 ;
587 fPreviousSftOrigin = G4ThreeVector(0., 0., 0.) ;
588
589 // reset looping counter -- for motion in field
590 fNoLooperTrials = 0;
591 // Must clear this state .. else it depends on last track's value
592 // --> a better solution would set this from state of suspended track TODO ?
593 // Was if( aTrack->GetCurrentStepNumber()==1 ) { .. }
594
595 // ChordFinder reset internal state
596 //
597 if (DoesGlobalFieldExist()) {
598 fFieldPropagator->ClearPropagatorState();
599 // Resets all state of field propagator class (ONLY)
600 // including safety values (in case of overlaps and to wipe for first track).
601
602 // G4ChordFinder* chordF= fFieldPropagator->GetChordFinder();
603 // if( chordF ) chordF->ResetStepEstimate();
604 }
605
606 // Make sure to clear the chord finders of all fields (ie managers)
607 static G4FieldManagerStore* fieldMgrStore = G4FieldManagerStore::GetInstance();
608 fieldMgrStore->ClearAllChordFindersState();
609
610 // Update the current touchable handle (from the track's)
611 //
612 fCurrentTouchableHandle = aTrack->GetTouchableHandle();
613}
void SwitchChordFinder(G4int val)
Switches chord finder between 1 - basf2 FullSim chord finder 2 - monopole chord finder Since monopole...
static G4MonopoleFieldSetup * GetMonopoleFieldSetup()
Returns G4MonopoleFieldSetup instance.
G4TouchableHandle fCurrentTouchableHandle
Current touchable handle.
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData)
G4VProcess::PostStepDoIt() implementation.
G4double fCandidateEndGlobalTime
The particle's state after this Step, Store for DoIt.
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
G4VProcess::AlongStepGetPhysicalInteractionLength() implementation.
G4int fThresholdTrials
Nubmer of trials for looping particles.
G4double fMaxEnergyKilled
Max of abandoned looping tracks energies.
G4bool fGeometryLimitedStep
Flag to determine whether a boundary was reached.
G4ThreeVector fTransportEndMomentumDir
The particle's state after this Step, Store for DoIt.
G4double fPreviousSafety
Remember last safety value.
virtual G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData)
G4VProcess::AlongStepDoIt() implementation, Proposes changes during step to fParticleChange of this c...
G4SafetyHelper * fpSafetyHelper
To pass it the safety value obtained.
G4double fSumEnergyKilled
Sum of abandoned looping tracks energies.
G4ThreeVector fPreviousSftOrigin
Remember last safety origin.
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *pForceCond)
G4VProcess::PostStepGetPhysicalInteractionLength() implementation.
G4ThreeVector fTransportEndSpin
The particle's state after this Step, Store for DoIt.
G4Navigator * fLinearNavigator
Propagator used to transport the particle.
G4double fThreshold_Important_Energy
Hesitate above this about looping particle for a certain no of trials.
G4ThreeVector fTransportEndPosition
The particle's state after this Step, Store for DoIt.
G4bool fEndGlobalTimeComputed
The particle's state after this Step, Store for DoIt.
G4bool fShortStepOptimisation
Whether to avoid calling G4Navigator for short step ( < safety) If using it, the safety estimate for ...
G4int fNoLooperTrials
Counter for steps in which particle reports 'looping', if it is above 'Important' Energy.
G4bool fMomentumChanged
The particle's state after this Step, Store for DoIt.
virtual void StartTracking(G4Track *aTrack)
Reset state for new (potentially resumed) track.
G4ParticleChangeForTransport fParticleChange
New ParticleChange.
G4MonopoleTransportation(const G4Monopole *mpl, G4int verb=1)
Constructor.
G4double fThreshold_Warning_Energy
Warn above this energy about looping particle.
G4bool DoesGlobalFieldExist()
Checks whether a field exists for the "global" field manager.
const G4Monopole * fParticleDef
Monopole definition for charge and mass reference.
G4bool fParticleIsLooping
Is the monopole stuck in looping.
G4double fThreshold_Trap_Energy
Assume monopoles below this can bound to material.
G4PropagatorInField * fFieldPropagator
Propagator used to transport the particle.
G4double fTransportEndKineticEnergy
The particle's state after this Step, Store for DoIt.
G4MonopoleFieldSetup * fMagSetup
Monpole field setup.
A class to hold monopole description as a particle.
Definition: G4Monopole.h:33
G4double MagneticCharge() const
Returns magnetic charge of the monopole.
Definition: G4Monopole.cc:47
Abstract base class for different kinds of events.
STL namespace.