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

#include <G4FieldManager.hh>

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
 

Detailed Description

Definition at line 83 of file G4FieldManager.hh.

Constructor & Destructor Documentation

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),
46  fDefault_Delta_One_Step_Value(0.01), // mm
47  fDefault_Delta_Intersection_Val(0.001), // mm
48  fEpsilonMin( fEpsilonMinDefault ),
49  fEpsilonMax( fEpsilonMaxDefault)
50 {
51  fDelta_One_Step_Value= fDefault_Delta_One_Step_Value;
52  fDelta_Intersection_Val= fDefault_Delta_Intersection_Val;
53  if ( detectorField )
54  fFieldChangesEnergy= detectorField->DoesFieldChangeEnergy();
55  else
56  fFieldChangesEnergy= fieldChangesEnergy;
57 
58  // Add to store
60 }
virtual G4bool DoesFieldChangeEnergy() const =0
static void Register(G4FieldManager *pVolume)

Here is the call graph for this function:

Here is the caller graph for this function:

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),
67  fDefault_Delta_One_Step_Value(0.01), // mm
68  fDefault_Delta_Intersection_Val(0.001), // mm
69  fEpsilonMin( fEpsilonMinDefault ),
70  fEpsilonMax( fEpsilonMaxDefault)
71 {
72  fChordFinder= new G4ChordFinder( detectorField );
73  fDelta_One_Step_Value= fDefault_Delta_One_Step_Value;
74  fDelta_Intersection_Val= fDefault_Delta_Intersection_Val;
75  // Add to store
77 }
static void Register(G4FieldManager *pVolume)

Here is the call graph for this function:

G4FieldManager::~G4FieldManager ( )
virtual

Definition at line 132 of file G4FieldManager.cc.

133 {
134  if( fAllocatedChordFinder ){
135  delete fChordFinder;
136  }
138 }
static void DeRegister(G4FieldManager *pVolume)

Here is the call graph for this function:

Member Function Documentation

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;
107  aFM->fDefault_Delta_Intersection_Val = this->fDefault_Delta_Intersection_Val;
108  aFM->fDefault_Delta_One_Step_Value = this->fDefault_Delta_One_Step_Value;
109  aFM->fDelta_Intersection_Val = this->fDelta_Intersection_Val;
110  aFM->fDelta_One_Step_Value = this->fDelta_One_Step_Value;
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
G4FieldManager(G4Field *detectorField=0, G4ChordFinder *pChordFinder=0, G4bool b=true)
void CreateChordFinder(G4MagneticField *detectorMagField)

Here is the call graph for this function:

Here is the caller graph for this function:

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:

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 }

Here is the caller graph for this function:

G4bool G4FieldManager::DoesFieldChangeEnergy ( ) const
inline

Here is the caller graph for this function:

G4bool G4FieldManager::DoesFieldExist ( ) const
inline

Here is the caller graph for this function:

G4ChordFinder* G4FieldManager::GetChordFinder ( )
inline

Here is the caller graph for this function:

const G4ChordFinder* G4FieldManager::GetChordFinder ( ) const
inline
G4double G4FieldManager::GetDeltaIntersection ( ) const
inline

Here is the caller graph for this function:

G4double G4FieldManager::GetDeltaOneStep ( ) const
inline

Here is the caller graph for this function:

const G4Field* G4FieldManager::GetDetectorField ( ) const
inline

Here is the caller graph for this function:

G4double G4FieldManager::GetMaximumEpsilonStep ( ) const
inline

Here is the caller graph for this function:

G4double G4FieldManager::GetMinimumEpsilonStep ( ) const
inline

Here is the caller graph for this function:

void G4FieldManager::SetAccuraciesWithDeltaOneStep ( G4double  valDeltaOneStep)
inline
void G4FieldManager::SetChordFinder ( G4ChordFinder aChordFinder)
inline

Here is the caller graph for this function:

void G4FieldManager::SetDeltaIntersection ( G4double  valueDintersection)
inline
void G4FieldManager::SetDeltaOneStep ( G4double  valueD1step)
inline
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 }

Here is the call graph for this function:

void G4FieldManager::SetFieldChangesEnergy ( G4bool  value)
inline
void G4FieldManager::SetMaximumEpsilonStep ( G4double  newEpsMax)
inline
void G4FieldManager::SetMinimumEpsilonStep ( G4double  newEpsMin)
inline

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