Geant4  10.02
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 84606 2014-10-17 07:50:04Z 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() ;
158  G4ThreeVector startPosition = track.GetPosition() ;
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
287  ->SetChargeMomentumMass( chargeState, // Was particleMagneticCharge - in Mev/c
288  momentumMagnitude, // Was particleElectricCharge
289  restMass ) ;
290  // SetChargeMomentumMass now passes both the electric and magnetic charge - in chargeState
291 
292  G4ThreeVector spin = track.GetPolarization() ;
293  G4FieldTrack aFieldTrack = G4FieldTrack( startPosition,
294  track.GetMomentumDirection(),
295  0.0,
296  track.GetKineticEnergy(),
297  restMass,
298  track.GetVelocity(),
299  track.GetGlobalTime(), // Lab.
300  track.GetProperTime(), // Part.
301  &spin ) ;
302  if( currentMinimumStep > 0 )
303  {
304  // Do the Transport in the field (non recti-linear)
305  //
306  lengthAlongCurve = fFieldPropagator->ComputeStep( aFieldTrack,
307  currentMinimumStep,
308  currentSafety,
309  track.GetVolume() ) ;
310  fGeometryLimitedStep= lengthAlongCurve < currentMinimumStep;
311  if( fGeometryLimitedStep ) {
312  geometryStepLength = lengthAlongCurve ;
313  } else {
314  geometryStepLength = currentMinimumStep ;
315  }
316  }
317  else
318  {
319  geometryStepLength = lengthAlongCurve= 0.0 ;
320  fGeometryLimitedStep = false ;
321  }
322 
323  // Remember last safety origin & value.
324  //
325  fPreviousSftOrigin = startPosition ;
326  fPreviousSafety = currentSafety ;
327  // fpSafetyHelper->SetCurrentSafety( newSafety, startPosition);
328 
329  // Get the End-Position and End-Momentum (Dir-ection)
330  //
331  fTransportEndPosition = aFieldTrack.GetPosition() ;
332 
333  // Momentum: Magnitude and direction can be changed too now ...
334  //
335  fMomentumChanged = true ;
336  fTransportEndMomentumDir = aFieldTrack.GetMomentumDir() ;
337 
339 
341  fEndGlobalTimeComputed = true;
342 
343  fTransportEndSpin = aFieldTrack.GetSpin();
345  endpointDistance = (fTransportEndPosition - startPosition).mag() ;
346  }
347 
348  // If we are asked to go a step length of 0, and we are on a boundary
349  // then a boundary will also limit the step -> we must flag this.
350  //
351  if( currentMinimumStep == 0.0 )
352  {
353  if( currentSafety == 0.0 ) fGeometryLimitedStep = true ;
354  }
355 
356  // Update the safety starting from the end-point,
357  // if it will become negative at the end-point.
358  //
359  if( currentSafety < endpointDistance )
360  {
361  // if( particleMagneticCharge == 0.0 )
362  // G4cout << " Avoiding call to ComputeSafety : charge = 0.0 " << G4endl;
363 
364  if( particleMagneticCharge != 0.0 ) {
365 
366  G4double endSafety =
368  currentSafety = endSafety ;
369  fPreviousSftOrigin = fTransportEndPosition ;
370  fPreviousSafety = currentSafety ;
372 
373  // Because the Stepping Manager assumes it is from the start point,
374  // add the StepLength
375  //
376  currentSafety += endpointDistance ;
377 
378 #ifdef G4DEBUG_TRANSPORT
379  G4cout.precision(12) ;
380  G4cout << "***G4MonopoleTransportation::AlongStepGPIL ** " << G4endl ;
381  G4cout << " Called Navigator->ComputeSafety at " << fTransportEndPosition
382  << " and it returned safety= " << endSafety << G4endl ;
383  G4cout << " Adding endpoint distance " << endpointDistance
384  << " to obtain pseudo-safety= " << currentSafety << G4endl ;
385 #endif
386  }
387  }
388 
389  fParticleChange.ProposeTrueStepLength(geometryStepLength) ;
390 
392  // change back to usual equation
393 
394  return geometryStepLength ;
395 }
396 
397 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
398 //
399 // Initialize ParticleChange (by setting all its members equal
400 // to corresponding members in G4Track)
401 
403  const G4Step& stepData )
404 {
405  static G4int noCalls=0;
406  static const G4ParticleDefinition* fOpticalPhoton =
408 
409  noCalls++;
410 
411  fParticleChange.Initialize(track) ;
412 
413  // Code for specific process
414  //
419 
421 
422  G4double deltaTime = 0.0 ;
423 
424  // Calculate Lab Time of Flight (ONLY if field Equations used it!)
425  // G4double endTime = fCandidateEndGlobalTime;
426  // G4double delta_time = endTime - startTime;
427 
428  G4double startTime = track.GetGlobalTime() ;
429 
431  {
432  // The time was not integrated .. make the best estimate possible
433  //
434  G4double finalVelocity = track.GetVelocity() ;
435  G4double initialVelocity = stepData.GetPreStepPoint()->GetVelocity() ;
436  G4double stepLength = track.GetStepLength() ;
437 
438  deltaTime= 0.0; // in case initialVelocity = 0
439  const G4DynamicParticle* fpDynamicParticle = track.GetDynamicParticle();
440  if (fpDynamicParticle->GetDefinition()== fOpticalPhoton)
441  {
442  // A photon is in the medium of the final point
443  // during the step, so it has the final velocity.
444  deltaTime = stepLength/finalVelocity ;
445  }
446  else if (finalVelocity > 0.0)
447  {
448  G4double meanInverseVelocity ;
449  // deltaTime = stepLength/finalVelocity ;
450  meanInverseVelocity = 0.5
451  * ( 1.0 / initialVelocity + 1.0 / finalVelocity ) ;
452  deltaTime = stepLength * meanInverseVelocity ;
453  }
454  else if( initialVelocity > 0.0 )
455  {
456  deltaTime = stepLength/initialVelocity ;
457  }
458  fCandidateEndGlobalTime = startTime + deltaTime ;
459  }
460  else
461  {
462  deltaTime = fCandidateEndGlobalTime - startTime ;
463  }
464 
466 
467  // Now Correct by Lorentz factor to get "proper" deltaTime
468 
469  G4double restMass = track.GetDynamicParticle()->GetMass() ;
470  G4double deltaProperTime = deltaTime*( restMass/track.GetTotalEnergy() ) ;
471 
472  fParticleChange.ProposeProperTime(track.GetProperTime() + deltaProperTime) ;
473  //fParticleChange. ProposeTrueStepLength( track.GetStepLength() ) ;
474 
475  // If the particle is caught looping or is stuck (in very difficult
476  // boundaries) in a magnetic field (doing many steps)
477  // THEN this kills it ...
478  //
479  if ( fParticleIsLooping )
480  {
482 
483  if( (endEnergy < fThreshold_Important_Energy)
485  // Kill the looping particle
486  //
488 
489  // 'Bare' statistics
490  fSumEnergyKilled += endEnergy;
491  if( endEnergy > fMaxEnergyKilled) { fMaxEnergyKilled= endEnergy; }
492 
493 #ifdef G4VERBOSE
494  if( (verboseLevel > 1) ||
495  ( endEnergy > fThreshold_Warning_Energy ) ) {
496  G4cout << " G4MonopoleTransportation is killing track that is looping or stuck "
497  << G4endl
498  << " This track has " << track.GetKineticEnergy() / MeV
499  << " MeV energy." << G4endl;
500  G4cout << " Number of trials = " << fNoLooperTrials
501  << " No of calls to AlongStepDoIt = " << noCalls
502  << G4endl;
503  }
504 #endif
505  fNoLooperTrials=0;
506  }
507  else{
508  fNoLooperTrials ++;
509 #ifdef G4VERBOSE
510  if( (verboseLevel > 2) ){
511  G4cout << " G4MonopoleTransportation::AlongStepDoIt(): Particle looping - "
512  << " Number of trials = " << fNoLooperTrials
513  << " No of calls to = " << noCalls
514  << G4endl;
515  }
516 #endif
517  }
518  }else{
519  fNoLooperTrials=0;
520  }
521 
522  // Another (sometimes better way) is to use a user-limit maximum Step size
523  // to alleviate this problem ..
524 
525  // Introduce smooth curved trajectories to particle-change
526  //
529 
530  return &fParticleChange ;
531 }
532 
533 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
534 //
535 // This ensures that the PostStep action is always called,
536 // so that it can do the relocation if it is needed.
537 //
538 
541  G4double, // previousStepSize
542  G4ForceCondition* pForceCond )
543 {
544  *pForceCond = Forced ;
545  return DBL_MAX ; // was kInfinity ; but convention now is DBL_MAX
546 }
547 
548 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
549 
551  const G4Step& )
552 {
553  G4TouchableHandle retCurrentTouchable ; // The one to return
554 
555  // Initialize ParticleChange (by setting all its members equal
556  // to corresponding members in G4Track)
557  // fParticleChange.Initialize(track) ; // To initialise TouchableChange
558 
560 
561  // If the Step was determined by the volume boundary,
562  // logically relocate the particle
563 
565  {
566  // fCurrentTouchable will now become the previous touchable,
567  // and what was the previous will be freed.
568  // (Needed because the preStepPoint can point to the previous touchable)
569 
572  LocateGlobalPointAndUpdateTouchableHandle( track.GetPosition(),
573  track.GetMomentumDirection(),
575  true ) ;
576  // Check whether the particle is out of the world volume
577  // If so it has exited and must be killed.
578  //
579  if( fCurrentTouchableHandle->GetVolume() == 0 )
580  {
582  }
583  retCurrentTouchable = fCurrentTouchableHandle ;
585  }
586  else // fGeometryLimitedStep is false
587  {
588  // This serves only to move the Navigator's location
589  //
591 
592  // The value of the track's current Touchable is retained.
593  // (and it must be correct because we must use it below to
594  // overwrite the (unset) one in particle change)
595  // It must be fCurrentTouchable too ??
596  //
598  retCurrentTouchable = track.GetTouchableHandle() ;
599  } // endif ( fGeometryLimitedStep )
600 
601  const G4VPhysicalVolume* pNewVol = retCurrentTouchable->GetVolume() ;
602  const G4Material* pNewMaterial = 0 ;
603  const G4VSensitiveDetector* pNewSensitiveDetector = 0 ;
604 
605  if( pNewVol != 0 )
606  {
607  pNewMaterial= pNewVol->GetLogicalVolume()->GetMaterial();
608  pNewSensitiveDetector= pNewVol->GetLogicalVolume()->GetSensitiveDetector();
609  }
610 
611  // ( <const_cast> pNewMaterial ) ;
612  // ( <const_cast> pNewSensitiveDetector) ;
613 
615  (G4Material *) pNewMaterial ) ;
617  (G4VSensitiveDetector *) pNewSensitiveDetector ) ;
618 
619  const G4MaterialCutsCouple* pNewMaterialCutsCouple = 0;
620  if( pNewVol != 0 )
621  {
622  pNewMaterialCutsCouple=pNewVol->GetLogicalVolume()->GetMaterialCutsCouple();
623  }
624 
625  if( pNewVol!=0 && pNewMaterialCutsCouple!=0 &&
626  pNewMaterialCutsCouple->GetMaterial()!=pNewMaterial )
627  {
628  // for parametrized volume
629  //
630  pNewMaterialCutsCouple =
632  ->GetMaterialCutsCouple(pNewMaterial,
633  pNewMaterialCutsCouple->GetProductionCuts());
634  }
635  fParticleChange.SetMaterialCutsCoupleInTouchable( pNewMaterialCutsCouple );
636 
637  // temporarily until Get/Set Material of ParticleChange,
638  // and StepPoint can be made const.
639  // Set the touchable in ParticleChange
640  // this must always be done because the particle change always
641  // uses this value to overwrite the current touchable pointer.
642  //
643  fParticleChange.SetTouchableHandle(retCurrentTouchable) ;
644 
645  return &fParticleChange ;
646 }
647 
648 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
649 
650 // New method takes over the responsibility to reset the state
651 // of G4MonopoleTransportation object at the start of a new track
652 // or the resumption of a suspended track.
653 
654 void
656 {
658 
659 // The actions here are those that were taken in AlongStepGPIL
660 // when track.GetCurrentStepNumber()==1
661 
662  // reset safety value and center
663  //
664  fPreviousSafety = 0.0 ;
665  fPreviousSftOrigin = G4ThreeVector(0.,0.,0.) ;
666 
667  // reset looping counter -- for motion in field
668  fNoLooperTrials= 0;
669  // Must clear this state .. else it depends on last track's value
670  // --> a better solution would set this from state of suspended track TODO ?
671  // Was if( aTrack->GetCurrentStepNumber()==1 ) { .. }
672 
673  // ChordFinder reset internal state
674  //
675  if( DoesGlobalFieldExist() ) {
677  // Resets all state of field propagator class (ONLY)
678  // including safety values (in case of overlaps and to wipe for first track).
679 
680  // G4ChordFinder* chordF= fFieldPropagator->GetChordFinder();
681  // if( chordF ) chordF->ResetStepEstimate();
682  }
683 
684  // Make sure to clear the chord finders of all fields (ie managers)
685  static G4FieldManagerStore* fieldMgrStore= G4FieldManagerStore::GetInstance();
686  fieldMgrStore->ClearAllChordFindersState();
687 
688  // Update the current touchable handle (from the track's)
689  //
691 }
692 
693 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
694 
void SetMaterialInTouchable(G4Material *fMaterial)
void SetTouchableHandle(const G4TouchableHandle &fTouchable)
const G4ThreeVector & GetPolarization() const
static const double MeV
Definition: G4SIunits.hh:211
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:744
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()
#define fPreviousSftOrigin
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()
#define fPreviousSafety
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
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:579
G4MonopoleFieldSetup * fMagSetup
G4GPILSelection
const G4Material * GetMaterial() const