BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSFieldMagDipoleQuadrupole.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 "BDSFieldMagDipoleQuadrupole.hh"
21#include "BDSFieldMagQuadrupole.hh"
22#include "BDSFieldMagDipole.hh"
23
24#include "globals.hh" // geant4 types / globals
25#include "G4ThreeVector.hh"
26
27
29 G4double const brho):
30 quad(new BDSFieldMagQuadrupole(strength, brho)),
31 dipole(new BDSFieldMagDipole(strength))
32{
33 finiteStrength = quad->FiniteStrength() || dipole->FiniteStrength();
34}
35
36BDSFieldMagDipoleQuadrupole::~BDSFieldMagDipoleQuadrupole()
37{
38 delete quad;
39 delete dipole;
40}
41
42G4ThreeVector BDSFieldMagDipoleQuadrupole::GetField(const G4ThreeVector &position,
43 const G4double t) const
44{
45 G4ThreeVector quadField = quad->GetFieldTransformed(position, t);
46 G4ThreeVector dipoleField = dipole->GetFieldTransformed(position, t);
47 return quadField + dipoleField;
48}
49
50void BDSFieldMagDipoleQuadrupole::SetTransform(const G4Transform3D& transformIn)
51{
52 quad->SetTransform(transformIn);
53 dipole->SetTransform(transformIn);
54}
virtual void SetTransform(const G4Transform3D &transformIn)
virtual G4ThreeVector GetField(const G4ThreeVector &position, const G4double t=0) const
Access the field value.
BDSFieldMagDipoleQuadrupole()
Private default constructor to force use of supplied constructor.
A uniform dipole field.
Class that provides the magnetic strength in a quadrupole.
virtual void SetTransform(const G4Transform3D &transformIn)
Definition: BDSFieldMag.cc:88
virtual G4ThreeVector GetFieldTransformed(const G4ThreeVector &position, const G4double t) const
Get the field value after applying transform for local offset.
Definition: BDSFieldMag.cc:44
Efficient storage of magnet strengths.