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

#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 *)
 
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 G4double GetMaxAcceptedEpsilon ()
 
static G4bool SetMaxAcceptedEpsilon (G4double maxEps, G4bool softFail=false)
 

Protected Member Functions

void ReportBadEpsilonValue (G4ExceptionDescription &erm, G4double value, 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

Definition at line 84 of file G4FieldManager.hh.

Constructor & Destructor Documentation

◆ G4FieldManager() [1/3]

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

Definition at line 49 of file G4FieldManager.cc.

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

◆ G4FieldManager() [2/3]

G4FieldManager::G4FieldManager ( G4MagneticField detectorMagneticField)

Definition at line 76 of file G4FieldManager.cc.

77 : fDetectorField(detectorField), fAllocatedChordFinder(true),
78 fDelta_One_Step_Value( fDefault_Delta_One_Step_Value ),
79 fDelta_Intersection_Val( fDefault_Delta_Intersection_Val ),
80 fEpsilonMin( fEpsilonMinDefault ),
81 fEpsilonMax( fEpsilonMaxDefault )
82{
83 fChordFinder = new G4ChordFinder( detectorField );
84
86 G4cout << "G4FieldManager/ctor#2 fEpsilon Min/Max: eps_min = " << fEpsilonMin << " eps_max=" << fEpsilonMax << G4endl;
87 // Add to store
88 //
90}

◆ ~G4FieldManager()

G4FieldManager::~G4FieldManager ( )
virtual

Definition at line 154 of file G4FieldManager.cc.

155{
156 if( fAllocatedChordFinder )
157 {
158 delete fChordFinder;
159 }
161}
static void DeRegister(G4FieldManager *pVolume)

◆ G4FieldManager() [3/3]

G4FieldManager::G4FieldManager ( const G4FieldManager )
delete

Member Function Documentation

◆ ChangeDetectorField()

void G4FieldManager::ChangeDetectorField ( G4Field detectorField)
inline

◆ Clone()

G4FieldManager * G4FieldManager::Clone ( ) const
virtual

Definition at line 92 of file G4FieldManager.cc.

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

Referenced by G4VUserDetectorConstruction::CloneF().

◆ ConfigureForTrack()

void G4FieldManager::ConfigureForTrack ( const G4Track )
virtual

Definition at line 149 of file G4FieldManager.cc.

150{
151 // Default is to do nothing!
152}

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

◆ CreateChordFinder()

void G4FieldManager::CreateChordFinder ( G4MagneticField detectorMagField)

Definition at line 164 of file G4FieldManager.cc.

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

Referenced by Clone().

◆ DoesFieldChangeEnergy()

◆ 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

◆ GetDeltaOneStep()

G4double G4FieldManager::GetDeltaOneStep ( ) const
inline

◆ GetDetectorField()

◆ GetMaxAcceptedEpsilon()

G4double G4FieldManager::GetMaxAcceptedEpsilon ( )
static

Definition at line 334 of file G4FieldManager.cc.

335{
336 return fMaxAcceptedEpsilon;
337}
static G4double fMaxAcceptedEpsilon

◆ GetMaximumEpsilonStep()

G4double G4FieldManager::GetMaximumEpsilonStep ( ) const
inline

◆ GetMinimumEpsilonStep()

G4double G4FieldManager::GetMinimumEpsilonStep ( ) const
inline

◆ operator=()

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

◆ ProposeDetectorField()

void G4FieldManager::ProposeDetectorField ( G4Field detectorField)
inline

◆ ReportBadEpsilonValue()

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

Definition at line 409 of file G4FieldManager.cc.

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

◆ SetChordFinder()

void G4FieldManager::SetChordFinder ( G4ChordFinder aChordFinder)
inline

◆ SetDeltaIntersection()

void G4FieldManager::SetDeltaIntersection ( G4double  valueDintersection)
inline

◆ SetDeltaOneStep()

void G4FieldManager::SetDeltaOneStep ( G4double  valueD1step)
inline

◆ SetDetectorField()

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

Definition at line 195 of file G4FieldManager.cc.

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

◆ SetFieldChangesEnergy()

void G4FieldManager::SetFieldChangesEnergy ( G4bool  value)
inline

◆ SetMaxAcceptedEpsilon()

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

Definition at line 341 of file G4FieldManager.cc.

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

◆ SetMaximumEpsilonStep()

G4bool G4FieldManager::SetMaximumEpsilonStep ( G4double  newEpsMax)

Definition at line 253 of file G4FieldManager.cc.

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

◆ SetMinimumEpsilonStep()

G4bool G4FieldManager::SetMinimumEpsilonStep ( G4double  newEpsMin)

Definition at line 295 of file G4FieldManager.cc.

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

Member Data Documentation

◆ fMaxAcceptedEpsilon

G4double G4FieldManager::fMaxAcceptedEpsilon = 0.01
staticprotected

◆ fMaxFinalEpsilon

constexpr G4double G4FieldManager::fMaxFinalEpsilon = 0.02
staticconstexprprotected

Definition at line 186 of file G4FieldManager.hh.

Referenced by SetMaxAcceptedEpsilon().

◆ fMaxWarningEpsilon

constexpr G4double G4FieldManager::fMaxWarningEpsilon = 0.001
staticconstexprprotected

Definition at line 185 of file G4FieldManager.hh.

Referenced by SetMaxAcceptedEpsilon().

◆ fMinAcceptedEpsilon

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

◆ fVerboseConstruction

G4bool G4FieldManager::fVerboseConstruction = false
staticprotected

Definition at line 188 of file G4FieldManager.hh.

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


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