Geant4  10.00.p01
G4MonopoleTransportation.cc
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 69705 2013-05-13 09:09:52Z 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 "G4ProductionCutsTable.hh"
48 #include "G4ParticleTable.hh"
49 #include "G4ChordFinder.hh"
50 #include "G4SafetyHelper.hh"
51 #include "G4FieldManagerStore.hh"
52 #include "G4Monopole.hh"
54 #include "G4SystemOfUnits.hh"
55 
57 
58 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
59 
61  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 ), fMaxEnergyKilled( 0.0 ),
76  fShortStepOptimisation(false), // Old default: true (=fast short steps)
77  fpSafetyHelper(0)
78 {
79  verboseLevel = verb;
80 
81  // set Process Sub Type
83 
85 
86  G4TransportationManager* transportMgr ;
87 
89 
90  fLinearNavigator = transportMgr->GetNavigatorForTracking() ;
91 
92  // fGlobalFieldMgr = transportMgr->GetFieldManager() ;
93 
94  fFieldPropagator = transportMgr->GetPropagatorInField() ;
95 
96  fpSafetyHelper = transportMgr->GetSafetyHelper(); // New
97 
98  // Cannot determine whether a field exists here,
99  // because it would only work if the field manager has informed
100  // about the detector's field before this transportation process
101  // is constructed.
102  // Instead later the method DoesGlobalFieldExist() is called
103 
104  static G4TouchableHandle nullTouchableHandle; // Points to (G4VTouchable*) 0
105  fCurrentTouchableHandle = nullTouchableHandle;
106 
107  fEndGlobalTimeComputed = false;
109 }
110 
111 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
112 
114 {
115  if( (verboseLevel > 0) && (fSumEnergyKilled > 0.0 ) ){
116  G4cout << " G4MonopoleTransportation: Statistics for looping particles "
117  << G4endl;
118  G4cout << " Sum of energy of loopers killed: " << fSumEnergyKilled << G4endl;
119  G4cout << " Max energy of loopers killed: " << fMaxEnergyKilled << G4endl;
120  }
121 }
122 
123 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
124 //
125 // Responsibilities:
126 // Find whether the geometry limits the Step, and to what length
127 // Calculate the new value of the safety and return it.
128 // Store the final time, position and momentum.
129 
132  G4double, // previousStepSize
133  G4double currentMinimumStep,
134  G4double& currentSafety,
135  G4GPILSelection* selection )
136 {
138  // change to monopole equation
139 
140  G4double geometryStepLength, newSafety ;
141  fParticleIsLooping = false ;
142 
143  // Initial actions moved to StartTrack()
144  // --------------------------------------
145  // Note: in case another process changes touchable handle
146  // it will be necessary to add here (for all steps)
147  // fCurrentTouchableHandle = aTrack->GetTouchableHandle();
148 
149  // GPILSelection is set to defaule value of CandidateForSelection
150  // It is a return value
151  //
152  *selection = CandidateForSelection ;
153 
154  // Get initial Energy/Momentum of the track
155  //
156  const G4DynamicParticle* pParticle = track.GetDynamicParticle() ;
157  G4ThreeVector startMomentumDir = pParticle->GetMomentumDirection() ;
159 
160  // G4double theTime = track.GetGlobalTime() ;
161 
162  // The Step Point safety can be limited by other geometries and/or the
163  // assumptions of any process - it's not always the geometrical safety.
164  // We calculate the starting point's isotropic safety here.
165  //
166  G4ThreeVector OriginShift = startPosition - fPreviousSftOrigin ;
167  G4double MagSqShift = OriginShift.mag2() ;
168  if( MagSqShift >= sqr(fPreviousSafety) )
169  {
170  currentSafety = 0.0 ;
171  }
172  else
173  {
174  currentSafety = fPreviousSafety - std::sqrt(MagSqShift) ;
175  }
176 
177  // Is the monopole charged ?
178  //
179  G4double particleMagneticCharge = fParticleDef->MagneticCharge() ;
180  G4double particleElectricCharge = pParticle->GetCharge();
181 
182  fGeometryLimitedStep = false ;
183  // fEndGlobalTimeComputed = 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 whether the particle have an (EM) field force exerting upon it
190  //
191  G4FieldManager* fieldMgr=0;
192  G4bool fieldExertsForce = false ;
193 
194  if( (particleMagneticCharge != 0.0) )
195  {
196  fieldMgr= fFieldPropagator->FindAndSetFieldManager( track.GetVolume() );
197  if (fieldMgr != 0) {
198  // Message the field Manager, to configure it for this track
199  fieldMgr->ConfigureForTrack( &track );
200  // Moved here, in order to allow a transition
201  // from a zero-field status (with fieldMgr->(field)0
202  // to a finite field status
203 
204  // If the field manager has no field, there is no field !
205  fieldExertsForce = (fieldMgr->GetDetectorField() != 0);
206  }
207  }
208 
209  // G4cout << " G4Transport: field exerts force= " << fieldExertsForce
210  // << " fieldMgr= " << fieldMgr << G4endl;
211 
212  // Choose the calculation of the transportation: Field or not
213  //
214  if( !fieldExertsForce )
215  {
216  G4double linearStepLength ;
217  if( fShortStepOptimisation && (currentMinimumStep <= currentSafety) )
218  {
219  // The Step is guaranteed to be taken
220  //
221  geometryStepLength = currentMinimumStep ;
222  fGeometryLimitedStep = false ;
223  }
224  else
225  {
226  // Find whether the straight path intersects a volume
227  //
228  linearStepLength = fLinearNavigator->ComputeStep( startPosition,
229  startMomentumDir,
230  currentMinimumStep,
231  newSafety) ;
232  // Remember last safety origin & value.
233  //
234  fPreviousSftOrigin = startPosition ;
235  fPreviousSafety = newSafety ;
236  // fpSafetyHelper->SetCurrentSafety( newSafety, startPosition);
237 
238  // The safety at the initial point has been re-calculated:
239  //
240  currentSafety = newSafety ;
241 
242  fGeometryLimitedStep= (linearStepLength <= currentMinimumStep);
244  {
245  // The geometry limits the Step size (an intersection was found.)
246  geometryStepLength = linearStepLength ;
247  }
248  else
249  {
250  // The full Step is taken.
251  geometryStepLength = currentMinimumStep ;
252  }
253  }
254  endpointDistance = geometryStepLength ;
255 
256  // Calculate final position
257  //
258  fTransportEndPosition = startPosition+geometryStepLength*startMomentumDir ;
259 
260  // Momentum direction, energy and polarisation are unchanged by transport
261  //
262  fTransportEndMomentumDir = startMomentumDir ;
265  fParticleIsLooping = false ;
266  fMomentumChanged = false ;
267  fEndGlobalTimeComputed = false ;
268  }
269  else // A field exerts force
270  {
271  G4double momentumMagnitude = pParticle->GetTotalMomentum() ;
272  G4ThreeVector EndUnitMomentum ;
273  G4double lengthAlongCurve ;
274  G4double restMass = fParticleDef->GetPDGMass() ;
275 
276  G4ChargeState chargeState(particleElectricCharge, // The charge can change (dynamic)
278  0, // Magnetic moment: pParticleDef->GetMagneticMoment(),
279  0, // Electric Dipole moment - not in Particle Definition
280  particleMagneticCharge ); // in Mev/c
281 
282  G4EquationOfMotion* equationOfMotion =
284  ->GetEquationOfMotion();
285 
286  equationOfMotion->SetChargeMomentumMass( chargeState, // Was particleMagneticCharge - in Mev/c
287  momentumMagnitude, // Was particleElectricCharge
288  restMass ) ;
289  // SetChargeMomentumMass now passes both the electric and magnetic charge - in chargeState
290 
291  G4ThreeVector spin = track.GetPolarization() ;
292  G4FieldTrack aFieldTrack = G4FieldTrack( startPosition,
293  track.GetMomentumDirection(),
294  0.0,
295  track.GetKineticEnergy(),
296  restMass,
297  track.GetVelocity(),
298  track.GetGlobalTime(), // Lab.
299  track.GetProperTime(), // Part.
300  &spin ) ;
301  if( currentMinimumStep > 0 )
302  {
303  // Do the Transport in the field (non recti-linear)
304  //
305  lengthAlongCurve = fFieldPropagator->ComputeStep( aFieldTrack,
306  currentMinimumStep,
307  currentSafety,
308  track.GetVolume() ) ;
309  fGeometryLimitedStep= lengthAlongCurve < currentMinimumStep;
310  if( fGeometryLimitedStep ) {
311  geometryStepLength = lengthAlongCurve ;
312  } else {
313  geometryStepLength = currentMinimumStep ;
314  }
315  }
316  else
317  {
318  geometryStepLength = lengthAlongCurve= 0.0 ;
319  fGeometryLimitedStep = false ;
320  }
321 
322  // Remember last safety origin & value.
323  //
324  fPreviousSftOrigin = startPosition ;
325  fPreviousSafety = currentSafety ;
326  // fpSafetyHelper->SetCurrentSafety( newSafety, startPosition);
327 
328  // Get the End-Position and End-Momentum (Dir-ection)
329  //
330  fTransportEndPosition = aFieldTrack.GetPosition() ;
331 
332  // Momentum: Magnitude and direction can be changed too now ...
333  //
334  fMomentumChanged = true ;
335  fTransportEndMomentumDir = aFieldTrack.GetMomentumDir() ;
336 
338 
340  fEndGlobalTimeComputed = true;
341 
342  fTransportEndSpin = aFieldTrack.GetSpin();
345  }
346 
347  // If we are asked to go a step length of 0, and we are on a boundary
348  // then a boundary will also limit the step -> we must flag this.
349  //
350  if( currentMinimumStep == 0.0 )
351  {
352  if( currentSafety == 0.0 ) fGeometryLimitedStep = true ;
353  }
354 
355  // Update the safety starting from the end-point,
356  // if it will become negative at the end-point.
357  //
358  if( currentSafety < endpointDistance )
359  {
360  // if( particleMagneticCharge == 0.0 )
361  // G4cout << " Avoiding call to ComputeSafety : charge = 0.0 " << G4endl;
362 
363  if( particleMagneticCharge != 0.0 ) {
364 
365  G4double endSafety =
367  currentSafety = endSafety ;
368  fPreviousSftOrigin = fTransportEndPosition ;
369  fPreviousSafety = currentSafety ;
371 
372  // Because the Stepping Manager assumes it is from the start point,
373  // add the StepLength
374  //
375  currentSafety += endpointDistance ;
376 
377 #ifdef G4DEBUG_TRANSPORT
378  G4cout.precision(12) ;
379  G4cout << "***G4MonopoleTransportation::AlongStepGPIL ** " << G4endl ;
380  G4cout << " Called Navigator->ComputeSafety at " << fTransportEndPosition
381  << " and it returned safety= " << endSafety << G4endl ;
382  G4cout << " Adding endpoint distance " << endpointDistance
383  << " to obtain pseudo-safety= " << currentSafety << G4endl ;
384 #endif
385  }
386  }
387 
388  fParticleChange.ProposeTrueStepLength(geometryStepLength) ;
389 
391  // change back to usual equation
392 
393  return geometryStepLength ;
394 }
395 
396 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
397 //
398 // Initialize ParticleChange (by setting all its members equal
399 // to corresponding members in G4Track)
400 
402  const G4Step& stepData )
403 {
404  static G4int noCalls=0;
405  static const G4ParticleDefinition* fOpticalPhoton =
407 
408  noCalls++;
409 
410  fParticleChange.Initialize(track) ;
411 
412  // Code for specific process
413  //
418 
420 
421  G4double deltaTime = 0.0 ;
422 
423  // Calculate Lab Time of Flight (ONLY if field Equations used it!)
424  // G4double endTime = fCandidateEndGlobalTime;
425  // G4double delta_time = endTime - startTime;
426 
427  G4double startTime = track.GetGlobalTime() ;
428 
430  {
431  // The time was not integrated .. make the best estimate possible
432  //
433  G4double finalVelocity = track.GetVelocity() ;
434  G4double initialVelocity = stepData.GetPreStepPoint()->GetVelocity() ;
435  G4double stepLength = track.GetStepLength() ;
436 
437  deltaTime= 0.0; // in case initialVelocity = 0
438  const G4DynamicParticle* fpDynamicParticle = track.GetDynamicParticle();
439  if (fpDynamicParticle->GetDefinition()== fOpticalPhoton)
440  {
441  // A photon is in the medium of the final point
442  // during the step, so it has the final velocity.
443  deltaTime = stepLength/finalVelocity ;
444  }
445  else if (finalVelocity > 0.0)
446  {
447  G4double meanInverseVelocity ;
448  // deltaTime = stepLength/finalVelocity ;
449  meanInverseVelocity = 0.5
450  * ( 1.0 / initialVelocity + 1.0 / finalVelocity ) ;
451  deltaTime = stepLength * meanInverseVelocity ;
452  }
453  else if( initialVelocity > 0.0 )
454  {
455  deltaTime = stepLength/initialVelocity ;
456  }
457  fCandidateEndGlobalTime = startTime + deltaTime ;
458  }
459  else
460  {
461  deltaTime = fCandidateEndGlobalTime - startTime ;
462  }
463 
465 
466  // Now Correct by Lorentz factor to get "proper" deltaTime
467 
468  G4double restMass = track.GetDynamicParticle()->GetMass() ;
469  G4double deltaProperTime = deltaTime*( restMass/track.GetTotalEnergy() ) ;
470 
471  fParticleChange.ProposeProperTime(track.GetProperTime() + deltaProperTime) ;
472  //fParticleChange. ProposeTrueStepLength( track.GetStepLength() ) ;
473 
474  // If the particle is caught looping or is stuck (in very difficult
475  // boundaries) in a magnetic field (doing many steps)
476  // THEN this kills it ...
477  //
478  if ( fParticleIsLooping )
479  {
481 
482  if( (endEnergy < fThreshold_Important_Energy)
484  // Kill the looping particle
485  //
487 
488  // 'Bare' statistics
489  fSumEnergyKilled += endEnergy;
490  if( endEnergy > fMaxEnergyKilled) { fMaxEnergyKilled= endEnergy; }
491 
492 #ifdef G4VERBOSE
493  if( (verboseLevel > 1) ||
494  ( endEnergy > fThreshold_Warning_Energy ) ) {
495  G4cout << " G4MonopoleTransportation is killing track that is looping or stuck "
496  << G4endl
497  << " This track has " << track.GetKineticEnergy() / MeV
498  << " MeV energy." << G4endl;
499  G4cout << " Number of trials = " << fNoLooperTrials
500  << " No of calls to AlongStepDoIt = " << noCalls
501  << G4endl;
502  }
503 #endif
504  fNoLooperTrials=0;
505  }
506  else{
507  fNoLooperTrials ++;
508 #ifdef G4VERBOSE
509  if( (verboseLevel > 2) ){
510  G4cout << " G4MonopoleTransportation::AlongStepDoIt(): Particle looping - "
511  << " Number of trials = " << fNoLooperTrials
512  << " No of calls to = " << noCalls
513  << G4endl;
514  }
515 #endif
516  }
517  }else{
518  fNoLooperTrials=0;
519  }
520 
521  // Another (sometimes better way) is to use a user-limit maximum Step size
522  // to alleviate this problem ..
523 
524  // Introduce smooth curved trajectories to particle-change
525  //
528 
529  return &fParticleChange ;
530 }
531 
532 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
533 //
534 // This ensures that the PostStep action is always called,
535 // so that it can do the relocation if it is needed.
536 //
537 
540  G4double, // previousStepSize
541  G4ForceCondition* pForceCond )
542 {
543  *pForceCond = Forced ;
544  return DBL_MAX ; // was kInfinity ; but convention now is DBL_MAX
545 }
546 
547 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
548 
550  const G4Step& )
551 {
552  G4TouchableHandle retCurrentTouchable ; // The one to return
553 
554  // Initialize ParticleChange (by setting all its members equal
555  // to corresponding members in G4Track)
556  // fParticleChange.Initialize(track) ; // To initialise TouchableChange
557 
559 
560  // If the Step was determined by the volume boundary,
561  // logically relocate the particle
562 
564  {
565  // fCurrentTouchable will now become the previous touchable,
566  // and what was the previous will be freed.
567  // (Needed because the preStepPoint can point to the previous touchable)
568 
571  LocateGlobalPointAndUpdateTouchableHandle( track.GetPosition(),
572  track.GetMomentumDirection(),
574  true ) ;
575  // Check whether the particle is out of the world volume
576  // If so it has exited and must be killed.
577  //
578  if( fCurrentTouchableHandle->GetVolume() == 0 )
579  {
581  }
582  retCurrentTouchable = fCurrentTouchableHandle ;
584  }
585  else // fGeometryLimitedStep is false
586  {
587  // This serves only to move the Navigator's location
588  //
590 
591  // The value of the track's current Touchable is retained.
592  // (and it must be correct because we must use it below to
593  // overwrite the (unset) one in particle change)
594  // It must be fCurrentTouchable too ??
595  //
597  retCurrentTouchable = track.GetTouchableHandle() ;
598  } // endif ( fGeometryLimitedStep )
599 
600  const G4VPhysicalVolume* pNewVol = retCurrentTouchable->GetVolume() ;
601  const G4Material* pNewMaterial = 0 ;
602  const G4VSensitiveDetector* pNewSensitiveDetector = 0 ;
603 
604  if( pNewVol != 0 )
605  {
606  pNewMaterial= pNewVol->GetLogicalVolume()->GetMaterial();
607  pNewSensitiveDetector= pNewVol->GetLogicalVolume()->GetSensitiveDetector();
608  }
609 
610  // ( <const_cast> pNewMaterial ) ;
611  // ( <const_cast> pNewSensitiveDetector) ;
612 
615 
616  const G4MaterialCutsCouple* pNewMaterialCutsCouple = 0;
617  if( pNewVol != 0 )
618  {
619  pNewMaterialCutsCouple=pNewVol->GetLogicalVolume()->GetMaterialCutsCouple();
620  }
621 
622  if( pNewVol!=0 && pNewMaterialCutsCouple!=0 && pNewMaterialCutsCouple->GetMaterial()!=pNewMaterial )
623  {
624  // for parametrized volume
625  //
626  pNewMaterialCutsCouple =
628  ->GetMaterialCutsCouple(pNewMaterial,
629  pNewMaterialCutsCouple->GetProductionCuts());
630  }
631  fParticleChange.SetMaterialCutsCoupleInTouchable( pNewMaterialCutsCouple );
632 
633  // temporarily until Get/Set Material of ParticleChange,
634  // and StepPoint can be made const.
635  // Set the touchable in ParticleChange
636  // this must always be done because the particle change always
637  // uses this value to overwrite the current touchable pointer.
638  //
639  fParticleChange.SetTouchableHandle(retCurrentTouchable) ;
640 
641  return &fParticleChange ;
642 }
643 
644 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
645 
646 // New method takes over the responsibility to reset the state
647 // of G4MonopoleTransportation object at the start of a new track
648 // or the resumption of a suspended track.
649 
650 void
652 {
654 
655 // The actions here are those that were taken in AlongStepGPIL
656 // when track.GetCurrentStepNumber()==1
657 
658  // reset safety value and center
659  //
660  fPreviousSafety = 0.0 ;
661  fPreviousSftOrigin = G4ThreeVector(0.,0.,0.) ;
662 
663  // reset looping counter -- for motion in field
664  fNoLooperTrials= 0;
665  // Must clear this state .. else it depends on last track's value
666  // --> a better solution would set this from state of suspended track TODO ?
667  // Was if( aTrack->GetCurrentStepNumber()==1 ) { .. }
668 
669  // ChordFinder reset internal state
670  //
671  if( DoesGlobalFieldExist() ) {
673  // Resets all state of field propagator class (ONLY)
674  // including safety values (in case of overlaps and to wipe for first track).
675 
676  // G4ChordFinder* chordF= fFieldPropagator->GetChordFinder();
677  // if( chordF ) chordF->ResetStepEstimate();
678  }
679 
680  // Make sure to clear the chord finders of all fields (ie managers)
681  static G4FieldManagerStore* fieldMgrStore= G4FieldManagerStore::GetInstance();
682  fieldMgrStore->ClearAllChordFindersState();
683 
684  // Update the current touchable handle (from the track's)
685  //
687 }
688 
689 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
690 
void SetMaterialInTouchable(G4Material *fMaterial)
void SetTouchableHandle(const G4TouchableHandle &fTouchable)
const G4ThreeVector & GetPolarization() const
static const double MeV
Definition: G4SIunits.hh:193
G4SafetyHelper * GetSafetyHelper() const
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
G4int verboseLevel
Definition: G4VProcess.hh:368
void SetMaterialCutsCoupleInTouchable(const G4MaterialCutsCouple *fMaterialCutsCouple)
virtual void SetChargeMomentumMass(G4ChargeState particleCharge, G4double MomentumXc, G4double MassXc2)=0
G4double GetProperTime() const
std::vector< G4ThreeVector > * GimmeTrajectoryVectorAndForgetIt() const
CLHEP::Hep3Vector G4ThreeVector
G4double GetVelocity() const
G4double GetKineticEnergy() const
const G4DynamicParticle * GetDynamicParticle() const
const G4MagIntegratorStepper * GetStepper() const
G4Material * GetMaterial() const
const G4ThreeVector & GetPosition() const
const G4ThreeVector & GetMomentumDir() const
G4TrackStatus GetTrackStatus() const
G4ThreeVector GetSpin() const
virtual G4double ComputeStep(const G4ThreeVector &pGlobalPoint, const G4ThreeVector &pDirection, const G4double pCurrentProposedStepLength, G4double &pNewSafety)
Definition: G4Navigator.cc:701
G4Navigator * GetNavigatorForTracking() const
virtual void StartTracking(G4Track *aTrack)
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData)
G4PropagatorInField * fFieldPropagator
void ProposePolarization(G4double Px, G4double Py, G4double Pz)
G4ParticleDefinition * GetDefinition() const
G4double GetVelocity() const
G4ParticleChangeForTransport fParticleChange
Definition of the G4MonopoleTransportation class.
int G4int
Definition: G4Types.hh:78
void ProposePosition(G4double x, G4double y, G4double z)
G4double GetTotalMomentum() const
G4StepPoint * GetPreStepPoint() const
virtual void Initialize(const G4Track &)
virtual void StartTracking(G4Track *)
Definition: G4VProcess.cc:101
virtual void ConfigureForTrack(const G4Track *)
Definition of the G4Monopole class.
void SetGeometricallyLimitedStep()
G4double GetKineticEnergy() const
G4ThreeVector GetPosition() const
G4GLOB_DLL std::ostream G4cout
G4double MagneticCharge() const
Definition: G4Monopole.cc:125
G4double GetMass() const
bool G4bool
Definition: G4Types.hh:79
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4ForceCondition *pForceCond)
const G4ThreeVector & GetMomentumDirection() const
G4double GetCharge() const
void ProposeTrueStepLength(G4double truePathLength)
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
Definition: G4Step.hh:76
static G4MonopoleFieldSetup * GetMonopoleFieldSetup()
G4double GetGlobalTime() const
G4FieldManager * FindAndSetFieldManager(G4VPhysicalVolume *pCurrentPhysVol)
const G4TouchableHandle & GetTouchableHandle() const
static G4TransportationManager * GetTransportationManager()
G4MonopoleTransportation(const G4Monopole *p, G4int verbosityLevel=1)
G4double ComputeStep(G4FieldTrack &pFieldTrack, G4double pCurrentProposedStepLength, G4double &pNewSafety, G4VPhysicalVolume *pPhysVol=0)
static G4ProductionCutsTable * GetProductionCutsTable()
const G4ThreeVector & GetMomentumDirection() const
void SetStepperAndChordFinder(G4int val)
G4LogicalVolume * GetLogicalVolume() const
void ProposeProperTime(G4double finalProperTime)
G4double GetPDGMass() const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
static G4FieldManagerStore * GetInstance()
void SetSensitiveDetectorInTouchable(G4VSensitiveDetector *fSensitiveDetector)
static G4ParticleTable * GetParticleTable()
G4bool IsParticleLooping() const
G4double GetLabTimeOfFlight() const
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:44
void SetPointerToVectorOfAuxiliaryPoints(std::vector< G4ThreeVector > *theNewVectorPointer)
void ProposeGlobalTime(G4double t)
G4ChordFinder * GetChordFinder()
virtual G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData)
void ProposeEnergy(G4double finalEnergy)
G4VPhysicalVolume * GetVolume() const
G4double GetTotalEnergy() const
G4double GetPDGSpin() const
void SetCurrentSafety(G4double val, const G4ThreeVector &pos)
#define G4endl
Definition: G4ios.hh:61
static char * startPosition
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
T sqr(const T &x)
Definition: templates.hh:145
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
const G4Field * GetDetectorField() const
virtual G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4double pProposedMaxLength=DBL_MAX, const G4bool keepState=true)
G4ForceCondition
G4MagInt_Driver * GetIntegrationDriver()
void SetMomentumChanged(G4bool b)
G4PropagatorInField * GetPropagatorInField() const
G4ProductionCuts * GetProductionCuts() const
#define DBL_MAX
Definition: templates.hh:83
G4VSensitiveDetector * GetSensitiveDetector() const
G4double GetStepLength() const
virtual void LocateGlobalPointWithinVolume(const G4ThreeVector &position)
Definition: G4Navigator.cc:550
G4MonopoleFieldSetup * fMagSetup
G4GPILSelection
const G4Material * GetMaterial() const