BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSFieldMagGradient.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 BDSFIELDMAGGRADIENT_H
20#define BDSFIELDMAGGRADIENT_H
21
22#include "globals.hh"
23
24#include <vector>
25
26class BDSFieldMag;
28class BDSFieldMagSkew;
29
43{
44public:
47
49 BDSMagnetStrength* CalculateMultipoles(BDSFieldMag* BField, G4int order, G4double Brho=4);
50
52 G4double GetBy(BDSFieldMag* BField, G4double x, G4double y=0) const;
53
57 std::vector<G4double> PrepareValues(BDSFieldMag* field,
58 G4int order,
59 G4double centreX,
60 G4double h,
61 G4int& centreIndex) const;
62
66 std::vector<std::vector<G4double> > PrepareSkewValues(BDSFieldMag* field,
67 G4int order,
68 G4double centreX,
69 G4double h,
70 G4int& centreIndex) const;
71
75 G4double Derivative(const std::vector<G4double>& data,
76 const G4int order,
77 const G4int startIndex,
78 const G4double h) const;
79
80};
81
82#endif
Find multipole components of fieldmaps by numerical differentiation.
G4double Derivative(const std::vector< G4double > &data, const G4int order, const G4int startIndex, const G4double h) const
std::vector< G4double > PrepareValues(BDSFieldMag *field, G4int order, G4double centreX, G4double h, G4int &centreIndex) const
BDSMagnetStrength * CalculateMultipoles(BDSFieldMag *BField, G4int order, G4double Brho=4)
Find the Magnet strengths of a multipole field to nth order.
std::vector< std::vector< G4double > > PrepareSkewValues(BDSFieldMag *field, G4int order, G4double centreX, G4double h, G4int &centreIndex) const
G4double GetBy(BDSFieldMag *BField, G4double x, G4double y=0) const
Query a point on x axis (by default) for By.
A wrapper class for BDSFieldMag that rotates it.
Interface for static magnetic fields that may or may not be local.
Definition: BDSFieldMag.hh:37
Efficient storage of magnet strengths.