Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 a Chord Finder object
54// "manually" and setting this 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, 10.03.97 - design and implementation
73// -------------------------------------------------------------------
74#ifndef G4FIELDMANAGER_HH
75#define G4FIELDMANAGER_HH 1
76
77#include "globals.hh"
78
79class G4Field;
80class G4MagneticField;
81class G4ChordFinder;
82class G4Track; // Forward reference for parameter configuration
83
85{
86 public: // with description
87
88 G4FieldManager(G4Field* detectorField = nullptr,
89 G4ChordFinder* pChordFinder = nullptr,
90 G4bool b = true ); // fieldChangesEnergy is taken from field
91 // General constructor for any field.
92 // -> Must be set with field and chordfinder for use.
93 G4FieldManager(G4MagneticField* detectorMagneticField);
94 // Creates ChordFinder
95 // -> Assumes pure magnetic field (so energy constant)
96
97 virtual ~G4FieldManager();
98
101
102 G4bool SetDetectorField(G4Field* detectorField, G4int failMode = 0);
103 // Pushes the field to the equation.
104 // Failure to push the field (due to absence of a chord finder, driver,
105 // stepper or equation) is
106 // - '0' = quiet : Do not complain if chordFinder == 0
107 // (It will still warn for other error.)
108 // - '1' = warn : a warning if anything is missing
109 // - '2'/else = FATAL : a fatal error for all other values.
110 // Returns success (true) or failure (false)
111
112 inline void ProposeDetectorField(G4Field* detectorField);
113 // Pushes the field to this class only -- no further.
114 // Should be used to initialise this field, only *before* creating
115 // the chord finder and its dependent classes.
116 // User is then responsible to ensure that:
117 // i) an equation, stepper, driver and chord finder are created
118 // ii) this field is used by the equation.
119
120 inline void ChangeDetectorField(G4Field* detectorField);
121 // Pushes the field to the equation ( & keeps its address )
122 // Can be used only once the equation, stepper, driver and chord finder
123 // have all been created. Else it is an error.
124
125 inline const G4Field* GetDetectorField() const;
126 inline G4bool DoesFieldExist() const;
127 // Set, get and check the field object
128
129 void CreateChordFinder(G4MagneticField* detectorMagField);
130 inline void SetChordFinder(G4ChordFinder* aChordFinder);
132 inline const G4ChordFinder* GetChordFinder() const;
133 // Create, set or get the associated Chord Finder
134
135 virtual void ConfigureForTrack( const G4Track * );
136 // Setup the choice of the configurable parameters
137 // relying on the current track's energy, particle identity, ..
138 // Note: in addition to the values of member variables,
139 // a user can use this to change the ChordFinder, the field, ...
140
141 public: // with description
142
144 // Accuracy for boundary intersection.
145
146 inline G4double GetDeltaOneStep() const;
147 // Accuracy for one tracking/physics step.
148
149 inline void SetAccuraciesWithDeltaOneStep(G4double valDeltaOneStep);
150 // Sets both accuracies, maintaining a fixed ratio for accuracies
151 // of volume Intersection and Integration (in One Step)
152
153 inline void SetDeltaOneStep(G4double valueD1step);
154 // Set accuracy for integration of one step. (only)
155 inline void SetDeltaIntersection(G4double valueDintersection);
156 // Set accuracy of intersection of a volume. (only)
157
160 // Minimum for Relative accuracy of a Step
161
164 // Maximum for Relative accuracy of a Step
165
167 inline void SetFieldChangesEnergy(G4bool value);
168 // For electric field this should be true
169 // For magnetic field this should be false
170
171 virtual G4FieldManager* Clone() const;
172 // Needed for multi-threading, create a clone of this object
173
174 public:
176 static G4bool SetMaxAcceptedEpsilon(G4double maxEps, G4bool softFail= false);
177 // Set value -- within limits.
178 // If it fails, with softFail=true it gives Warning, else FatalException
179
180 protected:
182 static constexpr G4double fMinAcceptedEpsilon= 1000.0 * std::numeric_limits<G4double>::epsilon();
183 // Epsilon_min/max values must be smaller than this - for robust integration
184
185 static constexpr G4double fMaxWarningEpsilon= 0.001; // Setting larger value will give warning.
186 static constexpr G4double fMaxFinalEpsilon= 0.02; // Will not accept larger values
187
189 // Control verbosity of constructors
190
191 private:
192
193 void InitialiseFieldChangesEnergy();
194 // Check whether field/equation change the energy,
195 // and sets the data member accordingly
196 // Note: does not handle special cases - this must be done
197 // separately (e.g. magnetic monopole in B field )
198
199 protected:
201 G4String& name) const;
202
203 private:
204
205 G4Field* fDetectorField = nullptr;
206 G4ChordFinder* fChordFinder = nullptr;
207 // Dependent objects -- with state that depends on tracking
208
209 G4bool fAllocatedChordFinder = false; // Did we used "new" to
210 // create fChordFinder ?
211 // INVARIANTS of tracking ---------------------------------------
212 //
213 // 1. 'CONSTANTS' - default values for accuracy parameters
214 //
215 const G4double fEpsilonMinDefault= 5.0e-5; // Expected: 5.0e-5 to 1.0e-10 ...
216 const G4double fEpsilonMaxDefault= 1.0e-3; // Expected: 1.0e-3 to 1.0e-8 ...
217
218 static G4double fDefault_Delta_One_Step_Value; // = 0.01 * millimeter;
219 static G4double fDefault_Delta_Intersection_Val; // = 0.001 * millimeter;
220 // Default values for accuracy parameters
221
222 // 2. CHARACTERISTIC of field
223 //
224 G4bool fFieldChangesEnergy = false;
225
226 // 3. PARAMETERS that determine the accuracy of integration or intersection
227 //
228 G4double fDelta_One_Step_Value; // for one tracking/physics step
229 G4double fDelta_Intersection_Val; // for boundary intersection
230 // Values for the required accuracies
231
232 G4double fEpsilonMin;
233 G4double fEpsilonMax;
234 // Values for the small possible relative accuracy of a step
235 // (corresponding to the greatest possible integration accuracy)
236};
237
238// Implementation of inline functions
239
240#include "G4FieldManager.icc"
241
242#endif
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
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 SetFieldChangesEnergy(G4bool value)
static G4double fMaxAcceptedEpsilon
void SetDeltaOneStep(G4double valueD1step)
void ReportBadEpsilonValue(G4ExceptionDescription &erm, G4double value, G4String &name) const
void SetChordFinder(G4ChordFinder *aChordFinder)
G4bool SetMinimumEpsilonStep(G4double newEpsMin)
virtual void ConfigureForTrack(const G4Track *)
void ProposeDetectorField(G4Field *detectorField)
G4double GetMinimumEpsilonStep() const
const G4ChordFinder * GetChordFinder() const
G4double GetMaximumEpsilonStep() const
static G4bool fVerboseConstruction
static constexpr G4double fMaxFinalEpsilon
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