Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4FieldManager.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// G4FieldManager implementation
27//
28// Author: John Apostolakis (CERN), 10.03.1997 - Design and implementation
29// -------------------------------------------------------------------
30
31#include "G4FieldManager.hh"
32#include "G4Field.hh"
33#include "G4MagneticField.hh"
34#include "G4ChordFinder.hh"
36#include "G4SystemOfUnits.hh"
37
38G4ThreadLocal G4FieldManager* G4FieldManager::fGlobalFieldManager = nullptr;
40
42// Legacy value. Future value = 0.001
43// Requesting a large epsilon (max) value provides poor accuracy for
44// every integration segment.
45// Problems occur because some methods (including DormandPrince(7)45)
46// the estimation of local error appears to be a substantial underestimate
47// at large epsilon values ( > 0.001 ). So the value for fMaxAcceptedEpsilon
48// is recommended to be 0.001 or below.
49
51 G4ChordFinder* pChordFinder,
52 G4bool fieldChangesEnergy )
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}
74
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}
88
90{
91 if( fAllocatedChordFinder )
92 {
93 delete fChordFinder;
94 }
96}
97
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}
155
157{
158 // Default is to do nothing!
159}
160
161void
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}
180
181void G4FieldManager::InitialiseFieldChangesEnergy()
182{
183 if ( fDetectorField != nullptr )
184 {
185 fFieldChangesEnergy = fDetectorField->DoesFieldChangeEnergy();
186 }
187 else
188 {
189 fFieldChangesEnergy = false; // No field, no change!
190 }
191}
192
194 G4int failMode )
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}
250
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}
293
294// -----------------------------------------------------------------------------
295
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}
333
334// -----------------------------------------------------------------------------
335
340
341// -----------------------------------------------------------------------------
342
344// Set value -- within limits
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}
412
413// -----------------------------------------------------------------------------
414
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}
G4ExceptionSeverity
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
double G4double
Definition G4Types.hh:83
long G4long
Definition G4Types.hh:87
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4ChordFinder is a class that provides Runge-Kutta integration of motion ODE and also has a method th...
G4EquationOfMotion is the abstract base class for the right hand size of the equation of motion of a ...
void SetFieldObj(G4Field *pField)
static void DeRegister(G4FieldManager *pFieldMan)
static void Register(G4FieldManager *pFieldMan)
G4FieldManager is a manager (store) for a pointer to the Field subclass that describes the field of a...
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)
void ReportBadEpsilonValue(G4ExceptionDescription &erm, G4double value, const G4String &name) const
static G4double fMaxAcceptedEpsilon
G4bool SetMinimumEpsilonStep(G4double newEpsMin)
G4FieldManager(G4Field *detectorField=nullptr, G4ChordFinder *pChordFinder=nullptr, G4bool b=true)
static G4bool fVerboseConstruction
static constexpr G4double fMaxFinalEpsilon
virtual void ConfigureForTrack(const G4Track *pTrack)
static G4bool SetMaxAcceptedEpsilon(G4double maxEps, G4bool softFail=false)
G4Field is the abstract class for any kind of field. It allows any kind of field (vector,...
Definition G4Field.hh:67
virtual G4bool DoesFieldChangeEnergy() const =0
virtual G4Field * Clone() const
Definition G4Field.cc:45
virtual G4EquationOfMotion * GetEquationOfMotion()=0
#define G4ThreadLocal
Definition tls.hh:77