Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
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//
27// $Id$
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// - 05.11.03 John Apostolakis, Added Min/MaximumEpsilonStep
69// - 20.06.03 John Apostolakis, Abstract & ability to ConfigureForTrack
70// - 10.03.97 John Apostolakis, design and implementation.
71// -------------------------------------------------------------------
72
73#ifndef G4FIELDMANAGER_HH
74#define G4FIELDMANAGER_HH 1
75
76#include "globals.hh"
77
78class G4Field;
79class G4MagneticField;
80class G4ChordFinder;
81class G4Track; // Forward reference for parameter configuration
82
84{
85 public: // with description
86 G4FieldManager(G4Field *detectorField=0,
87 G4ChordFinder *pChordFinder=0,
88 G4bool b=true ); // fieldChangesEnergy is taken from field
89 // General constructor for any field.
90 // -> Must be set with field and chordfinder for use.
91 G4FieldManager(G4MagneticField *detectorMagneticField);
92 // Creates ChordFinder
93 // - assumes pure magnetic field (so Energy constant)
94 virtual ~G4FieldManager();
95
96 G4bool SetDetectorField(G4Field *detectorField);
97 inline const G4Field* GetDetectorField() const;
98 inline G4bool DoesFieldExist() const;
99 // Set, get and check the field object
100
101 void CreateChordFinder(G4MagneticField *detectorMagField);
102 inline void SetChordFinder(G4ChordFinder *aChordFinder);
104 inline const G4ChordFinder* GetChordFinder() const;
105 // Create, set or get the associated Chord Finder
106
107 virtual void ConfigureForTrack( const G4Track * );
108 // Setup the choice of the configurable parameters
109 // relying on the current track's energy, particle identity, ..
110 // Note: In addition to the values of member variables,
111 // a user can use this to change the ChordFinder, the field, ...
112
113 public: // with description
114
115 inline G4double GetDeltaIntersection() const; // virtual ?
116 // Accuracy for boundary intersection.
117
118 inline G4double GetDeltaOneStep() const; // virtual ?
119 // Accuracy for one tracking/physics step.
120
121 inline void SetAccuraciesWithDeltaOneStep(G4double valDeltaOneStep);
122 // Sets both accuracies, maintaining a fixed ratio for accuracties
123 // of volume Intersection and Integration (in One Step)
124
125 inline void SetDeltaOneStep(G4double valueD1step);
126 // Set accuracy for integration of one step. (only)
127 inline void SetDeltaIntersection(G4double valueDintersection);
128 // Set accuracy of intersection of a volume. (only)
129
131 inline void SetMinimumEpsilonStep( G4double newEpsMin );
132 // Minimum for Relative accuracy of a Step
133
135 inline void SetMaximumEpsilonStep( G4double newEpsMax );
136 // Maximum for Relative accuracy of a Step
137
139 inline void SetFieldChangesEnergy(G4bool value);
140 // For electric field this should be true
141 // For magnetic field this should be false
142
143 private:
144
146 G4FieldManager& operator=(const G4FieldManager&);
147 // Private copy constructor and assignment operator.
148
149 private:
150 // Dependent objects -- with state that depends on tracking
151 G4Field* fDetectorField;
152 G4ChordFinder* fChordFinder;
153
154 G4bool fAllocatedChordFinder; // Did we used "new" to
155 // create fChordFinder ?
156 // INVARIANTS of tracking ---------------------------------------
157 //
158 // 1. CONSTANTS
159 const G4double fEpsilonMinDefault; // Can be 1.0e-5 to 1.0e-10 ...
160 const G4double fEpsilonMaxDefault; // Can be 1.0e-3 to 1.0e-8 ...
161
162 // 2. CHARACTERISTIC of field
163 G4bool fFieldChangesEnergy;
164
165 // 3. PARAMETERS
166 //
167 // Values for the required accuracies
168 G4double fDelta_One_Step_Value; // for one tracking/physics step
169 G4double fDelta_Intersection_Val; // for boundary intersection
170
171 G4double fDefault_Delta_One_Step_Value; // = 0.25 * mm;
172 G4double fDefault_Delta_Intersection_Val; // = 0.1 * mm;
173
174 // Values for the small possible relative accuracy of a step
175 // (corresponding to the greatest possible integration accuracy)
176 G4double fEpsilonMin;
177 G4double fEpsilonMax;
178
179};
180
181// Our current design and implementation expect that a particular
182// geometrical region has a Field manager.
183// By default a Field Manager is created for the world volume, and
184// will be utilised for all volumes unless it is overridden by a 'local'
185// field manager.
186
187// Note also that a region with both electric E and magnetic B field will
188// have these treated as one field.
189// Similarly it could be extended to treat other fields as additional components
190// of a single field type.
191
192
193// Implementation of inline functions
194
195#include "G4FieldManager.icc"
196
197#endif /* G4FIELDMANAGER_HH */
double G4double
Definition: G4Types.hh:64
bool G4bool
Definition: G4Types.hh:67
void SetAccuraciesWithDeltaOneStep(G4double valDeltaOneStep)
G4bool SetDetectorField(G4Field *detectorField)
G4bool DoesFieldChangeEnergy() const
virtual ~G4FieldManager()
void CreateChordFinder(G4MagneticField *detectorMagField)
void SetFieldChangesEnergy(G4bool value)
void SetDeltaOneStep(G4double valueD1step)
void SetChordFinder(G4ChordFinder *aChordFinder)
virtual void ConfigureForTrack(const G4Track *)
G4double GetMinimumEpsilonStep() const
const G4ChordFinder * GetChordFinder() const
void SetMinimumEpsilonStep(G4double newEpsMin)
G4double GetMaximumEpsilonStep() const
G4double GetDeltaOneStep() const
const G4Field * GetDetectorField() const
void SetDeltaIntersection(G4double valueDintersection)
G4ChordFinder * GetChordFinder()
void SetMaximumEpsilonStep(G4double newEpsMax)
G4double GetDeltaIntersection() const
G4bool DoesFieldExist() const