Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4MonopoleTransportation Class Reference

#include <G4MonopoleTransportation.hh>

Inheritance diagram for G4MonopoleTransportation:
Collaboration diagram for G4MonopoleTransportation:

Public Member Functions

 G4MonopoleTransportation (const G4Monopole *p, G4int verbosityLevel=1)
 
 ~G4MonopoleTransportation ()
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &, G4double previousStepSize, G4ForceCondition *pForceCond)
 
G4PropagatorInFieldGetPropagatorInField ()
 
void SetPropagatorInField (G4PropagatorInField *pFieldPropagator)
 
G4double GetThresholdWarningEnergy () const
 
G4double GetThresholdImportantEnergy () const
 
G4int GetThresholdTrials () const
 
void SetThresholdWarningEnergy (G4double newEnWarn)
 
void SetThresholdImportantEnergy (G4double newEnImp)
 
void SetThresholdTrials (G4int newMaxTrials)
 
G4double GetMaxEnergyKilled () const
 
G4double GetSumEnergyKilled () const
 
void ResetKilledStatistics (G4int report=1)
 
void EnableShortStepOptimisation (G4bool optimise=true)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual void StartTracking (G4Track *aTrack)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4int operator== (const G4VProcess &right) const
 
G4int operator!= (const G4VProcess &right) const
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Protected Member Functions

G4bool DoesGlobalFieldExist ()
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChangepParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 

Detailed Description

Definition at line 66 of file G4MonopoleTransportation.hh.

Constructor & Destructor Documentation

G4MonopoleTransportation::G4MonopoleTransportation ( const G4Monopole p,
G4int  verbosityLevel = 1 
)

Definition at line 60 of file G4MonopoleTransportation.cc.

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;
108  fCandidateEndGlobalTime = 0;
109 }
G4SafetyHelper * GetSafetyHelper() const
G4int verboseLevel
Definition: G4VProcess.hh:368
G4VProcess(const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
Definition: G4VProcess.cc:52
G4Navigator * GetNavigatorForTracking() const
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
static G4MonopoleFieldSetup * GetMonopoleFieldSetup()
static G4TransportationManager * GetTransportationManager()
static constexpr double MeV
Definition: G4SIunits.hh:214
G4PropagatorInField * GetPropagatorInField() const

Here is the call graph for this function:

G4MonopoleTransportation::~G4MonopoleTransportation ( )

Definition at line 113 of file G4MonopoleTransportation.cc.

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 }
G4int verboseLevel
Definition: G4VProcess.hh:368
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61

Member Function Documentation

G4VParticleChange * G4MonopoleTransportation::AlongStepDoIt ( const G4Track track,
const G4Step stepData 
)
virtual

Implements G4VProcess.

Definition at line 402 of file G4MonopoleTransportation.cc.

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  //
415  fParticleChange.ProposePosition(fTransportEndPosition) ;
416  fParticleChange.ProposeMomentumDirection(fTransportEndMomentumDir) ;
417  fParticleChange.ProposeEnergy(fTransportEndKineticEnergy) ;
418  fParticleChange.SetMomentumChanged(fMomentumChanged) ;
419 
420  fParticleChange.ProposePolarization(fTransportEndSpin);
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 
430  if (!fEndGlobalTimeComputed)
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 
465  fParticleChange.ProposeGlobalTime( fCandidateEndGlobalTime ) ;
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  {
481  G4double endEnergy= fTransportEndKineticEnergy;
482 
483  if( (endEnergy < fThreshold_Important_Energy)
484  || (fNoLooperTrials >= fThresholdTrials ) ){
485  // Kill the looping particle
486  //
487  fParticleChange.ProposeTrackStatus( fStopAndKill ) ;
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  //
528  (fFieldPropagator->GimmeTrajectoryVectorAndForgetIt() );
529 
530  return &fParticleChange ;
531 }
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
G4int verboseLevel
Definition: G4VProcess.hh:368
G4double GetProperTime() const
std::vector< G4ThreeVector > * GimmeTrajectoryVectorAndForgetIt() const
G4double GetVelocity() const
const G4DynamicParticle * GetDynamicParticle() const
void ProposePolarization(G4double Px, G4double Py, G4double Pz)
G4ParticleDefinition * GetDefinition() const
G4double GetVelocity() const
int G4int
Definition: G4Types.hh:78
void ProposePosition(G4double x, G4double y, G4double z)
G4StepPoint * GetPreStepPoint() const
virtual void Initialize(const G4Track &)
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
G4double GetMass() const
G4double GetGlobalTime() const
void ProposeProperTime(G4double finalProperTime)
static G4ParticleTable * GetParticleTable()
void SetPointerToVectorOfAuxiliaryPoints(std::vector< G4ThreeVector > *theNewVectorPointer)
void ProposeGlobalTime(G4double t)
void ProposeEnergy(G4double finalEnergy)
G4double GetTotalEnergy() const
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
void SetMomentumChanged(G4bool b)
G4double GetStepLength() const

Here is the call graph for this function:

G4double G4MonopoleTransportation::AlongStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4double  currentMinimumStep,
G4double currentSafety,
G4GPILSelection selection 
)
virtual

