BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSIntegratorQuadrupole.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 BDSINTEGRATORQUADRUPOLE_H
20#define BDSINTEGRATORQUADRUPOLE_H
21
22#include "BDSIntegratorMag.hh"
23
24#include "globals.hh"
25
26class G4Mag_EqRhs;
28class BDSStep;
29
43{
44public:
46 G4double brho,
47 G4Mag_EqRhs* eqOfMIn,
48 G4double minimumRadiusOfCurvatureIn);
49
50 virtual ~BDSIntegratorQuadrupole(){;}
51
55 virtual void Stepper(const G4double y[],
56 const G4double dydx[],
57 const G4double h,
58 G4double yOut[],
59 G4double yErr[]);
60
61private:
64
66 G4double bPrime;
67
72};
73
74#endif
Common functionality to BDSIM integrators.
Integrator that ignores the field and uses the analytical solution to a quadrupole.
virtual void Stepper(const G4double y[], const G4double dydx[], const G4double h, G4double yOut[], G4double yErr[])
BDSIntegratorQuadrupole()
Private default constructor to enforce use of supplied constructor.
G4double bPrime
B Field Gradient.
Efficient storage of magnet strengths.
A simple class to represent the positions of a step.
Definition: BDSStep.hh:33