Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4FieldManager.hh
이 파일의 문서화 페이지로 가기
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 //
27 // $Id: G4FieldManager.hh 104525 2017-06-02 07:22:58Z gcosmo $
28 //
29 //
30 // class G4FieldManager
31 //
32 // Class description:
33 //
34 // A class to manage (Store) a pointer to the Field subclass that
35 // describes the field of a detector (magnetic, electric or other).
36 // Also stores a reference to the chord finder.
37 //
38 // The G4FieldManager class exists to allow the user program to specify
39 // the electric, magnetic and/or other field(s) of the detector.
40 //
41 // A field manager can be set to a logical volume (or to more than one),
42 // in order to vary its field from that of the world. In this manner
43 // a zero or constant field can override a global field, a more or
44 // less exact version can override the external approximation, lower
45 // or higher precision for tracking can be specified, a different
46 // stepper can be chosen for different volumes, ...
47 //
48 // It also stores a pointer to the ChordFinder object that can do the
49 // propagation in this field. All geometrical track "advancement"
50 // in the field is handled by this ChordFinder object.
51 //
52 // G4FieldManager allows the other classes/object (of the MagneticField
53 // & other class categories) to find out whether a detector field object
54 // exists and what that object is.
55 //
56 // The Chord Finder must be created either by calling CreateChordFinder
57 // for a Magnetic Field or by the user creating a a Chord Finder object
58 // "manually" and setting this pointer.
59 //
60 // A default FieldManager is created by the singleton class
61 // G4NavigatorForTracking and exists before main is called.
62 // However a new one can be created and given to G4NavigatorForTracking.
63 //
64 // Our current design envisions that one Field manager is
65 // valid for each region detector.
66 
67 // History:
68 // - 09.06.15 John Apostolakis, Fix to push G4FieldManager* to equation
69 // - 05.11.03 John Apostolakis, Added Min/MaximumEpsilonStep
70 // - 20.06.03 John Apostolakis, Abstract & ability to ConfigureForTrack
71 // - 10.03.97 John Apostolakis, design and implementation.
72 // -------------------------------------------------------------------
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 
85 {
86  public: // with description
87  G4FieldManager(G4Field *detectorField=0,
88  G4ChordFinder *pChordFinder=0,
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  virtual ~G4FieldManager();
96 
97  G4bool SetDetectorField(G4Field *detectorField, int failMode= 0); // =1 is for Debugging ## Was =0
98  // Pushes the field to the equation.
99  // ( New behaviour June 2015 - to avoid the simplest user confusion. )
100  // Failure to push the field ( due to absence of a chord finder, driver,
101  // stepper or equation ) is
102  // - '0' = quiet : Do not complain if chordFinder == 0
103  // (It will still warn for other error.)
104  // - '1' = warn : a warning if anything is missing
105  // - '2'/else = FATAL : a fatal error for all other values.
106  // Returns success (true) or failure (false)
107 
108  inline void ProposeDetectorField(G4Field *detectorField);
109  // Pushes the field to this class only -- no further.
110  // Should be used to initialise this field, only *before* creating
111  // the chord finder and its dependent classes.
112  // User is then responsible to ensure that:
113  // i) an equation, stepper, driver and chord finder are created
114  // ii) this field is used by the equation.
115 
116  inline void ChangeDetectorField(G4Field *detectorField);
117  // Pushes the field to the equation ( & keeps its address )
118  // Can be used only once the equation, stepper, driver and chord finder
119  // have all been created. Else it is an error.
120 
121  inline const G4Field* GetDetectorField() const;
122  inline G4bool DoesFieldExist() const;
123  // Set, get and check the field object
124 
125  void CreateChordFinder(G4MagneticField *detectorMagField);
126  inline void SetChordFinder(G4ChordFinder *aChordFinder);
127  inline G4ChordFinder* GetChordFinder();
128  inline const G4ChordFinder* GetChordFinder() const;
129  // Create, set or get the associated Chord Finder
130 
131  virtual void ConfigureForTrack( const G4Track * );
132  // Setup the choice of the configurable parameters
133  // relying on the current track's energy, particle identity, ..
134  // Note: In addition to the values of member variables,
135  // a user can use this to change the ChordFinder, the field, ...
136 
137  public: // with description
138 
139  inline G4double GetDeltaIntersection() const; // virtual ?
140  // Accuracy for boundary intersection.
141 
142  inline G4double GetDeltaOneStep() const; // virtual ?
143  // Accuracy for one tracking/physics step.
144 
145  inline void SetAccuraciesWithDeltaOneStep(G4double valDeltaOneStep);
146  // Sets both accuracies, maintaining a fixed ratio for accuracties
147  // of volume Intersection and Integration (in One Step)
148 
149  inline void SetDeltaOneStep(G4double valueD1step);
150  // Set accuracy for integration of one step. (only)
151  inline void SetDeltaIntersection(G4double valueDintersection);
152  // Set accuracy of intersection of a volume. (only)
153 
154  inline G4double GetMinimumEpsilonStep() const;
155  inline void SetMinimumEpsilonStep( G4double newEpsMin );
156  // Minimum for Relative accuracy of a Step
157 
158  inline G4double GetMaximumEpsilonStep() const;
159  inline void SetMaximumEpsilonStep( G4double newEpsMax );
160  // Maximum for Relative accuracy of a Step
161 
162  inline G4bool DoesFieldChangeEnergy() const;
163  inline void SetFieldChangesEnergy(G4bool value);
164  // For electric field this should be true
165  // For magnetic field this should be false
166 
167  virtual G4FieldManager* Clone() const;
168  //Needed for multi-threading, create a clone of this object
169 
170  private:
171 
174  // Private copy constructor and assignment operator.
175 
177  // Check whether field/equation change the energy,
178  // and sets the data member accordingly
179  // Note: does not handle special cases - this must be done
180  // separately (e.g. magnetic monopole in B field )
181 
182  private:
183  // Dependent objects -- with state that depends on tracking
186 
187  G4bool fAllocatedChordFinder; // Did we used "new" to
188  // create fChordFinder ?
189  // INVARIANTS of tracking ---------------------------------------
190  //
191  // 1. CONSTANTS
192  const G4double fEpsilonMinDefault; // Can be 1.0e-5 to 1.0e-10 ...
193  const G4double fEpsilonMaxDefault; // Can be 1.0e-3 to 1.0e-8 ...
194 
195  // 2. CHARACTERISTIC of field
197 
198  // 3. PARAMETERS
199  //
200  // Values for the required accuracies
201  G4double fDelta_One_Step_Value; // for one tracking/physics step
202  G4double fDelta_Intersection_Val; // for boundary intersection
203 
206 
207  // Values for the small possible relative accuracy of a step
208  // (corresponding to the greatest possible integration accuracy)
211 
212 };
213 
214 // Our current design and implementation expect that a particular
215 // geometrical region has a Field manager.
216 // By default a Field Manager is created for the world volume, and
217 // will be utilised for all volumes unless it is overridden by a 'local'
218 // field manager.
219 
220 // Note also that a region with both electric E and magnetic B field will
221 // have these treated as one field.
222 // Similarly it could be extended to treat other fields as additional components
223 // of a single field type.
224 
225 
226 // Implementation of inline functions
227 
228 #include "G4FieldManager.icc"
229 
230 #endif /* G4FIELDMANAGER_HH */
virtual void ConfigureForTrack(const G4Track *)
void SetMinimumEpsilonStep(G4double newEpsMin)
G4bool DoesFieldChangeEnergy() const
G4double fDefault_Delta_One_Step_Value
G4double fEpsilonMin
G4double GetDeltaOneStep() const
G4FieldManager & operator=(const G4FieldManager &)
G4bool DoesFieldExist() const
G4double GetDeltaIntersection() const
virtual ~G4FieldManager()
G4FieldManager(G4Field *detectorField=0, G4ChordFinder *pChordFinder=0, G4bool b=true)
void SetAccuraciesWithDeltaOneStep(G4double valDeltaOneStep)
G4bool SetDetectorField(G4Field *detectorField, int failMode=0)
void SetFieldChangesEnergy(G4bool value)
void CreateChordFinder(G4MagneticField *detectorMagField)
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
const XML_Char int const XML_Char * value
Definition: expat.h:331
const G4double fEpsilonMaxDefault
G4double GetMaximumEpsilonStep() const
G4double fEpsilonMax
G4double fDelta_One_Step_Value
virtual G4FieldManager * Clone() const
void InitialiseFieldChangesEnergy()
const G4Field * GetDetectorField() const
void ChangeDetectorField(G4Field *detectorField)
G4double fDefault_Delta_Intersection_Val
G4ChordFinder * fChordFinder
G4double fDelta_Intersection_Val
void SetChordFinder(G4ChordFinder *aChordFinder)
G4Field * fDetectorField
const G4double fEpsilonMinDefault
void ProposeDetectorField(G4Field *detectorField)
void SetDeltaIntersection(G4double valueDintersection)
G4bool fAllocatedChordFinder
G4double GetMinimumEpsilonStep() const
G4bool fFieldChangesEnergy
void SetMaximumEpsilonStep(G4double newEpsMax)
void SetDeltaOneStep(G4double valueD1step)
G4ChordFinder * GetChordFinder()