Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4FieldManager Class Reference

G4FieldManager is a manager (store) for a pointer to the Field subclass that describes the field of a detector (magnetic, electric or other). It also stores a reference to the chord finder. A field manager can be set to a logical volume (or to more than one), in order to vary its field from that of the world volume. In this manner a zero or constant field can override a global field, a more or less exact version can override the external approximation, lower or higher precision for tracking can be specified, a different stepper can be chosen for different volumes, etc... The Chord Finder must be created either by calling CreateChordFinder() for a Magnetic Field or by the user creating a Chord Finder object "manually" and setting the pointer. The current design envisions that one Field manager is valid for each detector region. It is expected that a particular geometrical region has a Field manager. By default a Field Manager is created for the world volume, and will be utilised for all volumes unless it is overridden by a 'local' field manager. Note also that a region with both electric E and magnetic B field will have these treated as one field. Similarly it could be extended to treat other fields as additional components of a single field type. More...

#include <G4FieldManager.hh>

Public Member Functions

 G4FieldManager (G4Field *detectorField=nullptr, G4ChordFinder *pChordFinder=nullptr, G4bool b=true)
 G4FieldManager (G4MagneticField *detectorMagneticField)
virtual ~G4FieldManager ()
 G4FieldManager (const G4FieldManager &)=delete
G4FieldManageroperator= (const G4FieldManager &)=delete
G4bool SetDetectorField (G4Field *detectorField, G4int failMode=0)
void ProposeDetectorField (G4Field *detectorField)
void ChangeDetectorField (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 *pTrack)
G4double GetDeltaIntersection () const
G4double GetDeltaOneStep () const
void SetAccuraciesWithDeltaOneStep (G4double valDeltaOneStep)
void SetDeltaOneStep (G4double valueD1step)
void SetDeltaIntersection (G4double valueDintersection)
G4double GetMinimumEpsilonStep () const
G4bool SetMinimumEpsilonStep (G4double newEpsMin)
G4double GetMaximumEpsilonStep () const
G4bool SetMaximumEpsilonStep (G4double newEpsMax)
G4bool DoesFieldChangeEnergy () const
void SetFieldChangesEnergy (G4bool value)
virtual G4FieldManagerClone () const

Static Public Member Functions

static void SetGlobalFieldManager (G4FieldManager *fieldManager)
static G4FieldManagerGetGlobalFieldManager ()
static G4double GetMaxAcceptedEpsilon ()
static G4bool SetMaxAcceptedEpsilon (G4double maxEps, G4bool softFail=false)

Protected Member Functions

void ReportBadEpsilonValue (G4ExceptionDescription &erm, G4double value, const G4String &name) const

Static Protected Attributes

static G4double fMaxAcceptedEpsilon = 0.01
static constexpr G4double fMinAcceptedEpsilon = 1000.0 * std::numeric_limits<G4double>::epsilon()
static constexpr G4double fMaxWarningEpsilon = 0.001
static constexpr G4double fMaxFinalEpsilon = 0.02
static G4bool fVerboseConstruction = false

Detailed Description

G4FieldManager is a manager (store) for a pointer to the Field subclass that describes the field of a detector (magnetic, electric or other). It also stores a reference to the chord finder. A field manager can be set to a logical volume (or to more than one), in order to vary its field from that of the world volume. In this manner a zero or constant field can override a global field, a more or less exact version can override the external approximation, lower or higher precision for tracking can be specified, a different stepper can be chosen for different volumes, etc... The Chord Finder must be created either by calling CreateChordFinder() for a Magnetic Field or by the user creating a Chord Finder object "manually" and setting the pointer. The current design envisions that one Field manager is valid for each detector region. It is expected that a particular geometrical region has a Field manager. By default a Field Manager is created for the world volume, and will be utilised for all volumes unless it is overridden by a 'local' field manager. Note also that a region with both electric E and magnetic B field will have these treated as one field. Similarly it could be extended to treat other fields as additional components of a single field type.

