BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSIntegratorMultipoleThin.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 BDSINTEGRATORMULTIPOLETHIN_H
20#define BDSINTEGRATORMULTIPOLETHIN_H
21
22#include "BDSIntegratorMag.hh"
23
24#include "globals.hh"
25#include <list>
26#include <vector>
27
28class G4Mag_EqRhs;
30
38{
39public:
41 G4double brhoIn,
42 G4Mag_EqRhs* eqOfMIn);
43
45
49 virtual void Stepper(const G4double yIn[],
50 const G4double dydx[],
51 const G4double h,
52 G4double yOut[],
53 G4double yErr[]);
54
58 void OneStep(const G4ThreeVector& posIn,
59 const G4ThreeVector& momUIn, // assumed unit momentum of momIn
60 const G4double momIn,
61 G4ThreeVector& posOut,
62 G4ThreeVector& momOut,
63 G4double h) const;
64
65private:
68
70 G4int Factorial(G4int n);
71
73 G4double brho;
74
76 G4double b0l;
78 std::list<double> bnl;
79 std::list<double> bsl;
80 std::vector<G4int> nfact;
82};
83
84#endif
Common functionality to BDSIM integrators.
Integrator that ignores the field and uses the analytical solution to a multipole.
std::vector< G4int > nfact
Higher order components.
G4int Factorial(G4int n)
Calculate the factorial of n.
void OneStep(const G4ThreeVector &posIn, const G4ThreeVector &momUIn, const G4double momIn, G4ThreeVector &posOut, G4ThreeVector &momOut, G4double h) const
BDSIntegratorMultipoleThin()=delete
Private default constructor to enforce use of supplied constructor.
std::list< double > bsl
Higher order components.
G4double brho
Magnetic rigidity for momentum scaling.
virtual void Stepper(const G4double yIn[], const G4double dydx[], const G4double h, G4double yOut[], G4double yErr[])
std::list< double > bnl
Higher order components.
Efficient storage of magnet strengths.