BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSIntegratorDipoleRodrigues.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 BDSINTEGRATORDIPOLERODRIGUES_H
20#define BDSINTEGRATORDIPOLERODRIGUES_H
21
22#include "BDSIntegratorMag.hh"
23
24#include "globals.hh"
25#include "G4ThreeVector.hh"
26
27#include <utility>
28
29class G4Mag_EqRhs;
31
44{
45public:
47 G4double brho,
48 G4Mag_EqRhs* eqOfMIn);
49
51
54 virtual void Stepper(const G4double yIn[],
55 const G4double dydx[],
56 const G4double h,
57 G4double yOut[],
58 G4double yErr[]);
59
60protected:
62 void AdvanceHelix(const G4double yIn[],
63 const G4double dydx[],
64 G4double h,
65 G4double yOut[],
66 G4double yErr[]);
67
69 G4double cOverGeV;
70
72 const G4double angle;
73
75 const G4double length;
76
78 G4double bField;
79
82
83 std::pair<G4ThreeVector,G4ThreeVector> UpdatePandR(G4double rho,
84 G4double h,
85 G4ThreeVector localPos,
86 G4ThreeVector localMomUnit);
87
88private:
90 G4ThreeVector yInitial, yMidPoint, yFinal;
91
94};
95
96#endif
Stepper that calculates trajectory through uniform magnetic field.
void AdvanceHelix(const G4double yIn[], const G4double dydx[], G4double h, G4double yOut[], G4double yErr[])
Calculate the new particle coordinates.
const G4double angle
Angle that the dipole induces in the reference trajectory.
virtual void Stepper(const G4double yIn[], const G4double dydx[], const G4double h, G4double yOut[], G4double yErr[])
const G4double minimumRadiusOfCurvature
Minimum tolerable radius of curvature - used to prevent spiraling particles.
G4double bField
Uniform magnetic field in global Y direction.
G4double cOverGeV
Scaling factor in brho calculation.
BDSMagnetStrength const * strength
Magnet strength object needed by curvilinear transforms.
G4ThreeVector yInitial
Data stored in order to find the chord.
const G4double length
Arc length of the magnetic field.
Common functionality to BDSIM integrators.
Efficient storage of magnet strengths.