Implements G4VProcess.

Definition at line 131 of file G4MonopoleTransportation.cc.

136 {
137  fMagSetup->SetStepperAndChordFinder(1);
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);
243  if( fGeometryLimitedStep )
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 ;
263  fTransportEndKineticEnergy = track.GetKineticEnergy() ;
264  fTransportEndSpin = track.GetPolarization();
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)
277  fParticleDef->GetPDGSpin(),
278  0, // Magnetic moment: pParticleDef->GetMagneticMoment(),
279  0, // Electric Dipole moment - not in Particle Definition
280  particleMagneticCharge ); // in Mev/c
281 
282  G4EquationOfMotion* equationOfMotion =
283  (fFieldPropagator->GetChordFinder()->GetIntegrationDriver()->GetStepper())
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 
338  fTransportEndKineticEnergy = aFieldTrack.GetKineticEnergy() ;
339 
340  fCandidateEndGlobalTime = aFieldTrack.GetLabTimeOfFlight();
341  fEndGlobalTimeComputed = true;
342 
343  fTransportEndSpin = aFieldTrack.GetSpin();
344  fParticleIsLooping = fFieldPropagator->IsParticleLooping() ;
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 =
367  fLinearNavigator->ComputeSafety( fTransportEndPosition) ;
368  currentSafety = endSafety ;
369  fPreviousSftOrigin = fTransportEndPosition ;
370  fPreviousSafety = currentSafety ;
371  fpSafetyHelper->SetCurrentSafety( currentSafety, fTransportEndPosition);
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 
391  fMagSetup->SetStepperAndChordFinder(0);
392  // change back to usual equation
393 
394  return geometryStepLength ;
395 }
const G4ThreeVector & GetPolarization() const
virtual void SetChargeMomentumMass(G4ChargeState particleCharge, G4double MomentumXc, G4double MassXc2)=0
G4double GetProperTime() const
G4double GetVelocity() const
G4double GetKineticEnergy() const
const G4DynamicParticle * GetDynamicParticle() const
const G4MagIntegratorStepper * GetStepper() const
const G4ThreeVector & GetPosition() const
const G4ThreeVector & GetMomentumDir() const
G4ThreeVector GetSpin() const
virtual G4double ComputeStep(const G4ThreeVector &pGlobalPoint, const G4ThreeVector &pDirection, const G4double pCurrentProposedStepLength, G4double &pNewSafety)
Definition: G4Navigator.cc:747
G4double GetTotalMomentum() const
virtual void ConfigureForTrack(const G4Track *)
G4double GetKineticEnergy() const
G4ThreeVector GetPosition() const
G4GLOB_DLL std::ostream G4cout
G4double MagneticCharge() const
Definition: G4Monopole.cc:125
bool G4bool
Definition: G4Types.hh:79
const G4ThreeVector & GetMomentumDirection() const
G4double GetCharge() const
void ProposeTrueStepLength(G4double truePathLength)
G4double GetGlobalTime() const
G4FieldManager * FindAndSetFieldManager(G4VPhysicalVolume *pCurrentPhysVol)
G4double ComputeStep(G4FieldTrack &pFieldTrack, G4double pCurrentProposedStepLength, G4double &pNewSafety, G4VPhysicalVolume *pPhysVol=0)
const G4ThreeVector & GetMomentumDirection() const
void SetStepperAndChordFinder(G4int val)
G4double GetPDGMass() const
G4bool IsParticleLooping() const
G4double GetLabTimeOfFlight() const
G4ChordFinder * GetChordFinder()
double mag2() const
G4VPhysicalVolume * GetVolume() const
G4double GetPDGSpin() const
void SetCurrentSafety(G4double val, const G4ThreeVector &pos)
#define G4endl
Definition: G4ios.hh:61
T sqr(const T &x)
Definition: templates.hh:145
double G4double
Definition: G4Types.hh:76
const G4Field * GetDetectorField() const
virtual G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4double pProposedMaxLength=DBL_MAX, const G4bool keepState=true)
G4MagInt_Driver * GetIntegrationDriver()

Here is the call graph for this function:

virtual G4VParticleChange* G4MonopoleTransportation::AtRestDoIt ( const G4Track ,
const G4Step  
)
inlinevirtual

Implements G4VProcess.

Definition at line 135 of file G4MonopoleTransportation.hh.

138  {return 0;};
virtual G4double G4MonopoleTransportation::AtRestGetPhysicalInteractionLength ( const G4Track ,
G4ForceCondition  
)
inlinevirtual

Implements G4VProcess.

Definition at line 129 of file G4MonopoleTransportation.hh.

132  { return -1.0; };
G4bool G4MonopoleTransportation::DoesGlobalFieldExist ( )
protected

Here is the caller graph for this function:

