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