Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4ErrorMagFieldLimitProcess Class Reference

#include <G4ErrorMagFieldLimitProcess.hh>

Inheritance diagram for G4ErrorMagFieldLimitProcess:
Collaboration diagram for G4ErrorMagFieldLimitProcess:

Public Member Functions

 G4ErrorMagFieldLimitProcess (const G4String &processName="G4ErrorMagFieldLimit")
 
 ~G4ErrorMagFieldLimitProcess ()
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
- Public Member Functions inherited from G4VErrorLimitProcess
 G4VErrorLimitProcess (const G4String &processName)
 
 ~G4VErrorLimitProcess ()
 
virtual G4double GetMeanFreePath (const class G4Track &, G4double, enum G4ForceCondition *)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
G4double GetStepLimit () const
 
void SetStepLimit (G4double val)
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
- 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 StartTracking (G4Track *)
 
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 &)
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Member Functions inherited from G4VDiscreteProcess
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- Protected Attributes inherited from G4VErrorLimitProcess
G4double theStepLimit
 
G4double theStepLength
 
G4VParticleChange theParticleChange
 
- 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 54 of file G4ErrorMagFieldLimitProcess.hh.

Constructor & Destructor Documentation

G4ErrorMagFieldLimitProcess::G4ErrorMagFieldLimitProcess ( const G4String processName = "G4ErrorMagFieldLimit")

Definition at line 46 of file G4ErrorMagFieldLimitProcess.cc.

47  : G4VErrorLimitProcess(processName)
48 {
50 }
static const G4double kInfinity
Definition: geomdefs.hh:42
G4VErrorLimitProcess(const G4String &processName)
G4ErrorMagFieldLimitProcess::~G4ErrorMagFieldLimitProcess ( )

Definition at line 54 of file G4ErrorMagFieldLimitProcess.cc.

55 { }

Member Function Documentation

G4double G4ErrorMagFieldLimitProcess::PostStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
virtual

Implements G4VErrorLimitProcess.

Definition at line 60 of file G4ErrorMagFieldLimitProcess.cc.

62 {
64  const G4Field* field =
66  ->GetDetectorField();
67 
69  if( field != 0 ) {
70  G4ThreeVector trkPosi = aTrack.GetPosition();
71  G4double pos1[3];
72  pos1[0] = trkPosi.x(); pos1[1] = trkPosi.y(); pos1[2] = trkPosi.z();
73 
74  G4double h1[3];
75  field->GetFieldValue( pos1, h1 );
76 
77  G4ThreeVector BVec(h1[0],h1[1],h1[2]);
78  G4double pmag = aTrack.GetMomentum().mag();
79  G4double BPerpMom = BVec.cross( G4ThreeVector(pmag,0.,0.) ).mag() / pmag;
80 
81  theStepLength = theStepLimit * pmag / BPerpMom;
82 #ifdef G4VERBOSE
83  if(G4ErrorPropagatorData::verbose() >= 3 ) {
84  G4cout << "G4ErrorMagFieldLimitProcess:: stepLength "
85  << theStepLength << " B " << BPerpMom << " BVec " << BVec
86  << " pmag " << pmag << G4endl;
87  }
88 #endif
89  }
90 
91  return theStepLength;
92 }
G4double condition(const G4ErrorSymMatrix &m)
static const G4double kInfinity
Definition: geomdefs.hh:42
CLHEP::Hep3Vector G4ThreeVector
virtual void GetFieldValue(const double Point[4], double *fieldArr) const =0
double x() const
double z() const
G4GLOB_DLL std::ostream G4cout
static G4TransportationManager * GetTransportationManager()
G4FieldManager * GetFieldManager() const
double y() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
const G4Field * GetDetectorField() const

Here is the call graph for this function:


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