BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSIntegratorDipoleRodrigues2.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 BDSINTEGRATORDIPOLERODRIGUES2_H
20#define BDSINTEGRATORDIPOLERODRIGUES2_H
21
22#include "BDSIntegratorDrift.hh"
23#include "BDSMagnetStrength.hh"
24
25#include "globals.hh"
26#include "G4MagHelicalStepper.hh"
27#include "G4ThreeVector.hh"
28
29class G4MagIntegratorStepper;
30class G4Mag_EqRhs;
31
53class BDSIntegratorDipoleRodrigues2: public G4MagHelicalStepper, public BDSIntegratorDrift
54{
55public:
56 BDSIntegratorDipoleRodrigues2(G4Mag_EqRhs* eqOfMIn,
58
60
63 virtual void DumbStepper(const G4double yIn[6],
64 G4ThreeVector field,
65 G4double stepLength,
66 G4double yOut[6]);
67
71 void SingleStep(const G4double yIn[6],
72 const G4double& h,
73 G4double yOut[6]);
74
78 virtual void Stepper(const G4double yIn[6],
79 const G4double dydx[],
80 G4double h,
81 G4double yOut[6],
82 G4double yErr[6]);
83
84 virtual G4int IntegratorOrder() const {return 1;}
85
88 void AdvanceHelixForSpiralling(const G4double yIn[6],
89 const G4ThreeVector& field,
90 const G4double& h,
91 G4double yOut[6],
92 G4double yErr[6]);
93
95 inline G4double RadiusOfHelix() const {return GetRadHelix();}
96
97protected:
101
103 G4Mag_EqRhs* eqOfM;
104
105private:
108
112
113};
114
115#endif
Exact helix through pure dipole field.
void AdvanceHelixForSpiralling(const G4double yIn[6], const G4ThreeVector &field, const G4double &h, G4double yOut[6], G4double yErr[6])
G4Mag_EqRhs * eqOfM
Cache of equation of motion. This class does not own it.
void SingleStep(const G4double yIn[6], const G4double &h, G4double yOut[6])
virtual void DumbStepper(const G4double yIn[6], G4ThreeVector field, G4double stepLength, G4double yOut[6])
G4double RadiusOfHelix() const
Public accessor for protected variable in base class.
virtual void Stepper(const G4double yIn[6], const G4double dydx[], G4double h, G4double yOut[6], G4double yErr[6])
Output error array.
BDSIntegratorDipoleRodrigues2()=delete
Private default constructor to force use of provided one.
Routine for drift algorithm.