|
||||
File indexing completed on 2025-01-18 09:58:17
0001 // 0002 // ******************************************************************** 0003 // * License and Disclaimer * 0004 // * * 0005 // * The Geant4 software is copyright of the Copyright Holders of * 0006 // * the Geant4 Collaboration. It is provided under the terms and * 0007 // * conditions of the Geant4 Software License, included in the file * 0008 // * LICENSE and available at http://cern.ch/geant4/license . These * 0009 // * include a list of copyright holders. * 0010 // * * 0011 // * Neither the authors of this software system, nor their employing * 0012 // * institutes,nor the agencies providing financial support for this * 0013 // * work make any representation or warranty, express or implied, * 0014 // * regarding this software system or assume any liability for its * 0015 // * use. Please see the license in the file LICENSE and URL above * 0016 // * for the full disclaimer and the limitation of liability. * 0017 // * * 0018 // * This code implementation is the result of the scientific and * 0019 // * technical work of the GEANT4 collaboration. * 0020 // * By using, copying, modifying or distributing the software (or * 0021 // * any work based on the software) you agree to acknowledge its * 0022 // * use in resulting scientific publications, and indicate your * 0023 // * acceptance of all terms of the Geant4 Software license. * 0024 // ******************************************************************** 0025 // 0026 // G4FieldManager 0027 // 0028 // Class description: 0029 // 0030 // A class to manage (Store) a pointer to the Field subclass that 0031 // describes the field of a detector (magnetic, electric or other). 0032 // Also stores a reference to the chord finder. 0033 // 0034 // The G4FieldManager class exists to allow the user program to specify 0035 // the electric, magnetic and/or other field(s) of the detector. 0036 // 0037 // A field manager can be set to a logical volume (or to more than one), 0038 // in order to vary its field from that of the world. In this manner 0039 // a zero or constant field can override a global field, a more or 0040 // less exact version can override the external approximation, lower 0041 // or higher precision for tracking can be specified, a different 0042 // stepper can be chosen for different volumes, ... 0043 // 0044 // It also stores a pointer to the ChordFinder object that can do the 0045 // propagation in this field. All geometrical track "advancement" 0046 // in the field is handled by this ChordFinder object. 0047 // 0048 // G4FieldManager allows the other classes/object (of the MagneticField 0049 // & other class categories) to find out whether a detector field object 0050 // exists and what that object is. 0051 // 0052 // The Chord Finder must be created either by calling CreateChordFinder 0053 // for a Magnetic Field or by the user creating a a Chord Finder object 0054 // "manually" and setting this pointer. 0055 // 0056 // A default FieldManager is created by the singleton class 0057 // G4NavigatorForTracking and exists before main is called. 0058 // However a new one can be created and given to G4NavigatorForTracking. 0059 // 0060 // Our current design envisions that one Field manager is 0061 // valid for each region detector. 0062 // 0063 // It is expected that a particular geometrical region has a Field manager. 0064 // By default a Field Manager is created for the world volume, and 0065 // will be utilised for all volumes unless it is overridden by a 'local' 0066 // field manager. 0067 // Note also that a region with both electric E and magnetic B field will 0068 // have these treated as one field. 0069 // Similarly it could be extended to treat other fields as additional 0070 // components of a single field type. 0071 0072 // Author: John Apostolakis, 10.03.97 - design and implementation 0073 // ------------------------------------------------------------------- 0074 #ifndef G4FIELDMANAGER_HH 0075 #define G4FIELDMANAGER_HH 1 0076 0077 #include "globals.hh" 0078 0079 class G4Field; 0080 class G4MagneticField; 0081 class G4ChordFinder; 0082 class G4Track; // Forward reference for parameter configuration 0083 0084 class G4FieldManager 0085 { 0086 public: // with description 0087 0088 G4FieldManager(G4Field* detectorField = nullptr, 0089 G4ChordFinder* pChordFinder = nullptr, 0090 G4bool b = true ); // fieldChangesEnergy is taken from field 0091 // General constructor for any field. 0092 // -> Must be set with field and chordfinder for use. 0093 G4FieldManager(G4MagneticField* detectorMagneticField); 0094 // Creates ChordFinder 0095 // -> Assumes pure magnetic field (so energy constant) 0096 0097 virtual ~G4FieldManager(); 0098 0099 G4FieldManager(const G4FieldManager&) = delete; 0100 G4FieldManager& operator=(const G4FieldManager&) = delete; 0101 0102 G4bool SetDetectorField(G4Field* detectorField, G4int failMode = 0); 0103 // Pushes the field to the equation. 0104 // Failure to push the field (due to absence of a chord finder, driver, 0105 // stepper or equation) is 0106 // - '0' = quiet : Do not complain if chordFinder == 0 0107 // (It will still warn for other error.) 0108 // - '1' = warn : a warning if anything is missing 0109 // - '2'/else = FATAL : a fatal error for all other values. 0110 // Returns success (true) or failure (false) 0111 0112 inline void ProposeDetectorField(G4Field* detectorField); 0113 // Pushes the field to this class only -- no further. 0114 // Should be used to initialise this field, only *before* creating 0115 // the chord finder and its dependent classes. 0116 // User is then responsible to ensure that: 0117 // i) an equation, stepper, driver and chord finder are created 0118 // ii) this field is used by the equation. 0119 0120 inline void ChangeDetectorField(G4Field* detectorField); 0121 // Pushes the field to the equation ( & keeps its address ) 0122 // Can be used only once the equation, stepper, driver and chord finder 0123 // have all been created. Else it is an error. 0124 0125 inline const G4Field* GetDetectorField() const; 0126 inline G4bool DoesFieldExist() const; 0127 // Set, get and check the field object 0128 0129 void CreateChordFinder(G4MagneticField* detectorMagField); 0130 inline void SetChordFinder(G4ChordFinder* aChordFinder); 0131 inline G4ChordFinder* GetChordFinder(); 0132 inline const G4ChordFinder* GetChordFinder() const; 0133 // Create, set or get the associated Chord Finder 0134 0135 virtual void ConfigureForTrack( const G4Track * ); 0136 // Setup the choice of the configurable parameters 0137 // relying on the current track's energy, particle identity, .. 0138 // Note: in addition to the values of member variables, 0139 // a user can use this to change the ChordFinder, the field, ... 0140 0141 public: // with description 0142 0143 inline G4double GetDeltaIntersection() const; 0144 // Accuracy for boundary intersection. 0145 0146 inline G4double GetDeltaOneStep() const; 0147 // Accuracy for one tracking/physics step. 0148 0149 inline void SetAccuraciesWithDeltaOneStep(G4double valDeltaOneStep); 0150 // Sets both accuracies, maintaining a fixed ratio for accuracies 0151 // of volume Intersection and Integration (in One Step) 0152 0153 inline void SetDeltaOneStep(G4double valueD1step); 0154 // Set accuracy for integration of one step. (only) 0155 inline void SetDeltaIntersection(G4double valueDintersection); 0156 // Set accuracy of intersection of a volume. (only) 0157 0158 inline G4double GetMinimumEpsilonStep() const; 0159 G4bool SetMinimumEpsilonStep( G4double newEpsMin ); 0160 // Minimum for Relative accuracy of a Step 0161 0162 inline G4double GetMaximumEpsilonStep() const; 0163 G4bool SetMaximumEpsilonStep( G4double newEpsMax ); 0164 // Maximum for Relative accuracy of a Step 0165 0166 inline G4bool DoesFieldChangeEnergy() const; 0167 inline void SetFieldChangesEnergy(G4bool value); 0168 // For electric field this should be true 0169 // For magnetic field this should be false 0170 0171 virtual G4FieldManager* Clone() const; 0172 // Needed for multi-threading, create a clone of this object 0173 0174 public: 0175 static G4double GetMaxAcceptedEpsilon(); 0176 static G4bool SetMaxAcceptedEpsilon(G4double maxEps, G4bool softFail= false); 0177 // Set value -- within limits. 0178 // If it fails, with softFail=true it gives Warning, else FatalException 0179 0180 protected: 0181 static G4double fMaxAcceptedEpsilon; 0182 static constexpr G4double fMinAcceptedEpsilon= 1000.0 * std::numeric_limits<G4double>::epsilon(); 0183 // Epsilon_min/max values must be smaller than this - for robust integration 0184 0185 static constexpr G4double fMaxWarningEpsilon= 0.001; // Setting larger value will give warning. 0186 static constexpr G4double fMaxFinalEpsilon= 0.02; // Will not accept larger values 0187 0188 static G4bool fVerboseConstruction; 0189 // Control verbosity of constructors 0190 0191 private: 0192 0193 void InitialiseFieldChangesEnergy(); 0194 // Check whether field/equation change the energy, 0195 // and sets the data member accordingly 0196 // Note: does not handle special cases - this must be done 0197 // separately (e.g. magnetic monopole in B field ) 0198 0199 protected: 0200 void ReportBadEpsilonValue(G4ExceptionDescription& erm, G4double value, 0201 G4String& name) const; 0202 0203 private: 0204 0205 G4Field* fDetectorField = nullptr; 0206 G4ChordFinder* fChordFinder = nullptr; 0207 // Dependent objects -- with state that depends on tracking 0208 0209 G4bool fAllocatedChordFinder = false; // Did we used "new" to 0210 // create fChordFinder ? 0211 // INVARIANTS of tracking --------------------------------------- 0212 // 0213 // 1. 'CONSTANTS' - default values for accuracy parameters 0214 // 0215 const G4double fEpsilonMinDefault= 5.0e-5; // Expected: 5.0e-5 to 1.0e-10 ... 0216 const G4double fEpsilonMaxDefault= 1.0e-3; // Expected: 1.0e-3 to 1.0e-8 ... 0217 0218 static G4double fDefault_Delta_One_Step_Value; // = 0.01 * millimeter; 0219 static G4double fDefault_Delta_Intersection_Val; // = 0.001 * millimeter; 0220 // Default values for accuracy parameters 0221 0222 // 2. CHARACTERISTIC of field 0223 // 0224 G4bool fFieldChangesEnergy = false; 0225 0226 // 3. PARAMETERS that determine the accuracy of integration or intersection 0227 // 0228 G4double fDelta_One_Step_Value; // for one tracking/physics step 0229 G4double fDelta_Intersection_Val; // for boundary intersection 0230 // Values for the required accuracies 0231 0232 G4double fEpsilonMin; 0233 G4double fEpsilonMax; 0234 // Values for the small possible relative accuracy of a step 0235 // (corresponding to the greatest possible integration accuracy) 0236 }; 0237 0238 // Implementation of inline functions 0239 0240 #include "G4FieldManager.icc" 0241 0242 #endif
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |