Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4FieldManager.hh
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
27//
28// Class description:
29//
30// A class to manage (Store) a pointer to the Field subclass that
31// describes the field of a detector (magnetic, electric or other).
32// Also stores a reference to the chord finder.
33//
34// The G4FieldManager class exists to allow the user program to specify
35// the electric, magnetic and/or other field(s) of the detector.
36//
37// A field manager can be set to a logical volume (or to more than one),
38// in order to vary its field from that of the world. In this manner
39// a zero or constant field can override a global field, a more or
40// less exact version can override the external approximation, lower
41// or higher precision for tracking can be specified, a different
42// stepper can be chosen for different volumes, ...
43//
44// It also stores a pointer to the ChordFinder object that can do the
45// propagation in this field. All geometrical track "advancement"
46// in the field is handled by this ChordFinder object.
47//
48// G4FieldManager allows the other classes/object (of the MagneticField
49// & other class categories) to find out whether a detector field object
50// exists and what that object is.
51//
52// The Chord Finder must be created either by calling CreateChordFinder
53// for a Magnetic Field or by the user creating a Chord Finder object
54// "manually" and setting the pointer.
55//
56// A default FieldManager is created by the singleton class
57// G4NavigatorForTracking and exists before main is called.
58// However a new one can be created and given to G4NavigatorForTracking.
59//
60// Our current design envisions that one Field manager is
61// valid for each region detector.
62//
63// It is expected that a particular geometrical region has a Field manager.
64// By default a Field Manager is created for the world volume, and
65// will be utilised for all volumes unless it is overridden by a 'local'
66// field manager.
67// Note also that a region with both electric E and magnetic B field will
68// have these treated as one field.
69// Similarly it could be extended to treat other fields as additional
70// components of a single field type.
71
72// Author: John Apostolakis (CERN), 10.03.1997 - Design and implementation
73// -------------------------------------------------------------------
74#ifndef G4FIELDMANAGER_HH
75#define G4FIELDMANAGER_HH
76
77#include "G4FieldParameters.hh"
78#include "globals.hh"
79
80class G4Field;
81class G4MagneticField;
82class G4ChordFinder;
83class G4Track; // Forward reference for parameter configuration
84
85/**
86 * @brief G4FieldManager is a manager (store) for a pointer to the Field
87 * subclass that describes the field of a detector (magnetic, electric or
88 * other). It also stores a reference to the chord finder.
89 * A field manager can be set to a logical volume (or to more than one),
90 * in order to vary its field from that of the world volume. In this manner
91 * a zero or constant field can override a global field, a more or less exact
92 * version can override the external approximation, lower or higher precision
93 * for tracking can be specified, a different stepper can be chosen for
94 * different volumes, etc...
95 * The Chord Finder must be created either by calling CreateChordFinder()
96 * for a Magnetic Field or by the user creating a Chord Finder object
97 * "manually" and setting the pointer.
98 * The current design envisions that one Field manager is valid for each
99 * detector region. It is expected that a particular geometrical region has
100 * a Field manager. By default a Field Manager is created for the world volume,
101 * and will be utilised for all volumes unless it is overridden by a 'local'
102 * field manager.
103 * Note also that a region with both electric E and magnetic B field will
104 * have these treated as one field. Similarly it could be extended to treat
105 * other fields as additional components of a single field type.
106 */
107
109{
110 public:
111
112 /**
113 * General Constructor for any field. Must be set with field and chord finder
114 * for use.
115 * @param[in] detectorField Pointer to the field.
116 * @param[in] pChordFinder Pointer to the chord finder object.
117 * @param[in] b Flag to indicate if the field changes the energy; it is
118 * taken from the provided field, if specified.
119 */
120 G4FieldManager(G4Field* detectorField = nullptr,
121 G4ChordFinder* pChordFinder = nullptr,
122 G4bool b = true ); // fieldChangesEnergy is taken from field
123
124 /**
125 * Constructor creating the chord finder. It assumes pure magnetic field,
126 * so energy constant.
127 * @param[in] detectorMagneticField Pointer to the magnetic field.
128 */
129 G4FieldManager(G4MagneticField* detectorMagneticField);
130
131 /**
132 * Virtual Destructor.
133 */
134 virtual ~G4FieldManager();
135
136 /**
137 * Copy constructor and assignment operator not allowed.
138 */
141
142 /**
143 * Pushes the field to the equation. Failure to push the field (due to
144 * absence of a chord finder, driver, stepper or equation) is
145 * - '0' = quiet : Do not complain if chordFinder == 0
146 * (It will still warn for other error);
147 * - '1' = warn : a warning if anything is missing;
148 * - '2'/else = FATAL : a fatal error for all other values.
149 * @param[in] detectorField Pointer to the field.
150 * @param[in] failMode Flag (0/1/2) for selected failure mode.
151 * @returns Success (true) or failure (false).
152 */
153 G4bool SetDetectorField(G4Field* detectorField, G4int failMode = 0);
154
155 /**
156 * Pushes the field to this class only -- no further.
157 * Should be used to initialise this field, only *before* creating
158 * the chord finder and its dependent classes.
159 * User is then responsible to ensure that:
160 * i) an equation, stepper, driver and chord finder are created;
161 * ii) this field is used by the equation.
162 * @param[in] detectorField Pointer to the field.
163 */
164 inline void ProposeDetectorField(G4Field* detectorField);
165
166 /**
167 * Pushes the field to the equation and keeps its address.
168 * Can be used only once the equation, stepper, driver and chord finder
169 * have all been created; else it is an error.
170 * @param[in] detectorField Pointer to the field.
171 */
172 inline void ChangeDetectorField(G4Field* detectorField);
173
174 /**
175 * Methods to get and check (existance of) the field object.
176 */
177 inline const G4Field* GetDetectorField() const;
178 inline G4bool DoesFieldExist() const;
179
180 /**
181 * Methods to create, set or get the associated Chord Finder.
182 */
183 void CreateChordFinder(G4MagneticField* detectorMagField);
184 inline void SetChordFinder(G4ChordFinder* aChordFinder);
186 inline const G4ChordFinder* GetChordFinder() const;
187
188 /**
189 * Setups the choice of the configurable parameters, relying on the
190 * current track's energy, particle identity...
191 * Note: in addition to the values of member variables, a user can use
192 * this to change the ChordFinder, the field, etc.
193 * @param[in] pTrack Pointer to a track.
194 */
195 virtual void ConfigureForTrack( const G4Track* pTrack );
196
197 /**
198 * Static methods to set/get the global field.
199 */
200 static void SetGlobalFieldManager(G4FieldManager* fieldManager);
202
203 /**
204 * Returns the accuracy for boundary intersection.
205 */
207
208 /**
209 * Returns the accuracy for one tracking/physics step.
210 */
211 inline G4double GetDeltaOneStep() const;
212
213 /**
214 * Sets both accuracies, maintaining a fixed ratio for accuracies
215 * of volume Intersection and Integration (in One Step).
216 */
217 inline void SetAccuraciesWithDeltaOneStep(G4double valDeltaOneStep);
218
219 /**
220 * Sets the accuracy for integration of one step (only).
221 */
222 inline void SetDeltaOneStep(G4double valueD1step);
223
224 /**
225 * Sets the accuracy of intersection of a volume (only).
226 */
227 inline void SetDeltaIntersection(G4double valueDintersection);
228
229 /**
230 * Methods to set/get the minimum for Relative accuracy of a Step.
231 */
234
235 /**
236 * Methods to set/get the maximum for Relative accuracy of a Step.
237 */
240
241 /**
242 * Methods to set/get flag for field changing energy.
243 * For electric field this should be true; for magnetic field this
244 * should be false.
245 */
247 inline void SetFieldChangesEnergy(G4bool value);
248
249 /**
250 * Needed for multi-threading, create and returns an allocated clone
251 * of this object.
252 */
253 virtual G4FieldManager* Clone() const;
254
255 /**
256 * Static methods to set/get the maximum accepted epsilon.
257 * If setting fails, with softFail=true it gives Warning, else
258 * a FatalException.
259 */
261 static G4bool SetMaxAcceptedEpsilon(G4double maxEps, G4bool softFail= false);
262
263 protected:
264
265 /**
266 * Logger for reporting on correctness of the proposed epsilon value.
267 */
269 const G4String& name) const;
270
271 /** Epsilon_min/max values must be smaller than this for robust integration. */
273 static constexpr G4double fMinAcceptedEpsilon = 1000.0 * std::numeric_limits<G4double>::epsilon();
274
275 /** Setting larger value will give warning. */
276 static constexpr G4double fMaxWarningEpsilon = 0.001;
277
278 /** Will not accept larger values. */
279 static constexpr G4double fMaxFinalEpsilon = 0.02;
280
281 /** Controls verbosity of constructors. */
283
284 private:
285
286 /**
287 * Checks whether the field/equation changes the energy and sets the data
288 * member accordingly. Does not handle special cases - this must be done
289 * separately (e.g. magnetic monopole in B field).
290 */
291 void InitialiseFieldChangesEnergy();
292
293 private:
294
295 /** Dependent objects -- with state that depends on tracking. */
296 G4Field* fDetectorField = nullptr;
297 G4ChordFinder* fChordFinder = nullptr;
298
299 /** Flag to indicate if "new" was used to create the Chord Finder. */
300 G4bool fAllocatedChordFinder = false; //
301
302 // 1. CHARACTERISTIC of field
303
304 G4bool fFieldChangesEnergy = false;
305
306 // 2. PARAMETERS that determine the accuracy of integration or intersection
307
308 /** Value for the required accuracies for one tracking/physics step. */
309 G4double fDelta_One_Step_Value = G4FieldDefaults::kDeltaOneStep;
310
311 /** Value for the required accuracies for boundary intersection. */
312 G4double fDelta_Intersection_Val = G4FieldDefaults::kDeltaIntersection;
313
314 /** Values for the small possible relative accuracy of a step
315 (corresponding to the greatest possible integration accuracy). */
318
319 /** Global field manager set by G4TransportationManager to allow accessing
320 the global field without dependency on navigation. */
321 static G4ThreadLocal G4FieldManager* fGlobalFieldManager;
322};
323
324// Implementation of inline functions
325
326#include "G4FieldManager.icc"
327
328#endif
std::ostringstream G4ExceptionDescription
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4ChordFinder is a class that provides Runge-Kutta integration of motion ODE and also has a method th...
void SetAccuraciesWithDeltaOneStep(G4double valDeltaOneStep)
virtual G4FieldManager * Clone() const
static G4double GetMaxAcceptedEpsilon()
G4bool DoesFieldChangeEnergy() const
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
void SetFieldChangesEnergy(G4bool value)
static G4double fMaxAcceptedEpsilon
static G4FieldManager * GetGlobalFieldManager()
void SetDeltaOneStep(G4double valueD1step)
void SetChordFinder(G4ChordFinder *aChordFinder)
G4bool SetMinimumEpsilonStep(G4double newEpsMin)
void ProposeDetectorField(G4Field *detectorField)
G4double GetMinimumEpsilonStep() const
const G4ChordFinder * GetChordFinder() const
G4FieldManager(G4Field *detectorField=nullptr, G4ChordFinder *pChordFinder=nullptr, G4bool b=true)
G4double GetMaximumEpsilonStep() const
static G4bool fVerboseConstruction
static constexpr G4double fMaxFinalEpsilon
static void SetGlobalFieldManager(G4FieldManager *fieldManager)
virtual void ConfigureForTrack(const G4Track *pTrack)
G4FieldManager & operator=(const G4FieldManager &)=delete
G4double GetDeltaOneStep() const
static G4bool SetMaxAcceptedEpsilon(G4double maxEps, G4bool softFail=false)
const G4Field * GetDetectorField() const
void ChangeDetectorField(G4Field *detectorField)
void SetDeltaIntersection(G4double valueDintersection)
G4FieldManager(const G4FieldManager &)=delete
G4ChordFinder * GetChordFinder()
G4double GetDeltaIntersection() const
G4bool DoesFieldExist() const
G4Field is the abstract class for any kind of field. It allows any kind of field (vector,...
Definition G4Field.hh:67
constexpr G4double kDeltaOneStep
Default delta one step in global field manager.
constexpr G4double kMaximumEpsilonStep
Default maximum epsilon step in global field manager.
constexpr G4double kDeltaIntersection
Delta intersection in global field manager.
constexpr G4double kMinimumEpsilonStep
Default minimum epsilon step in global field manager.
#define G4ThreadLocal
Definition tls.hh:77