Geant4 VMC Version 6.6
Loading...
Searching...
No Matches
G4MonopoleTransportation.cxx
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
28//
29// $Id: G4MonopoleTransportation.cc 111448 2018-08-10 07:54:47Z gcosmo $
30//
31//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
32//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
33//
34// This class is a process responsible for the transportation of
35// magnetic monopoles, ie the geometrical propagation that encounters the
36// geometrical sub-volumes of the detectors.
37//
38// For monopoles, uses a different equation of motion and ignores energy
39// conservation.
40//
41
42// =======================================================================
43// Created: 3 May 2010, J. Apostolakis, B. Bozsogi
44// =======================================================================
45
47#include "G4ChordFinder.hh"
48#include "G4FieldManagerStore.hh"
49#include "G4Monopole.hh"
50#include "G4ParticleTable.hh"
51#include "G4ProductionCutsTable.hh"
52#include "G4SafetyHelper.hh"
53#include "G4SystemOfUnits.hh"
54#include "G4TransportationProcessType.hh"
55
57
58//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
59
61 const G4Monopole* mpl, G4int verb)
62 : G4VProcess(G4String("MonopoleTransportation"), fTransportation),
63 fParticleDef(mpl),
64 fMagSetup(0),
65 fLinearNavigator(0),
66 fFieldPropagator(0),
67 fParticleIsLooping(false),
68 fPreviousSftOrigin(0., 0., 0.),
69 fPreviousSafety(0.0),
70 fThreshold_Warning_Energy(100 * MeV),
71 fThreshold_Important_Energy(250 * MeV),
72 fThresholdTrials(10),
73 // fUnimportant_Energy( 1 * MeV ),
74 fNoLooperTrials(0),
75 fSumEnergyKilled(0.0),
76 fMaxEnergyKilled(0.0),
77 fShortStepOptimisation(false), // Old default: true (=fast short steps)
78 fpSafetyHelper(0),
79 noCalls(0)
80{
81 verboseLevel = verb;
82
83 // set Process Sub Type
84 SetProcessSubType(TRANSPORTATION);
85
86 // Do not finalize the G4MonopoleTransportation class
87 if (G4Threading::IsMasterThread() &&
88 G4Threading::IsMultithreadedApplication()) {
89 return;
90 }
91
93
94 G4TransportationManager* transportMgr =
95 G4TransportationManager::GetTransportationManager();
96
97 fLinearNavigator = transportMgr->GetNavigatorForTracking();
98
99 fFieldPropagator = transportMgr->GetPropagatorInField();
100 fpSafetyHelper = transportMgr->GetSafetyHelper();
101
102 // New
103
104 // Cannot determine whether a field exists here,
105 // because it would only work if the field manager has informed
106 // about the detector's field before this transportation process
107 // is constructed.
108 // Instead later the method DoesGlobalFieldExist() is called
109 fCurrentTouchableHandle = nullptr;
110
113}
114
115//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
116
118{
119 if ((verboseLevel > 0) && (fSumEnergyKilled > 0.0)) {
120 G4cout << " G4MonopoleTransportation: Statistics for looping particles "
121 << G4endl;
122 G4cout << " Sum of energy of loopers killed: " << fSumEnergyKilled
123 << G4endl;
124 G4cout << " Max energy of loopers killed: " << fMaxEnergyKilled << G4endl;
125 }
126}
127
128//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
129//
130// Responsibilities:
131// Find whether the geometry limits the Step, and to what length
132// Calculate the new value of the safety and return it.
133// Store the final time, position and momentum.
134
136 const G4Track& track,
137 G4double, // previousStepSize
138 G4double currentMinimumStep, G4double& currentSafety,
139 G4GPILSelection* selection)
140{
142 // change to monopole equation
143
144 G4double geometryStepLength, newSafety;
145 fParticleIsLooping = false;
146
147 // Initial actions moved to StartTrack()
148 // --------------------------------------
149 // Note: in case another process changes touchable handle
150 // it will be necessary to add here (for all steps)
151 // fCurrentTouchableHandle = aTrack->GetTouchableHandle();
152
153 // GPILSelection is set to defaule value of CandidateForSelection
154 // It is a return value
155 //
156 *selection = CandidateForSelection;
157
158 // Get initial Energy/Momentum of the track
159 //
160 const G4DynamicParticle* pParticle = track.GetDynamicParticle();
161 G4ThreeVector startMomentumDir = pParticle->GetMomentumDirection();
162 G4ThreeVector startPosition = track.GetPosition();
163
164 // G4double theTime = track.GetGlobalTime() ;
165
166 // The Step Point safety can be limited by other geometries and/or the
167 // assumptions of any process - it's not always the geometrical safety.
168 // We calculate the starting point's isotropic safety here.
169 //
170 G4ThreeVector OriginShift = startPosition - fPreviousSftOrigin;
171 G4double MagSqShift = OriginShift.mag2();
172 if (MagSqShift >= sqr(fPreviousSafety)) {
173 currentSafety = 0.0;
174 }
175 else {
176 currentSafety = fPreviousSafety - std::sqrt(MagSqShift);
177 }
178
179 // Is the monopole charged ?
180 //
181 G4double particleMagneticCharge = fParticleDef->MagneticCharge();
182 G4double particleElectricCharge = pParticle->GetCharge();
183
184 fGeometryLimitedStep = false;
185 // fEndGlobalTimeComputed = false ;
186
187 // There is no need to locate the current volume. It is Done elsewhere:
188 // On track construction
189 // By the tracking, after all AlongStepDoIts, in "Relocation"
190
191 // Check whether the particle have an (EM) field force exerting upon it
192 //
193 G4FieldManager* fieldMgr = 0;
194 G4bool fieldExertsForce = false;
195
196 if ((particleMagneticCharge != 0.0)) {
197 fieldMgr = fFieldPropagator->FindAndSetFieldManager(track.GetVolume());
198 if (fieldMgr != 0) {
199 // Message the field Manager, to configure it for this track
200 fieldMgr->ConfigureForTrack(&track);
201 // Moved here, in order to allow a transition
202 // from a zero-field status (with fieldMgr->(field)0
203 // to a finite field status
204
205 // If the field manager has no field, there is no field !
206 fieldExertsForce = (fieldMgr->GetDetectorField() != 0);
207 }
208 }
209
210 // G4cout << " G4Transport: field exerts force= " << fieldExertsForce
211 // << " fieldMgr= " << fieldMgr << G4endl;
212
213 // Choose the calculation of the transportation: Field or not
214 //
215 if (!fieldExertsForce) {
216 G4double linearStepLength;
217 if (fShortStepOptimisation && (currentMinimumStep <= currentSafety)) {
218 // The Step is guaranteed to be taken
219 //
220 geometryStepLength = currentMinimumStep;
221 fGeometryLimitedStep = false;
222 }
223 else {
224 // Find whether the straight path intersects a volume
225 //
226 linearStepLength = fLinearNavigator->ComputeStep(
227 startPosition, startMomentumDir, currentMinimumStep, newSafety);
228 // Remember last safety origin & value.
229 //
230 fPreviousSftOrigin = startPosition;
231 fPreviousSafety = newSafety;
232 // fpSafetyHelper->SetCurrentSafety( newSafety, startPosition);
233
234 // The safety at the initial point has been re-calculated:
235 //
236 currentSafety = newSafety;
237
238 fGeometryLimitedStep = (linearStepLength <= currentMinimumStep);
240 // The geometry limits the Step size (an intersection was found.)
241 geometryStepLength = linearStepLength;
242 }
243 else {
244 // The full Step is taken.
245 geometryStepLength = currentMinimumStep;
246 }
247 }
248 endpointDistance = geometryStepLength;
249
250 // Calculate final position
251 //
253 startPosition + geometryStepLength * startMomentumDir;
254
255 // Momentum direction, energy and polarisation are unchanged by transport
256 //
257 fTransportEndMomentumDir = startMomentumDir;
258 fTransportEndKineticEnergy = track.GetKineticEnergy();
259 fTransportEndSpin = track.GetPolarization();
260 fParticleIsLooping = false;
261 fMomentumChanged = false;
263 }
264 else // A field exerts force
265 {
266 G4double momentumMagnitude = pParticle->GetTotalMomentum();
267 G4ThreeVector EndUnitMomentum;
268 G4double lengthAlongCurve;
269 G4double restMass = fParticleDef->GetPDGMass();
270
271 G4ChargeState chargeState(particleElectricCharge, // The charge can change
272 fParticleDef->GetPDGSpin(),
273 0, // Magnetic moment:
274 0, // Electric Dipole moment
275 particleMagneticCharge); // in Mev/c
276
277 G4EquationOfMotion* equationOfMotion = fFieldPropagator->GetChordFinder()
278 ->GetIntegrationDriver()
279 ->GetEquationOfMotion();
280
281 equationOfMotion->SetChargeMomentumMass(
282 chargeState, momentumMagnitude, restMass);
283 // SetChargeMomentumMass now passes both the electric and magnetic
284 // charge in chargeState
285
286 G4ThreeVector spin = track.GetPolarization();
287 G4FieldTrack aFieldTrack =
288 G4FieldTrack(startPosition, track.GetMomentumDirection(), 0.0,
289 track.GetKineticEnergy(), restMass, track.GetVelocity(),
290 track.GetGlobalTime(), // Lab.
291 track.GetProperTime(), // Part.
292 &spin);
293 if (currentMinimumStep > 0) {
294 // Do the Transport in the field (non recti-linear)
295 //
296 lengthAlongCurve = fFieldPropagator->ComputeStep(
297 aFieldTrack, currentMinimumStep, currentSafety, track.GetVolume());
298 fGeometryLimitedStep = lengthAlongCurve < currentMinimumStep;
300 geometryStepLength = lengthAlongCurve;
301 }
302 else {
303 geometryStepLength = currentMinimumStep;
304 }
305 }
306 else {
307 geometryStepLength = lengthAlongCurve = 0.0;
308 fGeometryLimitedStep = false;
309 }
310
311 // Remember last safety origin & value.
312 //
313 fPreviousSftOrigin = startPosition;
314 fPreviousSafety = currentSafety;
315 // fpSafetyHelper->SetCurrentSafety( newSafety, startPosition);
316
317 // Get the End-Position and End-Momentum (Dir-ection)
318 //
319 fTransportEndPosition = aFieldTrack.GetPosition();
320
321 // Momentum: Magnitude and direction can be changed too now ...
322 //
323 fMomentumChanged = true;
324 fTransportEndMomentumDir = aFieldTrack.GetMomentumDir();
325
326 fTransportEndKineticEnergy = aFieldTrack.GetKineticEnergy();
327
328 fCandidateEndGlobalTime = aFieldTrack.GetLabTimeOfFlight();
330
331 fTransportEndSpin = aFieldTrack.GetSpin();
332 fParticleIsLooping = fFieldPropagator->IsParticleLooping();
333 endpointDistance = (fTransportEndPosition - startPosition).mag();
334 }
335
336 // If we are asked to go a step length of 0, and we are on a boundary
337 // then a boundary will also limit the step -> we must flag this.
338 //
339 if (currentMinimumStep == 0.0) {
340 if (currentSafety == 0.0) fGeometryLimitedStep = true;
341 }
342
343 // Update the safety starting from the end-point,
344 // if it will become negative at the end-point.
345 //
346 if (currentSafety < endpointDistance) {
347 // if( particleMagneticCharge == 0.0 )
348 // G4cout << " Avoiding call to ComputeSafety : charge = 0.0 " << G4endl;
349
350 if (particleMagneticCharge != 0.0) {
351
352 G4double endSafety =
354 currentSafety = endSafety;
356 fPreviousSafety = currentSafety;
357 fpSafetyHelper->SetCurrentSafety(currentSafety, fTransportEndPosition);
358
359 // Because the Stepping Manager assumes it is from the start point,
360 // add the StepLength
361 //
362 currentSafety += endpointDistance;
363
364#ifdef G4DEBUG_TRANSPORT
365 G4cout.precision(12);
366 G4cout << "***G4MonopoleTransportation::AlongStepGPIL ** " << G4endl;
367 G4cout << " Called Navigator->ComputeSafety at " << fTransportEndPosition
368 << " and it returned safety= " << endSafety << G4endl;
369 G4cout << " Adding endpoint distance " << endpointDistance
370 << " to obtain pseudo-safety= " << currentSafety << G4endl;
371#endif
372 }
373 }
374
375 fParticleChange.ProposeTrueStepLength(geometryStepLength);
376
378 // change back to usual equation
379
380 return geometryStepLength;
381}
382
383//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
384//
385// Initialize ParticleChange (by setting all its members equal
386// to corresponding members in G4Track)
387
389 const G4Track& track, const G4Step& stepData)
390{
391 const G4ParticleDefinition* fOpticalPhoton =
392 G4ParticleTable::GetParticleTable()->FindParticle("opticalphoton");
393
394 ++noCalls;
395
396 fParticleChange.Initialize(track);
397
398 // Code for specific process
399 //
400 fParticleChange.ProposePosition(fTransportEndPosition);
401 fParticleChange.ProposeMomentumDirection(fTransportEndMomentumDir);
403 fParticleChange.SetMomentumChanged(fMomentumChanged);
404
405 fParticleChange.ProposePolarization(fTransportEndSpin);
406
407 G4double deltaTime = 0.0;
408
409 // Calculate Lab Time of Flight (ONLY if field Equations used it!)
410
411 G4double startTime = track.GetGlobalTime();
412
414 // The time was not integrated .. make the best estimate possible
415 //
416 G4double finalVelocity = track.GetVelocity();
417 G4double initialVelocity = stepData.GetPreStepPoint()->GetVelocity();
418 G4double stepLength = track.GetStepLength();
419
420 deltaTime = 0.0; // in case initialVelocity = 0
421 const G4DynamicParticle* fpDynamicParticle = track.GetDynamicParticle();
422 if (fpDynamicParticle->GetDefinition() == fOpticalPhoton) {
423 // A photon is in the medium of the final point
424 // during the step, so it has the final velocity.
425 deltaTime = stepLength / finalVelocity;
426 }
427 else if (finalVelocity > 0.0) {
428 G4double meanInverseVelocity;
429 // deltaTime = stepLength/finalVelocity ;
430 meanInverseVelocity = 0.5 * (1.0 / initialVelocity + 1.0 / finalVelocity);
431 deltaTime = stepLength * meanInverseVelocity;
432
433 // G4cout << "finalVelocity > 0.0 gives deltaTime " << deltaTime <<
434 // G4endl;
435 }
436 else if (initialVelocity > 0.0) {
437 deltaTime = stepLength / initialVelocity;
438 }
439 fCandidateEndGlobalTime = startTime + deltaTime;
440 // G4cout << "not EndGlobalTimeComputed" << G4endl;
441 // G4cout << "deltaTime, fCandidateEndGlobalTime "
442 // << deltaTime << ", " << fCandidateEndGlobalTime << G4endl;
443 }
444 else {
445 deltaTime = fCandidateEndGlobalTime - startTime;
446 // G4cout << "is EndGlobalTimeComputed" << G4endl;
447 // G4cout << "deltaTime, fCandidateEndGlobalTime "
448 // << deltaTime << ", " << fCandidateEndGlobalTime << G4endl;
449 }
450
451 fParticleChange.ProposeGlobalTime(fCandidateEndGlobalTime);
452
453 // Now Correct by Lorentz factor to get "proper" deltaTime
454
455 G4double restMass = track.GetDynamicParticle()->GetMass();
456 G4double deltaProperTime = deltaTime * (restMass / track.GetTotalEnergy());
457
458 fParticleChange.ProposeProperTime(track.GetProperTime() + deltaProperTime);
459 // G4cout << "ProposeProperTime: " << track.GetProperTime() + deltaProperTime
460 // << G4endl;
461
462 // If the particle is caught looping or is stuck (in very difficult
463 // boundaries) in a magnetic field (doing many steps)
464 // THEN this kills it ...
465 //
466 if (fParticleIsLooping) {
467 G4double endEnergy = fTransportEndKineticEnergy;
468
469 if ((endEnergy < fThreshold_Important_Energy) ||
471 // Kill the looping particle
472 //
473 fParticleChange.ProposeTrackStatus(fStopAndKill);
474
475 // 'Bare' statistics
476 fSumEnergyKilled += endEnergy;
477 if (endEnergy > fMaxEnergyKilled) {
478 fMaxEnergyKilled = endEnergy;
479 }
480
481#ifdef G4VERBOSE
482 if ((verboseLevel > 1) || (endEnergy > fThreshold_Warning_Energy)) {
483 G4cout << " G4MonopoleTransportation is killing track "
484 << "that is looping or stuck " << G4endl << " This track has "
485 << track.GetKineticEnergy() / MeV << " MeV energy." << G4endl;
486 G4cout << " Number of trials = " << fNoLooperTrials
487 << " No of calls to AlongStepDoIt = " << noCalls << G4endl;
488 }
489#endif
490 fNoLooperTrials = 0;
491 }
492 else {
494#ifdef G4VERBOSE
495 if ((verboseLevel > 2)) {
496 G4cout << " G4MonopoleTransportation::AlongStepDoIt(): "
497 << "Particle looping - "
498 << " Number of trials = " << fNoLooperTrials
499 << " No of calls to = " << noCalls << G4endl;
500 }
501#endif
502 }
503 }
504 else {
505 fNoLooperTrials = 0;
506 }
507
508 // Another (sometimes better way) is to use a user-limit maximum Step size
509 // to alleviate this problem ..
510
511 // Introduce smooth curved trajectories to particle-change
512 //
513 fParticleChange.SetPointerToVectorOfAuxiliaryPoints(
514 fFieldPropagator->GimmeTrajectoryVectorAndForgetIt());
515
516 return &fParticleChange;
517}
518
519//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
520//
521// This ensures that the PostStep action is always called,
522// so that it can do the relocation if it is needed.
523//
524
526 const G4Track&,
527 G4double, // previousStepSize
528 G4ForceCondition* pForceCond)
529{
530 *pForceCond = Forced;
531 return DBL_MAX; // was kInfinity ; but convention now is DBL_MAX
532}
533
534//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
535
537 const G4Track& track, const G4Step&)
538{
539 G4TouchableHandle retCurrentTouchable; // The one to return
540
541 // Initialize ParticleChange (by setting all its members equal
542 // to corresponding members in G4Track)
543 // fParticleChange.Initialize(track) ; // To initialise TouchableChange
544
545 fParticleChange.ProposeTrackStatus(track.GetTrackStatus());
546
547 // If the Step was determined by the volume boundary,
548 // logically relocate the particle
549
551 // fCurrentTouchable will now become the previous touchable,
552 // and what was the previous will be freed.
553 // (Needed because the preStepPoint can point to the previous touchable)
554
555 fLinearNavigator->SetGeometricallyLimitedStep();
556 fLinearNavigator->LocateGlobalPointAndUpdateTouchableHandle(
557 track.GetPosition(), track.GetMomentumDirection(),
559 // Check whether the particle is out of the world volume
560 // If so it has exited and must be killed.
561 //
562 if (fCurrentTouchableHandle->GetVolume() == 0) {
563 fParticleChange.ProposeTrackStatus(fStopAndKill);
564 }
565 retCurrentTouchable = fCurrentTouchableHandle;
566 fParticleChange.SetTouchableHandle(fCurrentTouchableHandle);
567 }
568 else // fGeometryLimitedStep is false
569 {
570 // This serves only to move the Navigator's location
571 //
572 fLinearNavigator->LocateGlobalPointWithinVolume(track.GetPosition());
573
574 // The value of the track's current Touchable is retained.
575 // (and it must be correct because we must use it below to
576 // overwrite the (unset) one in particle change)
577 // It must be fCurrentTouchable too ??
578 //
579 fParticleChange.SetTouchableHandle(track.GetTouchableHandle());
580 retCurrentTouchable = track.GetTouchableHandle();
581 } // endif ( fGeometryLimitedStep )
582
583 const G4VPhysicalVolume* pNewVol = retCurrentTouchable->GetVolume();
584 const G4Material* pNewMaterial = 0;
585 const G4VSensitiveDetector* pNewSensitiveDetector = 0;
586 if (pNewVol != 0) {
587 pNewMaterial = pNewVol->GetLogicalVolume()->GetMaterial();
588 pNewSensitiveDetector = pNewVol->GetLogicalVolume()->GetSensitiveDetector();
589 }
590
591 // ( <const_cast> pNewMaterial ) ;
592 // ( <const_cast> pNewSensitiveDetector) ;
593
594 fParticleChange.SetMaterialInTouchable((G4Material*)pNewMaterial);
595 fParticleChange.SetSensitiveDetectorInTouchable(
596 (G4VSensitiveDetector*)pNewSensitiveDetector);
597
598 const G4MaterialCutsCouple* pNewMaterialCutsCouple = 0;
599 if (pNewVol != 0) {
600 pNewMaterialCutsCouple =
601 pNewVol->GetLogicalVolume()->GetMaterialCutsCouple();
602 }
603
604 if (pNewVol != 0 && pNewMaterialCutsCouple != 0 &&
605 pNewMaterialCutsCouple->GetMaterial() != pNewMaterial) {
606 // for parametrized volume
607 //
608 pNewMaterialCutsCouple =
609 G4ProductionCutsTable::GetProductionCutsTable()->GetMaterialCutsCouple(
610 pNewMaterial, pNewMaterialCutsCouple->GetProductionCuts());
611 }
612 fParticleChange.SetMaterialCutsCoupleInTouchable(pNewMaterialCutsCouple);
613
614 // temporarily until Get/Set Material of ParticleChange,
615 // and StepPoint can be made const.
616 // Set the touchable in ParticleChange
617 // this must always be done because the particle change always
618 // uses this value to overwrite the current touchable pointer.
619 //
620 fParticleChange.SetTouchableHandle(retCurrentTouchable);
621
622 return &fParticleChange;
623}
624
625//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
626
627// New method takes over the responsibility to reset the state
628// of G4MonopoleTransportation object at the start of a new track
629// or the resumption of a suspended track.
630
632{
633 G4VProcess::StartTracking(aTrack);
634
635 // The actions here are those that were taken in AlongStepGPIL
636 // when track.GetCurrentStepNumber()==1
637
638 // reset safety value and center
639 //
640 fPreviousSafety = 0.0;
641 fPreviousSftOrigin = G4ThreeVector(0., 0., 0.);
642
643 // reset looping counter -- for motion in field
644 fNoLooperTrials = 0;
645 // Must clear this state .. else it depends on last track's value
646 // --> a better solution would set this from state of suspended track TODO ?
647 // Was if( aTrack->GetCurrentStepNumber()==1 ) { .. }
648
649 // ChordFinder reset internal state
650 //
651 if (DoesGlobalFieldExist()) {
652 fFieldPropagator->ClearPropagatorState();
653 // Resets all state of field propagator class (ONLY)
654 // including safety values
655 // in case of overlaps and to wipe for first track).
656
657 // G4ChordFinder* chordF= fFieldPropagator->GetChordFinder();
658 // if( chordF ) chordF->ResetStepEstimate();
659 }
660
661 // Make sure to clear the chord finders of all fields (ie managers)
662 G4FieldManagerStore* fieldMgrStore = G4FieldManagerStore::GetInstance();
663 fieldMgrStore->ClearAllChordFindersState();
664
665 // Update the current touchable handle (from the track's)
666 //
667 fCurrentTouchableHandle = aTrack->GetTouchableHandle();
668}
669
670//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Definition of the G4MonopoleTransportation class.
Definition of the G4Monopole class.
void SetStepperAndChordFinder(G4int val)
static G4MonopoleFieldSetup * GetMonopoleFieldSetup()
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
G4ParticleChangeForTransport fParticleChange
virtual G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData)
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4ForceCondition *pForceCond)
G4MonopoleTransportation(const G4Monopole *p, G4int verbosityLevel=1)
G4PropagatorInField * fFieldPropagator
virtual void StartTracking(G4Track *aTrack)
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData)
G4double MagneticCharge() const