Geant4  10.02.p03
G4FieldManager Class Reference

#include <G4FieldManager.hh>

Collaboration diagram for G4FieldManager:

Public Member Functions

 G4FieldManager (G4Field *detectorField=0, G4ChordFinder *pChordFinder=0, G4bool b=true)
 
 G4FieldManager (G4MagneticField *detectorMagneticField)
 
virtual ~G4FieldManager ()
 
G4bool SetDetectorField (G4Field *detectorField)
 
const G4FieldGetDetectorField () const
 
G4bool DoesFieldExist () const
 
void CreateChordFinder (G4MagneticField *detectorMagField)
 
void SetChordFinder (G4ChordFinder *aChordFinder)
 
G4ChordFinderGetChordFinder ()
 
const G4ChordFinderGetChordFinder () const
 
virtual void ConfigureForTrack (const G4Track *)
 
G4double GetDeltaIntersection () const
 
G4double GetDeltaOneStep () const
 
void SetAccuraciesWithDeltaOneStep (G4double valDeltaOneStep)
 
void SetDeltaOneStep (G4double valueD1step)
 
void SetDeltaIntersection (G4double valueDintersection)
 
G4double GetMinimumEpsilonStep () const
 
void SetMinimumEpsilonStep (G4double newEpsMin)
 
G4double GetMaximumEpsilonStep () const
 
void SetMaximumEpsilonStep (G4double newEpsMax)
 
G4bool DoesFieldChangeEnergy () const
 
void SetFieldChangesEnergy (G4bool value)
 
virtual G4FieldManagerClone () const
 

Private Member Functions

 G4FieldManager (const G4FieldManager &)
 
G4FieldManageroperator= (const G4FieldManager &)
 

Private Attributes

G4FieldfDetectorField
 
G4ChordFinderfChordFinder
 
G4bool fAllocatedChordFinder
 
const G4double fEpsilonMinDefault
 
const G4double fEpsilonMaxDefault
 
G4bool fFieldChangesEnergy
 
G4double fDelta_One_Step_Value
 
G4double fDelta_Intersection_Val
 
G4double fDefault_Delta_One_Step_Value
 
G4double fDefault_Delta_Intersection_Val
 
G4double fEpsilonMin
 
G4double fEpsilonMax
 

Detailed Description

Definition at line 83 of file G4FieldManager.hh.

Constructor & Destructor Documentation

◆ G4FieldManager() [1/3]

G4FieldManager::G4FieldManager ( G4Field detectorField = 0,
G4ChordFinder pChordFinder = 0,
G4bool  b = true 
)

Definition at line 37 of file G4FieldManager.cc.

41  : fDetectorField(detectorField),
42  fChordFinder(pChordFinder),
43  fAllocatedChordFinder(false),
44  fEpsilonMinDefault(5.0e-5),
45  fEpsilonMaxDefault(0.001),
50 {
53  if ( detectorField )
55  else
56  fFieldChangesEnergy= fieldChangesEnergy;
57 
58  // Add to store
60 }
G4double fDelta_One_Step_Value
G4double fEpsilonMin
G4double fDefault_Delta_Intersection_Val
virtual G4bool DoesFieldChangeEnergy() const =0
G4double fDefault_Delta_One_Step_Value
const G4double fEpsilonMaxDefault
G4double fEpsilonMax
G4ChordFinder * fChordFinder
G4bool fFieldChangesEnergy
static void Register(G4FieldManager *pVolume)
G4Field * fDetectorField
G4double fDelta_Intersection_Val
G4bool fAllocatedChordFinder
const G4double fEpsilonMinDefault
Here is the call graph for this function:
Here is the caller graph for this function:

◆ G4FieldManager() [2/3]

G4FieldManager::G4FieldManager ( G4MagneticField detectorMagneticField)

Definition at line 62 of file G4FieldManager.cc.

63  : fDetectorField(detectorField), fAllocatedChordFinder(true),
64  fEpsilonMinDefault(5.0e-5),
65  fEpsilonMaxDefault(0.001),
66  fFieldChangesEnergy(false),
71 {
72  fChordFinder= new G4ChordFinder( detectorField );
75  // Add to store
77 }
G4double fDelta_One_Step_Value
G4double fEpsilonMin
G4double fDefault_Delta_Intersection_Val
G4double fDefault_Delta_One_Step_Value
const G4double fEpsilonMaxDefault
G4double fEpsilonMax
G4ChordFinder * fChordFinder
G4bool fFieldChangesEnergy
static void Register(G4FieldManager *pVolume)
G4Field * fDetectorField
G4double fDelta_Intersection_Val
G4bool fAllocatedChordFinder
const G4double fEpsilonMinDefault
Here is the call graph for this function:

