source: trunk/source/geometry/magneticfield/include/G4FieldManager.hh@ 1215

Last change on this file since 1215 was 921, checked in by garnier, 17 years ago

en test de gl2ps. Problemes de libraries

File size: 8.1 KB
RevLine 
[831]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,v 1.16 2006/06/29 18:22:15 gunter Exp $
[921]28// GEANT4 tag $Name: geant4-09-02-cand-01 $
[831]29//
30//
31// class G4FieldManager
32//
33// Class description:
34//
35// A class to manage (Store) a pointer to the Field subclass that
36// describes the field of a detector (magnetic, electric or other).
37// Also stores a reference to the chord finder.
38//
39// The G4FieldManager class exists to allow the user program to specify
40// the electric, magnetic and/or other field(s) of the detector.
41//
42// A field manager can be set to a logical volume (or to more than one),
43// in order to vary its field from that of the world. In this manner
44// a zero or constant field can override a global field, a more or
45// less exact version can override the external approximation, lower
46// or higher precision for tracking can be specified, a different
47// stepper can be chosen for different volumes, ...
48//
49// It also stores a pointer to the ChordFinder object that can do the
50// propagation in this field. All geometrical track "advancement"
51// in the field is handled by this ChordFinder object.
52//
53// G4FieldManager allows the other classes/object (of the MagneticField
54// & other class categories) to find out whether a detector field object
55// exists and what that object is.
56//
57// The Chord Finder must be created either by calling CreateChordFinder
58// for a Magnetic Field or by the user creating a a Chord Finder object
59// "manually" and setting this pointer.
60//
61// A default FieldManager is created by the singleton class
62// G4NavigatorForTracking and exists before main is called.
63// However a new one can be created and given to G4NavigatorForTracking.
64//
65// Our current design envisions that one Field manager is
66// valid for each region detector.
67
68// History:
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
79class G4Field;
80class G4MagneticField;
81class G4ChordFinder;
82class G4Track; // Forward reference for parameter configuration
83
84class G4FieldManager
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);
98 inline const G4Field* GetDetectorField() const;
99 inline G4bool DoesFieldExist() const;
100 // Set, get and check the field object
101
102 void CreateChordFinder(G4MagneticField *detectorMagField);
103 inline void SetChordFinder(G4ChordFinder *aChordFinder);
104 inline G4ChordFinder* GetChordFinder();
105 inline const G4ChordFinder* GetChordFinder() const;
106 // Create, set or get the associated Chord Finder
107
108 virtual void ConfigureForTrack( const G4Track * );
109 // Setup the choice of the configurable parameters
110 // relying on the current track's energy, particle identity, ..
111 // Note: In addition to the values of member variables,
112 // a user can use this to change the ChordFinder, the field, ...
113
114 public: // with description
115
116 inline G4double GetDeltaIntersection() const; // virtual ?
117 // Accuracy for boundary intersection.
118
119 inline G4double GetDeltaOneStep() const; // virtual ?
120 // Accuracy for one tracking/physics step.
121
122 inline void SetAccuraciesWithDeltaOneStep(G4double valDeltaOneStep);
123 // Sets both accuracies, maintaining a fixed ratio for accuracties
124 // of volume Intersection and Integration (in One Step)
125
126 inline void SetDeltaOneStep(G4double valueD1step);
127 // Set accuracy for integration of one step. (only)
128 inline void SetDeltaIntersection(G4double valueDintersection);
129 // Set accuracy of intersection of a volume. (only)
130
131 inline G4double GetMinimumEpsilonStep() const;
132 inline void SetMinimumEpsilonStep( G4double newEpsMin );
133 // Minimum for Relative accuracy of a Step
134
135 inline G4double GetMaximumEpsilonStep() const;
136 inline void SetMaximumEpsilonStep( G4double newEpsMax );
137 // Maximum for Relative accuracy of a Step
138
139 inline G4bool DoesFieldChangeEnergy() const;
140 inline void SetFieldChangesEnergy(G4bool value);
141 // For electric field this should be true
142 // For magnetic field this should be false
143
144 private:
145
146 G4FieldManager(const G4FieldManager&);
147 G4FieldManager& operator=(const G4FieldManager&);
148 // Private copy constructor and assignment operator.
149
150 private:
151
152 G4Field* fDetectorField;
153 G4ChordFinder* fChordFinder;
154
155 G4bool fAllocatedChordFinder; // Did we used "new" to
156 // create fChordFinder ?
157 G4bool fFieldChangesEnergy;
158
159 // Values for the required accuracies
160 //
161 G4double fDelta_One_Step_Value; // for one tracking/physics step
162 G4double fDelta_Intersection_Val; // for boundary intersection
163
164 G4double fDefault_Delta_One_Step_Value; // = 0.25 * mm;
165 G4double fDefault_Delta_Intersection_Val; // = 0.1 * mm;
166
167 // Values for the small possible relative accuracy of a step
168 // (corresponding to the greatest possible integration accuracy)
169
170 G4double fEpsilonMinDefault; // Can be 1.0e-5 to 1.0e-10 ...
171 G4double fEpsilonMaxDefault; // Can be 1.0e-3 to 1.0e-8 ...
172 G4double fEpsilonMin;
173 G4double fEpsilonMax;
174};
175
176// Our current design and implementation expect that a particular
177// geometrical region has a Field manager.
178// By default a Field Manager is created for the world volume, and
179// will be utilised for all volumes unless it is overridden by a 'local'
180// field manager.
181
182// Note also that a region with both electric E and magnetic B field will
183// have these treated as one field.
184// Similarly it could be extended to treat other fields as additional components
185// of a single field type.
186
187
188// Implementation of inline functions
189
190#include "G4FieldManager.icc"
191
192#endif /* G4FIELDMANAGER_HH */
Note: See TracBrowser for help on using the repository browser.