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