BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSFieldMagMultipoleOuter.hh
1/*
2Beam Delivery Simulation (BDSIM) Copyright (C) Royal Holloway,
3University of London 2001 - 2022.
4
5This file is part of BDSIM.
6
7BDSIM is free software: you can redistribute it and/or modify
8it under the terms of the GNU General Public License as published
9by the Free Software Foundation version 3 of the License.
10
11BDSIM is distributed in the hope that it will be useful, but
12WITHOUT ANY WARRANTY; without even the implied warranty of
13MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14GNU General Public License for more details.
15
16You should have received a copy of the GNU General Public License
17along with BDSIM. If not, see <http://www.gnu.org/licenses/>.
18*/
19#ifndef BDSFIELDMAGMULTIPOLEOUTER_H
20#define BDSFIELDMAGMULTIPOLEOUTER_H
21
22#include "BDSFieldMag.hh"
23
24#include "globals.hh" // geant4 types / globals
25#include "G4RotationMatrix.hh"
26#include "G4ThreeVector.hh"
27#include "G4TwoVector.hh"
28
29#include <vector>
30
32
56{
57public:
58 BDSFieldMagMultipoleOuter(G4int orderIn, // order must be > 0
59 G4double poleTipRadius,
60 const BDSFieldMag* innerFieldIn,
61 G4bool kPositive,
62 G4double bRho,
63 G4double arbitraryScaling = 1.0);
64
66
68 virtual G4ThreeVector GetField(const G4ThreeVector& position,
69 const double t = 0) const;
70
71private:
72 const G4int order;
73 G4double phiOffset;
74 G4double spatialLimit;
75 G4double normalisation;
78 G4double poleTipRadius;
79 std::vector<G4TwoVector> currents;
80 G4double maxField;
82};
83
84#endif
A simple parameterisation of N-Pole outer yoke magnetic field.
virtual G4ThreeVector GetField(const G4ThreeVector &position, const double t=0) const
Access the field value.
G4double spatialLimit
Radius from any current source within which the field is artificially saturated.
std::vector< G4TwoVector > currents
Locations of infinite wire current sources.
G4double phiOffset
Tilt in XY calculated from B vector of inner field. if B0=(0,1,0), phiOffset=0.
G4int poleNOffset
Offset for pole to start at - in effect this flips the sign of the field.
G4double maxField
Any field beyond this will curtailed to this value.
G4bool initialisationPhase
Need a way to control cludge normalisation behaviour during initial normalisation calculation.
G4double poleTipRadius
Radius of transition between inner and outer fields.
G4double normalisation
Storage of the overall normalisation factor.
G4bool positiveField
Sign of magnetic field.
Interface for static magnetic fields that may or may not be local.
Definition: BDSFieldMag.hh:37
Efficient storage of magnet strengths.