void G4MonopoleTransportation::EnableShortStepOptimisation ( G4bool  optimise = true)
inline
G4double G4MonopoleTransportation::GetMaxEnergyKilled ( ) const
inline
G4PropagatorInField* G4MonopoleTransportation::GetPropagatorInField ( )
G4double G4MonopoleTransportation::GetSumEnergyKilled ( ) const
inline
G4double G4MonopoleTransportation::GetThresholdImportantEnergy ( ) const
inline
G4int G4MonopoleTransportation::GetThresholdTrials ( ) const
inline
G4double G4MonopoleTransportation::GetThresholdWarningEnergy ( ) const
inline
G4VParticleChange * G4MonopoleTransportation::PostStepDoIt ( const G4Track track,
const G4Step stepData 
)
virtual

Implements G4VProcess.

Definition at line 550 of file G4MonopoleTransportation.cc.

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 
559  fParticleChange.ProposeTrackStatus(track.GetTrackStatus()) ;
560 
561  // If the Step was determined by the volume boundary,
562  // logically relocate the particle
563 
564  if(fGeometryLimitedStep)
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 
570  fLinearNavigator->SetGeometricallyLimitedStep() ;
571  fLinearNavigator->
572  LocateGlobalPointAndUpdateTouchableHandle( track.GetPosition(),
573  track.GetMomentumDirection(),
574  fCurrentTouchableHandle,
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  {
581  fParticleChange.ProposeTrackStatus( fStopAndKill ) ;
582  }
583  retCurrentTouchable = fCurrentTouchableHandle ;
584  fParticleChange.SetTouchableHandle( fCurrentTouchableHandle ) ;
585  }
586  else // fGeometryLimitedStep is false
587  {
588  // This serves only to move the Navigator's location
589  //
590  fLinearNavigator->LocateGlobalPointWithinVolume( track.GetPosition() ) ;
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  //
597  fParticleChange.SetTouchableHandle( track.GetTouchableHandle() ) ;
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 
614  fParticleChange.SetMaterialInTouchable(
615  (G4Material *) pNewMaterial ) ;
616  fParticleChange.SetSensitiveDetectorInTouchable(
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 }
void SetMaterialInTouchable(G4Material *fMaterial)
void SetTouchableHandle(const G4TouchableHandle &fTouchable)
G4Material * GetMaterial() const
void SetMaterialCutsCoupleInTouchable(const G4MaterialCutsCouple *fMaterialCutsCouple)
const G4ThreeVector & GetPosition() const
G4TrackStatus GetTrackStatus() const
void SetGeometricallyLimitedStep()
const G4TouchableHandle & GetTouchableHandle() const
static G4ProductionCutsTable * GetProductionCutsTable()
const G4ThreeVector & GetMomentumDirection() const
G4LogicalVolume * GetLogicalVolume() const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
void SetSensitiveDetectorInTouchable(G4VSensitiveDetector *fSensitiveDetector)
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:44
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
void ProposeTrackStatus(G4TrackStatus status)
G4ProductionCuts * GetProductionCuts() const
G4VSensitiveDetector * GetSensitiveDetector() const
virtual void LocateGlobalPointWithinVolume(const G4ThreeVector &position)
Definition: G4Navigator.cc:582
const G4Material * GetMaterial() const

Here is the call graph for this function:

G4double G4MonopoleTransportation::PostStepGetPhysicalInteractionLength ( const G4Track ,
G4double  previousStepSize,
G4ForceCondition pForceCond 
)
virtual

Implements G4VProcess.

Definition at line 540 of file G4MonopoleTransportation.cc.

543 {
544  *pForceCond = Forced ;
545  return DBL_MAX ; // was kInfinity ; but convention now is DBL_MAX
546 }
#define DBL_MAX
Definition: templates.hh:83
void G4MonopoleTransportation::ResetKilledStatistics ( G4int  report = 1)
inline
void G4MonopoleTransportation::SetPropagatorInField ( G4PropagatorInField pFieldPropagator)
void G4MonopoleTransportation::SetThresholdImportantEnergy ( G4double  newEnImp)
inline
void G4MonopoleTransportation::SetThresholdTrials ( G4int  newMaxTrials)
inline
void G4MonopoleTransportation::SetThresholdWarningEnergy ( G4double  newEnWarn)
inline
void G4MonopoleTransportation::StartTracking ( G4Track aTrack)
virtual

Reimplemented from G4VProcess.

Definition at line 655 of file G4MonopoleTransportation.cc.

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() ) {
676  fFieldPropagator->ClearPropagatorState();
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  //
690  fCurrentTouchableHandle = aTrack->GetTouchableHandle();
691 }
CLHEP::Hep3Vector G4ThreeVector
virtual void StartTracking(G4Track *)
Definition: G4VProcess.cc:101
const G4TouchableHandle & GetTouchableHandle() const
static G4FieldManagerStore * GetInstance()

Here is the call graph for this function:


The documentation for this class was generated from the following files: