39G4double G4FieldManager::fDefault_Delta_One_Step_Value= 0.01 * millimeter;
40G4double G4FieldManager::fDefault_Delta_Intersection_Val= 0.001 * millimeter;
52 G4bool fieldChangesEnergy )
53 : fDetectorField(detectorField),
54 fChordFinder(pChordFinder),
55 fDelta_One_Step_Value( fDefault_Delta_One_Step_Value ),
56 fDelta_Intersection_Val( fDefault_Delta_Intersection_Val ),
57 fEpsilonMin( fEpsilonMinDefault ),
58 fEpsilonMax( fEpsilonMaxDefault)
60 if ( detectorField !=
nullptr )
66 fFieldChangesEnergy = fieldChangesEnergy;
71 G4cout <<
"G4FieldManager/ctor#1 fEpsilon Min/Max: eps_min = " << fEpsilonMin <<
" eps_max=" << fEpsilonMax <<
G4endl;
80 : fDetectorField(detectorField), fAllocatedChordFinder(true),
81 fDelta_One_Step_Value( fDefault_Delta_One_Step_Value ),
82 fDelta_Intersection_Val( fDefault_Delta_Intersection_Val ),
83 fEpsilonMin( fEpsilonMinDefault ),
84 fEpsilonMax( fEpsilonMaxDefault )
90 G4cout <<
"G4FieldManager/ctor#2 fEpsilon Min/Max: eps_min = " << fEpsilonMin <<
" eps_max=" << fEpsilonMax <<
G4endl;
103 if ( fDetectorField !=
nullptr )
105 aField = fDetectorField->
Clone();
111 aFM =
new G4FieldManager( aField ,
nullptr , fFieldChangesEnergy );
116 if ( fAllocatedChordFinder )
127 aFM->fChordFinder = aCF;
132 aFM->fEpsilonMax = fEpsilonMax;
133 aFM->fEpsilonMin = fEpsilonMin;
134 aFM->fDelta_Intersection_Val = fDelta_Intersection_Val;
135 aFM->fDelta_One_Step_Value = fDelta_One_Step_Value;
150 G4cout <<
"G4FieldManager/clone fEpsilon Min/Max: eps_min = " << fEpsilonMin <<
" eps_max=" << fEpsilonMax <<
G4endl;
161 if( fAllocatedChordFinder )
171 if ( fAllocatedChordFinder )
175 fAllocatedChordFinder =
false;
177 if( detectorMagField !=
nullptr )
180 fAllocatedChordFinder =
true;
184 fChordFinder =
nullptr;
188void G4FieldManager::InitialiseFieldChangesEnergy()
190 if ( fDetectorField !=
nullptr )
196 fFieldChangesEnergy =
false;
208 fDetectorField = pDetectorField;
209 InitialiseFieldChangesEnergy();
213 if( fChordFinder !=
nullptr )
215 failMode= std::max( failMode, 1) ;
218 driver = fChordFinder->GetIntegrationDriver();
219 if( driver !=
nullptr )
226 if( equation !=
nullptr )
234 if( !ableToSet && (failMode > 0) )
239 msg <<
"Unable to set the field in the dependent objects of G4FieldManager"
241 msg <<
"All the dependent classes must be fully initialised,"
242 <<
"before it is possible to call this method." <<
G4endl;
243 msg <<
"The problem encountered was the following: " <<
G4endl;
244 if( fChordFinder ==
nullptr ) { msg <<
" No ChordFinder. " ; }
245 else if( driver ==
nullptr ) { msg <<
" No Integration Driver set. ";}
246 else if( equation ==
nullptr ) { msg <<
" No Equation found. " ; }
248 else { msg <<
" Can NOT find reason for failure. ";}
252 G4Exception(
"G4FieldManager::SetDetectorField",
"Geometry001",
264 if(newEpsMax >= fEpsilonMin)
266 fEpsilonMax = newEpsMax;
270 G4cout <<
"G4FieldManager/SetEpsMax : eps_max = " << std::setw(10) << fEpsilonMax
271 <<
" ( Note: unchanged eps_min=" << std::setw(10) << fEpsilonMin <<
" )" <<
G4endl;
277 erm <<
" Call to set eps_max = " << newEpsMax <<
" . The problem is that"
278 <<
" its value must be at larger or equal to eps_min= " << fEpsilonMin <<
G4endl;
279 erm <<
" Modifying both to the same value " << newEpsMax <<
" to ensure consistency."
281 <<
" To avoid this warning, please set eps_min first, and ensure that "
284 fEpsilonMax = newEpsMax;
285 fEpsilonMin = newEpsMax;
309 fEpsilonMin = newEpsMin;
315 G4cout <<
"G4FieldManager/SetEpsMin : eps_min = "
316 << std::setw(10) << fEpsilonMin <<
G4endl;
318 if( fEpsilonMax < fEpsilonMin )
322 erm <<
"Setting eps_min = " << newEpsMin
323 <<
" For consistency set eps_max= " << fEpsilonMin
324 <<
" ( Old value = " << fEpsilonMax <<
" )" <<
G4endl;
325 fEpsilonMax = fEpsilonMin;
366 G4cout <<
"G4FieldManager::" << __func__
375 erm <<
"Proposed value for maximum-accepted-epsilon = " << maxAcceptValue
378 <<
"This may impact the robustness of integration of tracks in field."
381 <<
" , but future releases are expected " <<
G4endl
382 <<
" to tighten the limit of acceptable values to "
384 <<
"Suggestion: If you need better performance investigate using "
385 <<
"alternative, low-order RK integration methods or " <<
G4endl
386 <<
" helix-based methods (for pure B-fields) for low(er) energy tracks, "
387 <<
" especially electrons if you need better performance." <<
G4endl;
393 erm <<
" Proposed value for maximum accepted epsilon " << maxAcceptValue
398 erm <<
" Using the latter value instead." <<
G4endl;
405 erm <<
" NOTE: you can accept the ceiling value and turn this into a "
406 <<
" warning by using a 2nd argument " <<
G4endl
407 <<
" in your call to SetMaxAcceptedEpsilon: softFailure = true ";
415 G4Exception(methodName.c_str(),
"Geometry003", severity, erm);
425 erm <<
"Incorrect proposed value of " << name <<
" = " << value <<
G4endl
426 <<
" Its value is outside the permitted range from "
428 <<
" Clarification: " <<
G4endl;
429 G4long oldPrec = erm.precision();
432 erm <<
" a) The value must be positive and enough larger than the accuracy limit"
433 <<
" of the (G4)double type - ("
435 <<
" i.e. std::numeric_limits<G4double>::epsilon()= "
436 << std::numeric_limits<G4double>::epsilon()
437 <<
" to ensure that integration " <<
G4endl
438 <<
" could potentially achieve this acccuracy." <<
G4endl
443 erm <<
" b) It must be smaller than (or equal) " << std::setw(8)
445 <<
" to ensure robustness of integration - ("
450 G4bool badRoundoff = (std::fabs(1.0+value) == 1.0);
451 erm <<
" Unknown ERROR case -- extra check: " <<
G4endl;
452 erm <<
" c) as a floating point number (of type G4double) the sum (1+" << name
453 <<
" ) must be > 1 , ("
454 << (badRoundoff ?
"FAILED" :
"OK" ) <<
")" <<
G4endl
455 <<
" Now 1+eps_min = " << std::setw(20)
456 << std::setprecision(17) << (1+value) <<
G4endl
457 <<
" and (1.0+" << name <<
") - 1.0 = " << std::setw(20)
458 << std::setprecision(9) << (1.0+value)-1.0;
460 erm.precision(oldPrec);
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
G4GLOB_DLL std::ostream G4cout
void SetFieldObj(G4Field *pField)
static void DeRegister(G4FieldManager *pVolume)
static void Register(G4FieldManager *pVolume)
virtual G4FieldManager * Clone() const
static G4double GetMaxAcceptedEpsilon()
virtual ~G4FieldManager()
G4bool SetDetectorField(G4Field *detectorField, G4int failMode=0)
static constexpr G4double fMinAcceptedEpsilon
static constexpr G4double fMaxWarningEpsilon
void CreateChordFinder(G4MagneticField *detectorMagField)
G4bool SetMaximumEpsilonStep(G4double newEpsMax)
static G4double fMaxAcceptedEpsilon
void ReportBadEpsilonValue(G4ExceptionDescription &erm, G4double value, G4String &name) const
G4bool SetMinimumEpsilonStep(G4double newEpsMin)
virtual void ConfigureForTrack(const G4Track *)
G4FieldManager(G4Field *detectorField=nullptr, G4ChordFinder *pChordFinder=nullptr, G4bool b=true)
static G4bool fVerboseConstruction
static constexpr G4double fMaxFinalEpsilon
static G4bool SetMaxAcceptedEpsilon(G4double maxEps, G4bool softFail=false)
virtual G4bool DoesFieldChangeEnergy() const =0
virtual G4Field * Clone() const
virtual G4EquationOfMotion * GetEquationOfMotion()=0