BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSFieldMagQuadrupole.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 "BDSFieldMagQuadrupole.hh"
21#include "BDSMagnetStrength.hh"
22#include "BDSUtilities.hh"
23
24#include "globals.hh" // geant4 types / globals
25#include "G4ThreeVector.hh"
26
27#include "CLHEP/Units/SystemOfUnits.h"
28
30 G4double const brho)
31{
32 // B' = dBy/dx = Brho * (1/Brho dBy/dx) = Brho * k1
33 bPrime = brho * (*strength)["k1"] / CLHEP::m2;
35#ifdef BDSDEBUG
36 G4cout << __METHOD_NAME__ << "B' = " << bPrime << G4endl;
37#endif
38}
39
40G4ThreeVector BDSFieldMagQuadrupole::GetField(const G4ThreeVector &position,
41 const G4double /*t*/) const
42{
43 G4ThreeVector field;
44 field[0] = position.y() * bPrime; // B_x = B' * y;
45 field[1] = position.x() * bPrime; // B_y = B' * x;
46 field[2] = 0; // B_z = 0
47
48 return field;
49}
BDSFieldMagQuadrupole()
Private default constructor to force use of supplied constructor.
virtual G4ThreeVector GetField(const G4ThreeVector &position, const G4double t=0) const
Access the field value.
G4double bPrime
B' - the field gradient - a constant for a quadrupole.
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())