BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSFieldE.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 "BDSFieldE.hh"
20#include "BDSModulator.hh"
21
22#include "globals.hh" // geant4 types / globals
23#include "G4ThreeVector.hh"
24#include "G4Transform3D.hh"
25
27 finiteStrength(true),
28 transform(G4Transform3D::Identity),
29 transformIsNotIdentity(false),
30 modulator(nullptr),
31 inverseTransform(G4Transform3D::Identity)
32{;}
33
34BDSFieldE::BDSFieldE(G4Transform3D transformIn):
35 finiteStrength(true),
36 transform(transformIn),
37 transformIsNotIdentity(transformIn != G4Transform3D::Identity),
38 modulator(nullptr),
39 inverseTransform(transformIn.inverse())
40{;}
41
42G4ThreeVector BDSFieldE::GetFieldTransformed(const G4ThreeVector& position,
43 const G4double t) const
44{
45 if (!finiteStrength)
46 {return G4ThreeVector();} // quicker than query
48 {
49 G4ThreeVector transformedPosition = inverseTransform * (HepGeom::Point3D<G4double>)position;
50 G4ThreeVector field = GetField(transformedPosition, t);
51 if (modulator)
52 {
53 G4double factor = modulator->Factor(transformedPosition, t);
54 field *= factor;
55 }
56 G4ThreeVector transformedField = transform * (HepGeom::Vector3D<G4double>)field;
57 return transformedField;
58 }
59 else
60 {
61 G4ThreeVector field = GetField(position,t);
62 if (modulator)
63 {
64 G4double factor = modulator->Factor(position, t);
65 field *= factor;
66 }
67 return field;
68 }
69}
70
71void BDSFieldE::GetFieldValue(const G4double point[4],
72 G4double* field) const
73{
74 G4ThreeVector fieldValue = GetFieldTransformed(G4ThreeVector(point[0], point[1],
75 point[2]), point[3]);
76 field[0] = 0; // B_x
77 field[1] = 0; // B_y
78 field[2] = 0; // B_z
79 field[3] = fieldValue[0]; // E_x
80 field[4] = fieldValue[1]; // E_y
81 field[5] = fieldValue[2]; // E_z
82}
83
84void BDSFieldE::SetTransform(const G4Transform3D& transformIn)
85{
86 transform = transformIn;
87 transformIsNotIdentity = transformIn != G4Transform3D::Identity;
88 inverseTransform = transformIn.inverse();
89}
G4Transform3D inverseTransform
The complimentary transform used to initially rotate the point of query.
Definition: BDSFieldE.hh:93
BDSModulator * modulator
Optional modulator;.
Definition: BDSFieldE.hh:89
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:42
virtual void GetFieldValue(const G4double point[4], G4double *field) const
Definition: BDSFieldE.cc:71
G4bool finiteStrength
Flag to cache whether finite nor not.
Definition: BDSFieldE.hh:83
G4Transform3D transform
Transform to apply for the field relative to the local coordinates of the geometry.
Definition: BDSFieldE.hh:86
G4bool transformIsNotIdentity
Cache of whether to use transform at all.
Definition: BDSFieldE.hh:88
virtual void SetTransform(const G4Transform3D &transformIn)
Definition: BDSFieldE.cc:84
virtual G4double Factor(const G4ThreeVector &xyz, G4double T) const =0