Definition at line 108 of file G4FieldManager.hh.

Constructor & Destructor Documentation

◆ G4FieldManager() [1/3]

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

General Constructor for any field. Must be set with field and chord finder for use.

Parameters
[in]detectorFieldPointer to the field.
[in]pChordFinderPointer to the chord finder object.
[in]bFlag to indicate if the field changes the energy; it is taken from the provided field, if specified.

Definition at line 50 of file G4FieldManager.cc.

53 : fDetectorField(detectorField),
54 fChordFinder(pChordFinder)
55{
56 if ( detectorField != nullptr )
57 {
58 fFieldChangesEnergy = detectorField->DoesFieldChangeEnergy();
59 }
60 else
61 {
62 fFieldChangesEnergy = fieldChangesEnergy;
63 }
64
66 {
67 G4cout << "G4FieldManager/ctor#1 fEpsilon Min/Max: eps_min = " << fEpsilonMin << " eps_max=" << fEpsilonMax << G4endl;
68 }
69
70 // Add to store
71 //
73}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
static void Register(G4FieldManager *pFieldMan)
static G4bool fVerboseConstruction
virtual G4bool DoesFieldChangeEnergy() const =0

Referenced by Clone(), G4FieldManager(), GetGlobalFieldManager(), operator=(), SetFieldChangesEnergy(), and SetGlobalFieldManager().

◆ G4FieldManager() [2/3]

G4FieldManager::G4FieldManager ( G4MagneticField * detectorMagneticField)

Constructor creating the chord finder. It assumes pure magnetic field, so energy constant.

Parameters
[in]detectorMagneticFieldPointer to the magnetic field.

Definition at line 75 of file G4FieldManager.cc.

76 : fDetectorField(detectorField), fAllocatedChordFinder(true)
77{
78 fChordFinder = new G4ChordFinder( detectorField );
79
81 {
82 G4cout << "G4FieldManager/ctor#2 fEpsilon Min/Max: eps_min = " << fEpsilonMin << " eps_max=" << fEpsilonMax << G4endl;
83 }
84 // Add to store
85 //
87}

◆ ~G4FieldManager()

G4FieldManager::~G4FieldManager ( )
virtual

Virtual Destructor.

Definition at line 89 of file G4FieldManager.cc.

90{
91 if( fAllocatedChordFinder )
92 {
93 delete fChordFinder;
94 }
96}
static void DeRegister(G4FieldManager *pFieldMan)

◆ G4FieldManager() [3/3]

G4FieldManager::G4FieldManager ( const G4FieldManager & )
delete

Copy constructor and assignment operator not allowed.

Member Function Documentation

◆ ChangeDetectorField()

void G4FieldManager::ChangeDetectorField ( G4Field * detectorField)
inline

Pushes the field to the equation and keeps its address. Can be used only once the equation, stepper, driver and chord finder have all been created; else it is an error.

Parameters
[in]detectorFieldPointer to the field.

◆ Clone()

G4FieldManager * G4FieldManager::Clone ( ) const
virtual

Needed for multi-threading, create and returns an allocated clone of this object.

Definition at line 98 of file G4FieldManager.cc.