◆ ~G4FieldManager()

G4FieldManager::~G4FieldManager ( )
virtual

Definition at line 132 of file G4FieldManager.cc.

133 {
134  if( fAllocatedChordFinder ){
135  delete fChordFinder;
136  }
138 }
G4ChordFinder * fChordFinder
G4bool fAllocatedChordFinder
static void DeRegister(G4FieldManager *pVolume)
Here is the call graph for this function:

◆ G4FieldManager() [3/3]

G4FieldManager::G4FieldManager ( const G4FieldManager )
private

Member Function Documentation

◆ Clone()

G4FieldManager * G4FieldManager::Clone ( ) const
virtual

Definition at line 79 of file G4FieldManager.cc.

80 {
81  G4Field* aField = 0;
82  G4FieldManager* aFM = 0;
83  G4ChordFinder* aCF = 0;
84  try {
85  if ( this->fDetectorField )
86  aField = this->fDetectorField->Clone();
87 
88  //Create a new field manager, note that we do not set any chordfinder now.
89  aFM = new G4FieldManager( aField , 0 , this->fFieldChangesEnergy );
90 
91  //Check if orignally we have the fAllocatedChordFinder variable set, in case, call chord
92  //constructor
93  if ( this->fAllocatedChordFinder )
94  {
95  aFM->CreateChordFinder( dynamic_cast<G4MagneticField*>(aField) );
96  }
97  else
98  {
99  //Chord was specified by user, should we clone?
100  //TODO: For the moment copy pointer, to be understood if cloning of ChordFinder is needed
101  aCF = this->fChordFinder;/*->Clone*/
102  aFM->fChordFinder = aCF;
103  }
104  //Copy values of other variables
105  aFM->fEpsilonMax = this->fEpsilonMax;
106  aFM->fEpsilonMin = this->fEpsilonMin;
111  //TODO: Should we really add to the store the cloned FM? Who will use this?
112  }
113  catch ( ... )
114  {
115  //Failed creating clone: probably user did not implement Clone method
116  //in derived classes?
117  //Perform clean-up after ourselves...
118  delete aField;
119  delete aFM;
120  delete aCF;
121  throw;
122  }
123  return aFM;
124 }
virtual G4Field * Clone() const
Definition: G4Field.cc:54
G4double fDelta_One_Step_Value
G4double fEpsilonMin
G4double fDefault_Delta_Intersection_Val
G4double fDefault_Delta_One_Step_Value
G4FieldManager(G4Field *detectorField=0, G4ChordFinder *pChordFinder=0, G4bool b=true)
G4double fEpsilonMax
G4ChordFinder * fChordFinder
G4bool fFieldChangesEnergy
G4Field * fDetectorField
void CreateChordFinder(G4MagneticField *detectorMagField)
G4double fDelta_Intersection_Val
G4bool fAllocatedChordFinder
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ConfigureForTrack()

void G4FieldManager::ConfigureForTrack ( const G4Track *  )
virtual

Definition at line 126 of file G4FieldManager.cc.

127 {
128  // Default is to do nothing!
129  ;
130 }
Here is the caller graph for this function:

◆ CreateChordFinder()

void G4FieldManager::CreateChordFinder ( G4MagneticField detectorMagField)

Definition at line 141 of file G4FieldManager.cc.

142 {
143  if ( fAllocatedChordFinder )
144  delete fChordFinder;
145  fChordFinder= new G4ChordFinder( detectorMagField );
146  fAllocatedChordFinder= true;
147 }
G4ChordFinder * fChordFinder
G4bool fAllocatedChordFinder
Here is the caller graph for this function:

◆ DoesFieldChangeEnergy()

G4bool G4FieldManager::DoesFieldChangeEnergy ( ) const
inline
Here is the caller graph for this function:

◆ DoesFieldExist()

G4bool G4FieldManager::DoesFieldExist ( ) const
inline
Here is the caller graph for this function:

◆ GetChordFinder() [1/2]

G4ChordFinder* G4FieldManager::GetChordFinder ( )
inline
Here is the caller graph for this function:

◆ GetChordFinder() [2/2]

const G4ChordFinder * G4FieldManager::GetChordFinder ( ) const
inline

Definition at line 47 of file pyG4FieldManager.cc.

48 {
49  return fChordFinder;
50 }
G4ChordFinder * fChordFinder

