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