99{
100 G4Field* aField = nullptr;
101 G4FieldManager* aFM = nullptr;
102 G4ChordFinder* aCF = nullptr;
103 try
104 {
105 if ( fDetectorField != nullptr )
106 {
107 aField = fDetectorField->Clone();
108 }
109
110 // Create a new field manager, note that we do not set
111 // any chordfinder now
112 //
113 aFM = new G4FieldManager( aField , nullptr , fFieldChangesEnergy );
114
115 // Check if originally we have the fAllocatedChordFinder variable
116 // set, in case, call chord constructor
117 //
118 if ( fAllocatedChordFinder )
119 {
120 aFM->CreateChordFinder( dynamic_cast<G4MagneticField*>(aField) );
121 }
122 else
123 {
124 // Chord was specified by user, should we clone?
125 // TODO: For the moment copy pointer, to be understood
126 // if cloning of ChordFinder is needed
127 //
128 aCF = fChordFinder; /*->Clone*/
129 aFM->fChordFinder = aCF;
130 }
131
132 // Copy values of other variables
133
134 aFM->fEpsilonMax = fEpsilonMax;
135 aFM->fEpsilonMin = fEpsilonMin;
136 aFM->fDelta_Intersection_Val = fDelta_Intersection_Val;
137 aFM->fDelta_One_Step_Value = fDelta_One_Step_Value;
138 // TODO: Should we really add to the store the cloned FM?
139 // Who will use this?
140 }
141 catch ( ... )
142 {
143 // Failed creating clone: probably user did not implement Clone method
144 // in derived classes?
145 // Perform clean-up after ourselves...
146 delete aField;
147 delete aFM;
148 delete aCF;
149 throw;
150 }
151
152 G4cout << "G4FieldManager/clone fEpsilon Min/Max: eps_min = " << fEpsilonMin << " eps_max=" << fEpsilonMax << G4endl;
153 return aFM;
154}
void CreateChordFinder(G4MagneticField *detectorMagField)
G4FieldManager(G4Field *detectorField=nullptr, G4ChordFinder *pChordFinder=nullptr, G4bool b=true)
virtual G4Field * Clone() const
Definition G4Field.cc:45

Referenced by G4VUserDetectorConstruction::CloneF().

◆ ConfigureForTrack()

void G4FieldManager::ConfigureForTrack ( const G4Track * pTrack)
virtual

Setups the choice of the configurable parameters, relying on the current track's energy, particle identity... Note: in addition to the values of member variables, a user can use this to change the ChordFinder, the field, etc.

Parameters
[in]pTrackPointer to a track.

Definition at line 156 of file G4FieldManager.cc.

157{
158 // Default is to do nothing!
159}

Referenced by G4CoupledTransportation::AlongStepGetPhysicalInteractionLength(), and G4ITTransportation::AlongStepGetPhysicalInteractionLength().

◆ CreateChordFinder()

void G4FieldManager::CreateChordFinder ( G4MagneticField * detectorMagField)

Methods to create, set or get the associated Chord Finder.

Definition at line 162 of file G4FieldManager.cc.

163{
164 if ( fAllocatedChordFinder )
165 {
166 delete fChordFinder;
167 }
168 fAllocatedChordFinder = false;
169
170 if( detectorMagField != nullptr )
171 {
172 fChordFinder = new G4ChordFinder( detectorMagField );
173 fAllocatedChordFinder = true;
174 }
175 else
176 {
177 fChordFinder = nullptr;
178 }
179}

Referenced by Clone().

◆ DoesFieldChangeEnergy()

G4bool G4FieldManager::DoesFieldChangeEnergy ( ) const
inline

Methods to set/get flag for field changing energy. For electric field this should be true; for magnetic field this should be false.

◆ DoesFieldExist()

G4bool G4FieldManager::DoesFieldExist ( ) const
inline

◆ GetChordFinder() [1/2]

G4ChordFinder * G4FieldManager::GetChordFinder ( )
inline

◆ GetChordFinder() [2/2]

const G4ChordFinder * G4FieldManager::GetChordFinder ( ) const
inline

◆ GetDeltaIntersection()

G4double G4FieldManager::GetDeltaIntersection ( ) const
inline

Returns the accuracy for boundary intersection.

◆ GetDeltaOneStep()

G4double G4FieldManager::GetDeltaOneStep ( ) const
inline

Returns the accuracy for one tracking/physics step.

◆ GetDetectorField()

◆ GetGlobalFieldManager()

G4FieldManager * G4FieldManager::GetGlobalFieldManager ( )
static

◆ GetMaxAcceptedEpsilon()

G4double G4FieldManager::GetMaxAcceptedEpsilon ( )
static

Static methods to set/get the maximum accepted epsilon. If setting fails, with softFail=true it gives Warning, else a FatalException.

Definition at line 336 of file G4FieldManager.cc.

337{
338 return fMaxAcceptedEpsilon;
339}
static G4double fMaxAcceptedEpsilon

◆ GetMaximumEpsilonStep()

G4double G4FieldManager::GetMaximumEpsilonStep ( ) const
inline

Methods to set/get the maximum for Relative accuracy of a Step.

◆ GetMinimumEpsilonStep()

G4double G4FieldManager::GetMinimumEpsilonStep ( ) const
inline

Methods to set/get the minimum for Relative accuracy of a Step.

◆ operator=()

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

◆ ProposeDetectorField()

void G4FieldManager::ProposeDetectorField ( G4Field * detectorField)
inline

Pushes the field to this class only – no further. Should be used to initialise this field, only before creating the chord finder and its dependent classes. User is then responsible to ensure that: i) an equation, stepper, driver and chord finder are created; ii) this field is used by the equation.

Parameters
[in]detectorFieldPointer to the field.

◆ ReportBadEpsilonValue()

void G4FieldManager::ReportBadEpsilonValue ( G4ExceptionDescription & erm,
G4double value,
const G4String & name ) const
protected

Logger for reporting on correctness of the proposed epsilon value.

Definition at line 415 of file G4FieldManager.cc.

417{
418 erm << "Incorrect proposed value of " << name << " = " << value << G4endl
419 << " Its value is outside the permitted range from "
421 << " Clarification: " << G4endl;
422 G4long oldPrec = erm.precision();
423 if(value < fMinAcceptedEpsilon )
424 {
425 erm << " a) The value must be positive and enough larger than the accuracy limit"
426 << " of the (G4)double type - ("
427 << (value < fMinAcceptedEpsilon ? "FAILED" : "OK" ) << ")" << G4endl
428 << " i.e. std::numeric_limits<G4double>::epsilon()= "
429 << std::numeric_limits<G4double>::epsilon()
430 << " to ensure that integration " << G4endl
431 << " could potentially achieve this acccuracy." << G4endl
432 << " Minimum accepted eps_min/max value = " << fMinAcceptedEpsilon << G4endl;
433 }
434 else if( value > fMaxAcceptedEpsilon)
435 {
436 erm << " b) It must be smaller than (or equal) " << std::setw(8)
437 << std::setprecision(4) << fMaxAcceptedEpsilon
438 << " to ensure robustness of integration - ("
439 << (( value < fMaxAcceptedEpsilon) ? "OK" : "FAILED" ) << ")" << G4endl;
440 }
441 else
442 {
443 G4bool badRoundoff = (std::fabs(1.0+value) == 1.0);
444 erm << " Unknown ERROR case -- extra check: " << G4endl;
445 erm << " c) as a floating point number (of type G4double) the sum (1+" << name
446 << " ) must be > 1 , ("
447 << (badRoundoff ? "FAILED" : "OK" ) << ")" << G4endl
448 << " Now 1+eps_min = " << std::setw(20)
449 << std::setprecision(17) << (1+value) << G4endl
450 << " and (1.0+" << name << ") - 1.0 = " << std::setw(20)
451 << std::setprecision(9) << (1.0+value)-1.0;
452 }
453 erm.precision(oldPrec);
454}
long G4long
Definition G4Types.hh:87
bool G4bool
Definition G4Types.hh:86
static constexpr G4double fMinAcceptedEpsilon
const char * name(G4int ptype)

Referenced by SetMaximumEpsilonStep(), and SetMinimumEpsilonStep().

◆ SetAccuraciesWithDeltaOneStep()

void G4FieldManager::SetAccuraciesWithDeltaOneStep ( G4double valDeltaOneStep)
inline

Sets both accuracies, maintaining a fixed ratio for accuracies of volume Intersection and Integration (in One Step).

◆ SetChordFinder()

void G4FieldManager::SetChordFinder ( G4ChordFinder * aChordFinder)
inline

◆ SetDeltaIntersection()

void G4FieldManager::SetDeltaIntersection ( G4double valueDintersection)
inline

Sets the accuracy of intersection of a volume (only).

◆ SetDeltaOneStep()

void G4FieldManager::SetDeltaOneStep ( G4double valueD1step)
inline

Sets the accuracy for integration of one step (only).

◆ SetDetectorField()

G4bool G4FieldManager::SetDetectorField ( G4Field * detectorField,
G4int failMode = 0 )

Pushes the field to the equation. Failure to push the field (due to absence of a chord finder, driver, stepper or equation) is

  • '0' = quiet : Do not complain if chordFinder == 0 (It will still warn for other error);
  • '1' = warn : a warning if anything is missing;
  • '2'/else = FATAL : a fatal error for all other values.
    Parameters
    [in]detectorFieldPointer to the field.
    [in]failModeFlag (0/1/2) for selected failure mode.
    Returns
    Success (true) or failure (false).

Definition at line 193 of file G4FieldManager.cc.

195{
196 G4VIntegrationDriver* driver = nullptr;
197 G4EquationOfMotion* equation = nullptr;
198 // G4bool compatibleField = false;
199 G4bool ableToSet = false;
200
201 fDetectorField = pDetectorField;
202 InitialiseFieldChangesEnergy();
203
204 // Must 'propagate' the field to the dependent classes
205 //
206 if( fChordFinder != nullptr )
207 {
208 failMode= std::max( failMode, 1) ;
209 // If a chord finder exists, warn in case of error!
210
211 driver = fChordFinder->GetIntegrationDriver();
212 if( driver != nullptr )
213 {
214 equation = driver->GetEquationOfMotion();
215
216 // Should check the compatibility between the
217 // field and the equation HERE
218
219 if( equation != nullptr )
220 {
221 equation->SetFieldObj(pDetectorField);
222 ableToSet = true;
223 }
224 }
225 }
226
227 if( !ableToSet && (failMode > 0) )
228 {
229 // If this fails, report the issue !
230
232 msg << "Unable to set the field in the dependent objects of G4FieldManager"
233 << G4endl;
234 msg << "All the dependent classes must be fully initialised,"
235 << "before it is possible to call this method." << G4endl;
236 msg << "The problem encountered was the following: " << G4endl;
237 if( fChordFinder == nullptr ) { msg << " No ChordFinder. " ; }
238 else if( driver == nullptr ) { msg << " No Integration Driver set. ";}
239 else if( equation == nullptr ) { msg << " No Equation found. " ; }
240 // else if( !compatibleField ) { msg << " Field not compatible. ";}
241 else { msg << " Can NOT find reason for failure. ";}
242 msg << G4endl;
243 G4ExceptionSeverity severity = (failMode != 1)
245 G4Exception("G4FieldManager::SetDetectorField", "Geometry001",
246 severity, msg);
247 }
248 return ableToSet;
249}
G4ExceptionSeverity
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
void SetFieldObj(G4Field *pField)
virtual G4EquationOfMotion * GetEquationOfMotion()=0

◆ SetFieldChangesEnergy()

void G4FieldManager::SetFieldChangesEnergy ( G4bool value)
inline

◆ SetGlobalFieldManager()

void G4FieldManager::SetGlobalFieldManager ( G4FieldManager * fieldManager)
static

Static methods to set/get the global field.

Referenced by G4TransportationManager::SetFieldManager().

◆ SetMaxAcceptedEpsilon()

G4bool G4FieldManager::SetMaxAcceptedEpsilon ( G4double maxEps,
G4bool softFail = false )
static

Definition at line 343 of file G4FieldManager.cc.

345{
346 G4bool success= false;
347 // Limit for warning and absolute limit chosen from experience in and
348 // investigation of integration with G4DormandPrince745 in HEP-type setups.
349 if( maxAcceptValue <= fMaxWarningEpsilon )
350 {
351 fMaxAcceptedEpsilon= maxAcceptValue;
352 success= true;
353 }
354 else
355 {
357 G4ExceptionSeverity severity;
358
359 G4cout << "G4FieldManager::" << __func__
360 << " Parameters: fMaxAcceptedEpsilon = " << fMaxAcceptedEpsilon
361 << " fMaxFinalEpsilon = " << fMaxFinalEpsilon << G4endl;
362
363 if( maxAcceptValue <= fMaxFinalEpsilon )
364 {
365 success= true;
366 fMaxAcceptedEpsilon = maxAcceptValue;
367 // Integration is poor, and robustness will likely suffer
368 erm << "Proposed value for maximum-accepted-epsilon = " << maxAcceptValue
369 << " is larger than the recommended = " << fMaxWarningEpsilon
370 << G4endl
371 << "This may impact the robustness of integration of tracks in field."
372 << G4endl
373 << "The request was accepted and the value = " << fMaxAcceptedEpsilon
374 << " , but future releases are expected " << G4endl
375 << " to tighten the limit of acceptable values to "
377 << "Suggestion: If you need better performance investigate using "
378 << "alternative, low-order RK integration methods or " << G4endl
379 << " helix-based methods (for pure B-fields) for low(er) energy tracks, "
380 << " especially electrons if you need better performance." << G4endl;
381 severity= JustWarning;
382 }
383 else
384 {
386 erm << " Proposed value for maximum accepted epsilon " << maxAcceptValue
387 << " is larger than the top of the range = " << fMaxFinalEpsilon
388 << G4endl;
389 if( softFailure )
390 {
391 erm << " Using the latter value instead." << G4endl;
392 }
393 erm << G4endl;
394 erm << " Please adjust to request maxAccepted <= " << fMaxFinalEpsilon
395 << G4endl << G4endl;
396 if( !softFailure )
397 {
398 erm << " NOTE: you can accept the ceiling value and turn this into a "
399 << " warning by using a 2nd argument " << G4endl
400 << " in your call to SetMaxAcceptedEpsilon: softFailure = true ";
401 }
402 severity = softFailure ? JustWarning : FatalException;
403 // if( softFailure ) severity= JustWarning;
404 // else severity= FatalException;
405 success = false;
406 }
407 G4String methodName = G4String("G4FieldManager::")+ G4String(__func__);
408 G4Exception(methodName.c_str(), "Geometry003", severity, erm);
409 }
410 return success;
411}
static constexpr G4double fMaxWarningEpsilon
static constexpr G4double fMaxFinalEpsilon

◆ SetMaximumEpsilonStep()

G4bool G4FieldManager::SetMaximumEpsilonStep ( G4double newEpsMax)

Definition at line 251 of file G4FieldManager.cc.

252{
253 G4bool succeeded= false;
254 if( (newEpsMax > 0.0) && ( newEpsMax <= fMaxAcceptedEpsilon)
255 && (fMinAcceptedEpsilon <= newEpsMax ) ) // (std::fabs(1.0+newEpsMax)>1.0) )
256 {
257 if(newEpsMax >= fEpsilonMin)
258 {
259 fEpsilonMax = newEpsMax;
260 succeeded = true;
262 {
263 G4cout << "G4FieldManager/SetEpsMax : eps_max = " << std::setw(10) << fEpsilonMax
264 << " ( Note: unchanged eps_min=" << std::setw(10) << fEpsilonMin << " )" << G4endl;
265 }
266 }
267 else
268 {
270 erm << " Call to set eps_max = " << newEpsMax << " . The problem is that"
271 << " its value must be at larger or equal to eps_min= " << fEpsilonMin << G4endl;
272 erm << " Modifying both to the same value " << newEpsMax << " to ensure consistency."
273 << G4endl
274 << " To avoid this warning, please set eps_min first, and ensure that "
275 << " 0 < eps_min <= eps_max <= " << fMaxAcceptedEpsilon << G4endl;
276
277 fEpsilonMax = newEpsMax;
278 fEpsilonMin = newEpsMax;
279 G4String methodName = G4String("G4FieldManager::")+ G4String(__func__);
280 G4Exception(methodName.c_str(), "Geometry003", JustWarning, erm);
281 }
282 }
283 else
284 {
286 G4String paramName("eps_max");
287 ReportBadEpsilonValue(erm, newEpsMax, paramName );
288 G4String methodName = G4String("G4FieldManager::")+ G4String(__func__);
289 G4Exception(methodName.c_str(), "Geometry001", FatalException, erm);
290 }
291 return succeeded;
292}
void ReportBadEpsilonValue(G4ExceptionDescription &erm, G4double value, const G4String &name) const

