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