◆ GetDeltaIntersection()

G4double G4FieldManager::GetDeltaIntersection ( ) const
inline
Here is the caller graph for this function:

◆ GetDeltaOneStep()

G4double G4FieldManager::GetDeltaOneStep ( ) const
inline
Here is the caller graph for this function:

◆ GetDetectorField()

const G4Field* G4FieldManager::GetDetectorField ( ) const
inline
Here is the caller graph for this function:

◆ GetMaximumEpsilonStep()

G4double G4FieldManager::GetMaximumEpsilonStep ( ) const
inline
Here is the caller graph for this function:

◆ GetMinimumEpsilonStep()

G4double G4FieldManager::GetMinimumEpsilonStep ( ) const
inline
Here is the caller graph for this function:

◆ operator=()

G4FieldManager& G4FieldManager::operator= ( const G4FieldManager )
private

◆ SetAccuraciesWithDeltaOneStep()

void G4FieldManager::SetAccuraciesWithDeltaOneStep ( G4double  valDeltaOneStep)
inline
Here is the caller graph for this function:

◆ SetChordFinder()

void G4FieldManager::SetChordFinder ( G4ChordFinder aChordFinder)
inline
Here is the caller graph for this function:

◆ SetDeltaIntersection()

void G4FieldManager::SetDeltaIntersection ( G4double  valueDintersection)
inline
Here is the caller graph for this function:

◆ SetDeltaOneStep()

void G4FieldManager::SetDeltaOneStep ( G4double  valueD1step)
inline
Here is the caller graph for this function:

◆ SetDetectorField()

G4bool G4FieldManager::SetDetectorField ( G4Field detectorField)

Definition at line 149 of file G4FieldManager.cc.

150 {
151  fDetectorField= pDetectorField;
152 
153  if ( pDetectorField )
154  fFieldChangesEnergy= pDetectorField->DoesFieldChangeEnergy();
155  else
156  fFieldChangesEnergy= false; // No field
157 
158  return false;
159 }
G4bool fFieldChangesEnergy
G4Field * fDetectorField
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetFieldChangesEnergy()

void G4FieldManager::SetFieldChangesEnergy ( G4bool  value)
inline
Here is the caller graph for this function:

◆ SetMaximumEpsilonStep()

void G4FieldManager::SetMaximumEpsilonStep ( G4double  newEpsMax)
inline
Here is the caller graph for this function:

◆ SetMinimumEpsilonStep()

void G4FieldManager::SetMinimumEpsilonStep ( G4double  newEpsMin)
inline
Here is the caller graph for this function:

Member Data Documentation

◆ fAllocatedChordFinder

G4bool G4FieldManager::fAllocatedChordFinder
private

Definition at line 157 of file G4FieldManager.hh.

◆ fChordFinder

G4ChordFinder* G4FieldManager::fChordFinder
private

Definition at line 155 of file G4FieldManager.hh.

◆ fDefault_Delta_Intersection_Val

G4double G4FieldManager::fDefault_Delta_Intersection_Val
private

Definition at line 175 of file G4FieldManager.hh.

◆ fDefault_Delta_One_Step_Value

G4double G4FieldManager::fDefault_Delta_One_Step_Value
private

Definition at line 174 of file G4FieldManager.hh.

◆ fDelta_Intersection_Val

G4double G4FieldManager::fDelta_Intersection_Val
private

Definition at line 172 of file G4FieldManager.hh.

◆ fDelta_One_Step_Value

G4double G4FieldManager::fDelta_One_Step_Value
private

Definition at line 171 of file G4FieldManager.hh.

◆ fDetectorField

G4Field* G4FieldManager::fDetectorField
private

Definition at line 154 of file G4FieldManager.hh.

◆ fEpsilonMax

G4double G4FieldManager::fEpsilonMax
private

Definition at line 180 of file G4FieldManager.hh.

◆ fEpsilonMaxDefault

const G4double G4FieldManager::fEpsilonMaxDefault
private

Definition at line 163 of file G4FieldManager.hh.

◆ fEpsilonMin

G4double G4FieldManager::fEpsilonMin
private

Definition at line 179 of file G4FieldManager.hh.

◆ fEpsilonMinDefault

const G4double G4FieldManager::fEpsilonMinDefault
private

Definition at line 162 of file G4FieldManager.hh.

◆ fFieldChangesEnergy

G4bool G4FieldManager::fFieldChangesEnergy
private

Definition at line 166 of file G4FieldManager.hh.


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