BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSIntegratorDipoleFringe.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 BDSINTEGRATORDIPOLEFRINGEBASE_H
20#define BDSINTEGRATORDIPOLEFRINGEBASE_H
21
22#include "BDSAuxiliaryNavigator.hh"
23#include "BDSIntegratorDipoleRodrigues2.hh"
24#include "BDSIntegratorDipoleQuadrupole.hh"
25#include "BDSIntegratorMultipoleThin.hh"
26
27#include "globals.hh"
28
29namespace BDS
30{
32 G4double FringeFieldCorrection(BDSMagnetStrength const* strength, G4bool entranceOrExit);
33
35 G4double SecondFringeFieldCorrection(BDSMagnetStrength const* strength, G4bool entranceOrExit);
36}
37
38class G4Mag_EqRhs;
40
48{
49public:
51 G4double brhoIn,
52 G4Mag_EqRhs* eqOfMIn,
54 const G4double& tiltIn = 0);
55
57
62
64 virtual void Stepper(const G4double yIn[6],
65 const G4double dydx[6],
66 const G4double h,
67 G4double yOut[6],
68 G4double yErr[6]);
69
73 void BaseStepper(const G4double yIn[6],
74 const G4double dydx[6],
75 const G4double& h,
76 G4double yOut[6],
77 G4double yErr[6],
78 const G4double& fcof,
79 const G4double& momScaling);
80
81
85 void OneStep(const G4ThreeVector& posIn,
86 const G4ThreeVector& momUIn, // assumed unit momentum of momIn
87 G4ThreeVector& posOut,
88 G4ThreeVector& momOut,
89 const G4double& bendingRadius) const;
90
93 void MultipoleStep(const G4double yIn[6],
94 G4double yMultipoleOut[6],
95 const G4double& h);
96
99 inline G4double GetPolefaceAngle() {return polefaceAngle;}
100 inline G4double GetFringeCorr() {return fringeCorr;}
101 inline G4double GetSecondFringeCorr() {return secondFringeCorr;}
102 inline G4double GetPolefaceCurv() {return polefaceCurvature;}
103
104protected:
108 G4double fringeCorr;
114 G4double rho;
115
116 G4ThreeVector unitField;
117 const G4double fieldArcLength;
118 const G4double fieldAngle;
119
120 const G4double tilt;
121 const G4bool finiteTilt;
122 G4double bx;
123 G4double by;
127
130 static G4double thinElementLength;
131
132 BDSIntegratorMultipoleThin* multipoleIntegrator;
133
134 G4bool isEntrance;
135
139
140private:
143};
144
145#endif
Extra G4Navigator to get coordinate transforms.
Integrator that ignores the field and uses the analytical solution for a dipole kick.
G4double secondFringeCorr
Second fringe field correction term.
void OneStep(const G4ThreeVector &posIn, const G4ThreeVector &momUIn, G4ThreeVector &posOut, G4ThreeVector &momOut, const G4double &bendingRadius) const
G4double rho
Nominal magnet bending radius.
G4double polefaceAngle
Poleface rotation angle.
G4double polefaceCurvature
Poleface curvature.
G4bool isEntrance
store if fringe is entrance or exit
BDSIntegratorDipoleFringe(BDSIntegratorDipoleFringe &)=delete
Assignment and copy constructor not implemented nor used.
G4double fringeCorr
Fringe field correction term.
const G4double fieldArcLength
Cache of the field arc length.
void BaseStepper(const G4double yIn[6], const G4double dydx[6], const G4double &h, G4double yOut[6], G4double yErr[6], const G4double &fcof, const G4double &momScaling)
BDSIntegratorDipoleFringe & operator=(const BDSIntegratorDipoleFringe &)=delete
Assignment and copy constructor not implemented nor used.
virtual void Stepper(const G4double yIn[6], const G4double dydx[6], const G4double h, G4double yOut[6], G4double yErr[6])
The stepper for integration. Calls base class stepper.
BDSIntegratorDipoleFringe()=delete
Private default constructor to enforce use of supplied constructor.
G4ThreeVector unitField
Cache of the unit field direction.
const G4double fieldAngle
Cache of the field angle.
void MultipoleStep(const G4double yIn[6], G4double yMultipoleOut[6], const G4double &h)
Exact helix through pure dipole field.
Integrator that ignores the field and uses the analytical solution to a multipole.
Efficient storage of magnet strengths.
Return either G4Tubs or G4CutTubs depending on flat face.
G4double SecondFringeFieldCorrection(BDSMagnetStrength const *strength, G4bool entranceOrExit)
Function to calculate the value of the second fringe field correction term.
G4double FringeFieldCorrection(BDSMagnetStrength const *strength, G4bool entranceOrExit)
Function to calculate the value of the fringe field correction term.