◆ SetMinimumEpsilonStep()

G4bool G4FieldManager::SetMinimumEpsilonStep ( G4double newEpsMin)

Definition at line 296 of file G4FieldManager.cc.

297{
298 G4bool succeeded= false;
299
300 if( fMinAcceptedEpsilon <= newEpsMin && newEpsMin <= fMaxAcceptedEpsilon )
301 {
302 fEpsilonMin = newEpsMin;
303 //*********
304 succeeded= true;
305
307 {
308 G4cout << "G4FieldManager/SetEpsMin : eps_min = "
309 << std::setw(10) << fEpsilonMin << G4endl;
310 }
311 if( fEpsilonMax < fEpsilonMin )
312 {
313 // Ensure consistency
315 erm << "Setting eps_min = " << newEpsMin
316 << " For consistency set eps_max= " << fEpsilonMin
317 << " ( Old value = " << fEpsilonMax << " )" << G4endl;
318 fEpsilonMax = fEpsilonMin;
319 G4String methodName = G4String("G4FieldManager::")+ G4String(__func__);
320 G4Exception(methodName.c_str(), "Geometry003", JustWarning, erm);
321 }
322 }
323 else
324 {
326 G4String paramName("eps_min");
327 ReportBadEpsilonValue(erm, newEpsMin, paramName );
328 G4String methodName = G4String("G4FieldManager::")+ G4String(__func__);
329 G4Exception(methodName.c_str(), "Geometry001", FatalException, erm);
330 }
331 return succeeded;
332}

Member Data Documentation

◆ fMaxAcceptedEpsilon

G4double G4FieldManager::fMaxAcceptedEpsilon = 0.01
staticprotected

Epsilon_min/max values must be smaller than this for robust integration.

Definition at line 272 of file G4FieldManager.hh.

Referenced by GetMaxAcceptedEpsilon(), ReportBadEpsilonValue(), SetMaxAcceptedEpsilon(), SetMaximumEpsilonStep(), and SetMinimumEpsilonStep().

◆ fMaxFinalEpsilon

G4double G4FieldManager::fMaxFinalEpsilon = 0.02
staticconstexprprotected

Will not accept larger values.

Definition at line 279 of file G4FieldManager.hh.

Referenced by SetMaxAcceptedEpsilon().

◆ fMaxWarningEpsilon

G4double G4FieldManager::fMaxWarningEpsilon = 0.001
staticconstexprprotected

Setting larger value will give warning.

Definition at line 276 of file G4FieldManager.hh.

Referenced by SetMaxAcceptedEpsilon().

◆ fMinAcceptedEpsilon

G4double G4FieldManager::fMinAcceptedEpsilon = 1000.0 * std::numeric_limits<G4double>::epsilon()
staticconstexprprotected

◆ fVerboseConstruction

G4bool G4FieldManager::fVerboseConstruction = false
staticprotected

Controls verbosity of constructors.

Definition at line 282 of file G4FieldManager.hh.

Referenced by G4FieldManager(), G4FieldManager(), SetMaximumEpsilonStep(), and SetMinimumEpsilonStep().


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