BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSFieldMag.cc
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#include "BDSFieldMag.hh"
20
21#include "globals.hh"
22#include "G4ThreeVector.hh"
23#include "G4Transform3D.hh"
24
25#include <cmath>
26
28 finiteStrength(true),
29 transform(G4Transform3D::Identity),
30 inverseTransform(G4Transform3D::Identity)
31{;}
32
33BDSFieldMag::BDSFieldMag(G4Transform3D transformIn):
34 finiteStrength(true),
35 transform(transformIn),
36 inverseTransform(transformIn.inverse())
37{;}
38
39G4ThreeVector BDSFieldMag::GetFieldTransformed(const G4ThreeVector& position,
40 const G4double t) const
41{
42 if (!finiteStrength)
43 {return G4ThreeVector();} // quicker than query
44 else if (transform != G4Transform3D::Identity)
45 {
46 G4ThreeVector transformedPosition = inverseTransform * (HepGeom::Point3D<G4double>)position;
47 G4ThreeVector field = GetField(transformedPosition, t);
48 G4ThreeVector transformedField = transform * (HepGeom::Vector3D<G4double>)field;
49 return transformedField;
50 }
51 else
52 {return GetField(position, t);}
53}
54
55void BDSFieldMag::GetFieldValue(const G4double point[4],
56 G4double* field) const
57{
58 G4double t = point[3];
59 if (std::isnan(t))
60 {t = 0;}
61 G4ThreeVector fieldValue = GetFieldTransformed(G4ThreeVector(point[0], point[1], point[2]), t);
62 field[0] = fieldValue[0]; // B_x
63 field[1] = fieldValue[1]; // B_y
64 field[2] = fieldValue[2]; // B_z
65 field[3] = 0; // E_x
66 field[4] = 0; // E_y
67 field[5] = 0; // E_z
68}
69
70void BDSFieldMag::SetTransform(const G4Transform3D& transformIn)
71{
72 transform = transformIn;
73 inverseTransform = transformIn.inverse();
74}
virtual void SetTransform(const G4Transform3D &transformIn)
Definition: BDSFieldMag.cc:70
virtual void GetFieldValue(const G4double point[4], G4double *field) const
Definition: BDSFieldMag.cc:55
virtual G4ThreeVector GetFieldTransformed(const G4ThreeVector &position, const G4double t) const
Get the field value after applying transform for local offset.
Definition: BDSFieldMag.cc:39
virtual G4ThreeVector GetField(const G4ThreeVector &position, const G4double t=0) const =0
G4bool finiteStrength
Flag to cache whether finite nor not.
Definition: BDSFieldMag.hh:74
G4Transform3D inverseTransform
The complimentary transform used to initially rotate the point of query.
Definition: BDSFieldMag.hh:81
G4Transform3D transform
Transform to apply for the field relative to the local coordinates of the geometry.
Definition: BDSFieldMag.hh:77