Geant4 11.2.2
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, 10.03.97 - 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
38G4double G4FieldManager::fDefault_Delta_One_Step_Value= 0.01 * millimeter;
39G4double G4FieldManager::fDefault_Delta_Intersection_Val= 0.001 * millimeter;
41
42G4double G4FieldManager::fMaxAcceptedEpsilon= 0.01; // 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 the estimation of local
46// error appears to be a substantial underestimate at large epsilon values ( > 0.001 ).
47// So the value for fMaxAcceptedEpsilon is recommended to be 0.001 or below.
48
50 G4ChordFinder* pChordFinder,
51 G4bool fieldChangesEnergy )
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 {
70 G4cout << "G4FieldManager/ctor#1 fEpsilon Min/Max: eps_min = " << fEpsilonMin << " eps_max=" << fEpsilonMax << G4endl;
71 }
72
73 // Add to store
74 //
76}
77
79 : fDetectorField(detectorField), fAllocatedChordFinder(true),
80 fDelta_One_Step_Value( fDefault_Delta_One_Step_Value ),
81 fDelta_Intersection_Val( fDefault_Delta_Intersection_Val ),
82 fEpsilonMin( fEpsilonMinDefault ),
83 fEpsilonMax( fEpsilonMaxDefault )
84{
85 fChordFinder = new G4ChordFinder( detectorField );
86
88 {
89 G4cout << "G4FieldManager/ctor#2 fEpsilon Min/Max: eps_min = " << fEpsilonMin << " eps_max=" << fEpsilonMax << G4endl;
90 }
91 // Add to store
92 //
94}
95
97{
98 G4Field* aField = nullptr;
99 G4FieldManager* aFM = nullptr;
100 G4ChordFinder* aCF = nullptr;
101 try {
102 if ( fDetectorField != nullptr )
103 {
104 aField = fDetectorField->Clone();
105 }
106
107 // Create a new field manager, note that we do not set
108 // any chordfinder now
109 //
110 aFM = new G4FieldManager( aField , nullptr , fFieldChangesEnergy );
111
112 // Check if originally we have the fAllocatedChordFinder variable
113 // set, in case, call chord constructor
114 //
115 if ( fAllocatedChordFinder )
116 {
117 aFM->CreateChordFinder( dynamic_cast<G4MagneticField*>(aField) );
118 }
119 else
120 {
121 // Chord was specified by user, should we clone?
122 // TODO: For the moment copy pointer, to be understood
123 // if cloning of ChordFinder is needed
124 //
125 aCF = fChordFinder; /*->Clone*/
126 aFM->fChordFinder = aCF;
127 }
128
129 // Copy values of other variables
130
131 aFM->fEpsilonMax = fEpsilonMax;
132 aFM->fEpsilonMin = fEpsilonMin;
133 aFM->fDelta_Intersection_Val = fDelta_Intersection_Val;
134 aFM->fDelta_One_Step_Value = fDelta_One_Step_Value;
135 // TODO: Should we really add to the store the cloned FM?
136 // Who will use this?
137 }
138 catch ( ... )
139 {
140 // Failed creating clone: probably user did not implement Clone method
141 // in derived classes?
142 // Perform clean-up after ourselves...
143 delete aField;
144 delete aFM;
145 delete aCF;
146 throw;
147 }
148
149 G4cout << "G4FieldManager/clone fEpsilon Min/Max: eps_min = " << fEpsilonMin << " eps_max=" << fEpsilonMax << G4endl;
150 return aFM;
151}
152
154{
155 // Default is to do nothing!
156}
157
159{
160 if( fAllocatedChordFinder )
161 {
162 delete fChordFinder;
163 }
165}
166
167void
169{
170 if ( fAllocatedChordFinder )
171 {
172 delete fChordFinder;
173 }
174 fAllocatedChordFinder = false;
175
176 if( detectorMagField != nullptr )
177 {
178 fChordFinder = new G4ChordFinder( detectorMagField );
179 fAllocatedChordFinder = true;
180 }
181 else
182 {
183 fChordFinder = nullptr;
184 }
185}
186
187void G4FieldManager::InitialiseFieldChangesEnergy()
188{
189 if ( fDetectorField != nullptr )
190 {
191 fFieldChangesEnergy = fDetectorField->DoesFieldChangeEnergy();
192 }
193 else
194 {
195 fFieldChangesEnergy = false; // No field, no change!
196 }
197}
198
200 G4int failMode )
201{
202 G4VIntegrationDriver* driver = nullptr;
203 G4EquationOfMotion* equation = nullptr;
204 // G4bool compatibleField = false;
205 G4bool ableToSet = false;
206
207 fDetectorField = pDetectorField;
208 InitialiseFieldChangesEnergy();
209
210 // Must 'propagate' the field to the dependent classes
211 //
212 if( fChordFinder != nullptr )
213 {
214 failMode= std::max( failMode, 1) ;
215 // If a chord finder exists, warn in case of error!
216
217 driver = fChordFinder->GetIntegrationDriver();
218 if( driver != nullptr )
219 {
220 equation = driver->GetEquationOfMotion();
221
222 // Should check the compatibility between the
223 // field and the equation HERE
224
225 if( equation != nullptr )
226 {
227 equation->SetFieldObj(pDetectorField);
228 ableToSet = true;
229 }
230 }
231 }
232
233 if( !ableToSet && (failMode > 0) )
234 {
235 // If this fails, report the issue !
236
238 msg << "Unable to set the field in the dependent objects of G4FieldManager"
239 << G4endl;
240 msg << "All the dependent classes must be fully initialised,"
241 << "before it is possible to call this method." << G4endl;
242 msg << "The problem encountered was the following: " << G4endl;
243 if( fChordFinder == nullptr ) { msg << " No ChordFinder. " ; }
244 else if( driver == nullptr ) { msg << " No Integration Driver set. ";}
245 else if( equation == nullptr ) { msg << " No Equation found. " ; }
246 // else if( !compatibleField ) { msg << " Field not compatible. ";}
247 else { msg << " Can NOT find reason for failure. ";}
248 msg << G4endl;
249 G4ExceptionSeverity severity = (failMode != 1)
251 G4Exception("G4FieldManager::SetDetectorField", "Geometry001",
252 severity, msg);
253 }
254 return ableToSet;
255}
256
258{
259 G4bool succeeded= false;
260 if( (newEpsMax > 0.0) && ( newEpsMax <= fMaxAcceptedEpsilon)
261 && (fMinAcceptedEpsilon <= newEpsMax ) ) // (std::fabs(1.0+newEpsMax)>1.0) )
262 {
263 if(newEpsMax >= fEpsilonMin)
264 {
265 fEpsilonMax = newEpsMax;
266 succeeded = true;
268 {
269 G4cout << "G4FieldManager/SetEpsMax : eps_max = " << std::setw(10) << fEpsilonMax
270 << " ( Note: unchanged eps_min=" << std::setw(10) << fEpsilonMin << " )" << G4endl;
271 }
272 }
273 else
274 {
276 erm << " Call to set eps_max = " << newEpsMax << " . The problem is that"
277 << " its value must be at larger or equal to eps_min= " << fEpsilonMin << G4endl;
278 erm << " Modifying both to the same value " << newEpsMax << " to ensure consistency."
279 << G4endl
280 << " To avoid this warning, please set eps_min first, and ensure that "
281 << " 0 < eps_min <= eps_max <= " << fMaxAcceptedEpsilon << G4endl;
282
283 fEpsilonMax = newEpsMax;
284 fEpsilonMin = newEpsMax;
285 G4String methodName = G4String("G4FieldManager::")+ G4String(__func__);
286 G4Exception(methodName.c_str(), "Geometry003", JustWarning, erm);
287 }
288 }
289 else
290 {
292 G4String paramName("eps_max");
293 ReportBadEpsilonValue(erm, newEpsMax, paramName );
294 G4String methodName = G4String("G4FieldManager::")+ G4String(__func__);
295 G4Exception(methodName.c_str(), "Geometry001", FatalException, erm);
296 }
297 return succeeded;
298}
299
300// -----------------------------------------------------------------------------
301
303{
304 G4bool succeeded= false;
305
306 if( fMinAcceptedEpsilon <= newEpsMin && newEpsMin <= fMaxAcceptedEpsilon )
307 {
308 fEpsilonMin = newEpsMin;
309 //*********
310 succeeded= true;
311
313 {
314 G4cout << "G4FieldManager/SetEpsMin : eps_min = "
315 << std::setw(10) << fEpsilonMin << G4endl;
316 }
317 if( fEpsilonMax < fEpsilonMin )
318 {
319 // Ensure consistency
321 erm << "Setting eps_min = " << newEpsMin
322 << " For consistency set eps_max= " << fEpsilonMin
323 << " ( Old value = " << fEpsilonMax << " )" << G4endl;
324 fEpsilonMax = fEpsilonMin;
325 G4String methodName = G4String("G4FieldManager::")+ G4String(__func__);
326 G4Exception(methodName.c_str(), "Geometry003", JustWarning, erm);
327 }
328 }
329 else
330 {
332 G4String paramName("eps_min");
333 ReportBadEpsilonValue(erm, newEpsMin, paramName );
334 G4String methodName = G4String("G4FieldManager::")+ G4String(__func__);
335 G4Exception(methodName.c_str(), "Geometry001", FatalException, erm);
336 }
337 return succeeded;
338}
339
340// -----------------------------------------------------------------------------
341
346
347// -----------------------------------------------------------------------------
348
350// Set value -- within limits
351{
352 G4bool success= false;
353 // Limit for warning and absolute limit chosen from experience in and
354 // investigation of integration with G4DormandPrince745 in HEP-type setups.
355 if( maxAcceptValue <= fMaxWarningEpsilon )
356 {
357 fMaxAcceptedEpsilon= maxAcceptValue;
358 success= true;
359 }
360 else
361 {
363 G4ExceptionSeverity severity;
364
365 G4cout << "G4FieldManager::" << __func__
366 << " Parameters: fMaxAcceptedEpsilon = " << fMaxAcceptedEpsilon
367 << " fMaxFinalEpsilon = " << fMaxFinalEpsilon << G4endl;
368
369 if( maxAcceptValue <= fMaxFinalEpsilon )
370 {
371 success= true;
372 fMaxAcceptedEpsilon = maxAcceptValue;
373 // Integration is poor, and robustness will likely suffer
374 erm << "Proposed value for maximum-accepted-epsilon = " << maxAcceptValue
375 << " is larger than the recommended = " << fMaxWarningEpsilon
376 << G4endl
377 << "This may impact the robustness of integration of tracks in field."
378 << G4endl
379 << "The request was accepted and the value = " << fMaxAcceptedEpsilon
380 << " , but future releases are expected " << G4endl
381 << " to tighten the limit of acceptable values to "
383 << "Suggestion: If you need better performance investigate using "
384 << "alternative, low-order RK integration methods or " << G4endl
385 << " helix-based methods (for pure B-fields) for low(er) energy tracks, "
386 << " especially electrons if you need better performance." << G4endl;
387 severity= JustWarning;
388 }
389 else
390 {
392 erm << " Proposed value for maximum accepted epsilon " << maxAcceptValue
393 << " is larger than the top of the range = " << fMaxFinalEpsilon
394 << G4endl;
395 if( softFailure )
396 {
397 erm << " Using the latter value instead." << G4endl;
398 }
399 erm << G4endl;
400 erm << " Please adjust to request maxAccepted <= " << fMaxFinalEpsilon
401 << G4endl << G4endl;
402 if( !softFailure )
403 {
404 erm << " NOTE: you can accept the ceiling value and turn this into a "
405 << " warning by using a 2nd argument " << G4endl
406 << " in your call to SetMaxAcceptedEpsilon: softFailure = true ";
407 }
408 severity = softFailure ? JustWarning : FatalException;
409 // if( softFailure ) severity= JustWarning;
410 // else severity= FatalException;
411 success = false;
412 }
413 G4String methodName = G4String("G4FieldManager::")+ G4String(__func__);
414 G4Exception(methodName.c_str(), "Geometry003", severity, erm);
415 }
416 return success;
417}
418
419// -----------------------------------------------------------------------------
420
423{
424 erm << "Incorrect proposed value of " << name << " = " << value << G4endl
425 << " Its value is outside the permitted range from "
427 << " Clarification: " << G4endl;
428 G4long oldPrec = erm.precision();
429 if(value < fMinAcceptedEpsilon )
430 {
431 erm << " a) The value must be positive and enough larger than the accuracy limit"
432 << " of the (G4)double type - ("
433 << (value < fMinAcceptedEpsilon ? "FAILED" : "OK" ) << ")" << G4endl
434 << " i.e. std::numeric_limits<G4double>::epsilon()= "
435 << std::numeric_limits<G4double>::epsilon()
436 << " to ensure that integration " << G4endl
437 << " could potentially achieve this acccuracy." << G4endl
438 << " Minimum accepted eps_min/max value = " << fMinAcceptedEpsilon << G4endl;
439 }
440 else if( value > fMaxAcceptedEpsilon)
441 {
442 erm << " b) It must be smaller than (or equal) " << std::setw(8)
443 << std::setprecision(4) << fMaxAcceptedEpsilon
444 << " to ensure robustness of integration - ("
445 << (( value < fMaxAcceptedEpsilon) ? "OK" : "FAILED" ) << ")" << G4endl;
446 }
447 else
448 {
449 G4bool badRoundoff = (std::fabs(1.0+value) == 1.0);
450 erm << " Unknown ERROR case -- extra check: " << G4endl;
451 erm << " c) as a floating point number (of type G4double) the sum (1+" << name
452 << " ) must be > 1 , ("
453 << (badRoundoff ? "FAILED" : "OK" ) << ")" << G4endl
454 << " Now 1+eps_min = " << std::setw(20)
455 << std::setprecision(17) << (1+value) << G4endl
456 << " and (1.0+" << name << ") - 1.0 = " << std::setw(20)
457 << std::setprecision(9) << (1.0+value)-1.0;
458 }
459 erm.precision(oldPrec);
460}
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
G4VIntegrationDriver * GetIntegrationDriver()
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
Definition G4Field.cc:49
virtual G4EquationOfMotion * GetEquationOfMotion()=0