Loading...
Searching...
No Matches
F04GlobalField.cc
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/// \file field/field04/src/F04GlobalField.cc
28/// \brief Implementation of the F04GlobalField class
29//
30
31#include <ctime>
32
33#include "Randomize.hh"
34#include "G4TransportationManager.hh"
35
36#include "G4ExplicitEuler.hh"
37#include "G4ImplicitEuler.hh"
38#include "G4SimpleRunge.hh"
39#include "G4SimpleHeum.hh"
40#include "G4ClassicalRK4.hh"
41#include "G4CashKarpRKF45.hh"
42#include "G4SystemOfUnits.hh"
43
44#include "F04GlobalField.hh"
45#include "F04SimpleSolenoid.hh"
46#include "F04FocusSolenoid.hh"
47
48//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
49
50G4ThreadLocal F04GlobalField* F04GlobalField::fObject = nullptr;
51
52//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
53
55 : fDetectorConstruction(det)
56{
57 fFieldMessenger = new F04FieldMessenger(this,det);
58
59 fFields = new FieldList();
60
61 // set object
62 fObject = this;
63
65}
66
67//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
68
70{
71 Clear();
72
73 delete fFields;
74
75 delete fFieldMessenger;
76
77 delete fEquation;
78 delete fStepper;
79 delete fChordFinder;
80}
81
82//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
83
85{
86 Clear();
87
88 // Construct equ. of motion of particles through B fields
89// fEquation = new G4Mag_EqRhs(this);
90 // Construct equ. of motion of particles through e.m. fields
91// fEquation = new G4EqMagElectricField(this);
92 // Construct equ. of motion of particles including spin through B fields
93// fEquation = new G4Mag_SpinEqRhs(this);
94 // Construct equ. of motion of particles including spin through e.m. fields
96
97 // Get transportation, field, and propagator managers
98 G4TransportationManager* transportManager =
99 G4TransportationManager::GetTransportationManager();
100
102
103 fFieldPropagator = transportManager->GetPropagatorInField();
104
105 // Need to SetFieldChangesEnergy to account for a time varying electric
106 // field (r.f. fields)
107 fFieldManager->SetFieldChangesEnergy(true);
108
109 // Set the field
110 fFieldManager->SetDetectorField(this);
111
112 // Choose a stepper for integration of the equation of motion
113 SetStepper();
114
115 // Create a cord finder providing the (global field, min step length,
116 // a pointer to the stepper)
118
119 // Set accuracy parameters
120 fChordFinder->SetDeltaChord( fDeltaChord );
121
122 fFieldManager->SetAccuraciesWithDeltaOneStep(fDeltaOneStep);
123
124 fFieldManager->SetDeltaIntersection(fDeltaIntersection);
125
126 fFieldPropagator->SetMinimumEpsilonStep(fEpsMin);
127 fFieldPropagator->SetMaximumEpsilonStep(fEpsMax);
128
129 G4cout << "Accuracy Parameters:" <<
130 " MinStep=" << fMinStep <<
131 " DeltaChord=" << fDeltaChord <<
132 " DeltaOneStep=" << fDeltaOneStep << G4endl;
133 G4cout << " " <<
134 " DeltaIntersection=" << fDeltaIntersection <<
135 " EpsMin=" << fEpsMin <<
136 " EpsMax=" << fEpsMax << G4endl;
137
138 fFieldManager->SetChordFinder(fChordFinder);
139
140 G4double l = 0.0;
143
145 G4ThreeVector captureMgntCenter =
147
148 auto focusSolenoid =
149 new F04FocusSolenoid(B1, B2, l, logicCaptureMgnt,captureMgntCenter);
150 focusSolenoid -> SetHalf(true);
151
153
154 G4LogicalVolume* logicTransferMgnt =
156 G4ThreeVector transferMgntCenter =
158
159 auto simpleSolenoid =
160 new F04SimpleSolenoid(B, l, logicTransferMgnt,transferMgntCenter);
161
162 simpleSolenoid->SetColor("1,0,1");
163 simpleSolenoid->SetColor("0,1,1");
164 simpleSolenoid->SetMaxStep(1.5*mm);
165 simpleSolenoid->SetMaxStep(2.5*mm);
166
167 if (fFields) {
168 if (fFields->size()>0) {
169 FieldList::iterator i;
170 for (i=fFields->begin(); i!=fFields->end(); ++i){
171 (*i)->Construct();
172 }
173 }
174 }
175}
176
177//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
178
184
185//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
186
188{
189 if (fObject) return fObject;
190 return nullptr;
191}
192
193//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
194
196{
197 delete fStepper;
198
199 switch ( fStepperType )
200 {
201 case 0:
202// fStepper = new G4ExplicitEuler( fEquation, 8 ); // no spin tracking
203 fStepper = new G4ExplicitEuler( fEquation, 12 ); // with spin tracking
204 G4cout << "G4ExplicitEuler is called" << G4endl;
205 break;
206 case 1:
207// fStepper = new G4ImplicitEuler( fEquation, 8 ); // no spin tracking
208 fStepper = new G4ImplicitEuler( fEquation, 12 ); // with spin tracking
209 G4cout << "G4ImplicitEuler is called" << G4endl;
210 break;
211 case 2:
212// fStepper = new G4SimpleRunge( fEquation, 8 ); // no spin tracking
213 fStepper = new G4SimpleRunge( fEquation, 12 ); // with spin tracking
214 G4cout << "G4SimpleRunge is called" << G4endl;
215 break;
216 case 3:
217// fStepper = new G4SimpleHeum( fEquation, 8 ); // no spin tracking
218 fStepper = new G4SimpleHeum( fEquation, 12 ); // with spin tracking
219 G4cout << "G4SimpleHeum is called" << G4endl;
220 break;
221 case 4:
222// fStepper = new G4ClassicalRK4( fEquation, 8 ); // no spin tracking
223 fStepper = new G4ClassicalRK4( fEquation, 12 ); // with spin tracking
224 G4cout << "G4ClassicalRK4 (default) is called" << G4endl;
225 break;
226 case 5:
227// fStepper = new G4CashKarpRKF45( fEquation, 8 ); // no spin tracking
228 fStepper = new G4CashKarpRKF45( fEquation, 12 ); // with spin tracking
229 G4cout << "G4CashKarpRKF45 is called" << G4endl;
230 break;
231 default: fStepper = nullptr;
232 }
233}
234
235//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
236
238{
239 return G4TransportationManager::GetTransportationManager()
240 ->GetFieldManager();
241}
242
243//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
244
245void F04GlobalField::GetFieldValue(const G4double* point, G4double* field) const
246{
247 // NOTE: this routine dominates the CPU time for tracking.
248 // Using the simple array fFp[] instead of fields[]
249 // directly sped it up
250
251 field[0] = field[1] = field[2] = field[3] = field[4] = field[5] = 0.0;
252
253 // protect against Geant4 bug that calls us with point[] NaN.
254 if(point[0] != point[0]) return;
255
256 // (can't use fNfp or fFp, as they may change)
257 if (fFirst) ((F04GlobalField*)this)->SetupArray(); // (cast away const)
258
259 for (int i=0; i<fNfp; ++i) {
260 const F04ElementField* p = fFp[i];
261 if (p->IsInBoundingBox(point)) {
262 p->AddFieldValue(point,field);
263 }
264 }
265
266}
267
268//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
269
271{
272 if (fFields) {
273 if (fFields->size()>0) {
274 FieldList::iterator i;
275 for (i=fFields->begin(); i!=fFields->end(); ++i) delete *i;
276 fFields->clear();
277 }
278 }
279
280 delete [] fFp;
281 fFirst = true;
282 fNfp = 0;
283 fFp = nullptr;
284}
285
286//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
287
289{
290 fFirst = false;
291 fNfp = fFields->size();
292 fFp = new const F04ElementField* [fNfp+1]; // add 1 so it's never 0
293 for (int i=0; i<fNfp; ++i) fFp[i] = (*fFields)[i];
294}
Definition of the F04FocusSolenoid class.
Definition of the F04GlobalField class.
std::vector< F04ElementField * > FieldList
Definition of the F04SimpleSolenoid class.
bool IsInBoundingBox(const G4double point[4]) const
IsInBoundingBox() returns true if the point is within the global bounding box - global coordinates.
virtual void AddFieldValue(const G4double point[4], G4double field[6]) const =0
AddFieldValue() will add the field value for this element to field[].
static G4ThreadLocal F04GlobalField * fObject
void ConstructField()
constructs all field tracking objects
~F04GlobalField() override
F04DetectorConstruction * fDetectorConstruction
G4FieldManager * fFieldManager
const F04ElementField ** fFp
G4ChordFinder * fChordFinder
G4PropagatorInField * fFieldPropagator
F04FieldMessenger * fFieldMessenger
void Clear()
Clear() removes all ElementField-s from the global object, and destroys them.
void GetFieldValue(const G4double *point, G4double *field) const override
GetFieldValue() returns the field value at a given point[].
G4EqEMFieldWithSpin * fEquation
FieldList * fFields
F04GlobalField(F04DetectorConstruction *const)
void SetStepper()
Set the Stepper.
static F04GlobalField * GetObject()
G4double fDeltaIntersection
G4FieldManager * GetGlobalFieldManager()
Get the global field manager.
G4MagIntegratorStepper * fStepper

Applications | User Support | Publications | Collaboration