BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSFieldE.hh
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#ifndef BDSFIELDE_H
20#define BDSFIELDE_H
21
22#include "globals.hh" // geant4 types / globals
23#include "G4ElectricField.hh"
24#include "G4ThreeVector.hh"
25#include "G4Transform3D.hh"
26
36class BDSFieldE: public G4ElectricField
37{
38public:
43 BDSFieldE();
44 explicit BDSFieldE(G4Transform3D transformIn);
45 virtual ~BDSFieldE(){;}
46
49 virtual G4ThreeVector GetField(const G4ThreeVector& position,
50 const G4double t = 0) const = 0;
51
56 virtual void GetFieldValue(const G4double point[4],
57 G4double* field) const;
58
60 virtual G4ThreeVector GetFieldTransformed(const G4ThreeVector& position,
61 const G4double t) const;
62
68 virtual void SetTransform(const G4Transform3D& transformIn);
69
71 inline G4bool FiniteStrength() const {return finiteStrength;}
72
73protected:
75
77 G4Transform3D transform;
78
79private:
80
82 G4Transform3D inverseTransform;
83};
84
85#endif
Interface for BDSIM electric fields that may or may not be local.
Definition: BDSFieldE.hh:37
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
G4bool FiniteStrength() const
Accessor.
Definition: BDSFieldE.hh:71
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