BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSFieldMagOctupole.cc
1/*
2Beam Delivery Simulation (BDSIM) Copyright (C) Royal Holloway,
3University of London 2001 - 2023.
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#include "BDSDebug.hh"
20#include "BDSFieldMagOctupole.hh"
21#include "BDSMagnetStrength.hh"
22#include "BDSUtilities.hh"
23
24#include "G4ThreeVector.hh"
25
26#include "CLHEP/Units/SystemOfUnits.h"
27
28#include <cmath>
29
31 G4double const brho)
32{
33 // B''' = d^3By/dx^3 = Brho * (1/Brho d^3By/dx^3) = Brho * k3
34 bTriplePrime = brho * (*strength)["k3"] / (CLHEP::m3*CLHEP::m);
37#ifdef BDSDEBUG
38 G4cout << __METHOD_NAME__ << "B''' = " << bTriplePrime << G4endl;
39#endif
40}
41
42G4ThreeVector BDSFieldMagOctupole::GetField(const G4ThreeVector &position,
43 const G4double /*t*/) const
44{
45 // B_x = (3x^2y - y^3) * (B'''/3!)
46 // B_y = (x^3-3xy^2) * (B'''/3!)
47 // B_z = 0
48
49 //shortcuts
50 G4double x = position.x();
51 G4double y = position.y();
52
53 G4ThreeVector localField;
54 localField[0] = (3 * std::pow(x,2) * y - std::pow(y,3)) * bTPNormed;
55 localField[1] = (std::pow(x,3) - 3* x * std::pow(y,2)) * bTPNormed;
56 localField[2] = 0;
57
58 return localField;
59}
BDSFieldMagOctupole()
Private default constructor to force use of supplied constructor.
G4double bTPNormed
bTriplePrime / 3! cached for simplicity
G4double bTriplePrime
3rd derivative of the magnetic field
virtual G4ThreeVector GetField(const G4ThreeVector &position, const G4double t=0) const
Access the field value.
G4bool finiteStrength
Flag to cache whether finite nor not.
Definition: BDSFieldMag.hh:83
Efficient storage of magnet strengths.
G4bool IsFinite(G4double value, G4double tolerance=std::numeric_limits< double >::epsilon())