BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSFieldEM.hh
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#ifndef BDSFIELDEM_H
20#define BDSFIELDEM_H
21
22#include "globals.hh" // geant4 types / globals
23#include "G4ElectroMagneticField.hh"
24#include "G4ThreeVector.hh"
25#include "G4Transform3D.hh"
26
27#include <utility>
28
29class BDSModulator;
30
40class BDSFieldEM: public G4ElectroMagneticField
41{
42public:
47 BDSFieldEM();
48 explicit BDSFieldEM(G4Transform3D transformIn);
49 virtual ~BDSFieldEM(){;}
50
53 virtual std::pair<G4ThreeVector,G4ThreeVector> GetField(const G4ThreeVector& position,
54 const G4double t = 0) const = 0;
55
58 virtual G4bool TimeVarying() const {return false;}
59
64 virtual void GetFieldValue(const G4double point[4],
65 G4double* field) const;
66
68 virtual std::pair<G4ThreeVector,G4ThreeVector> GetFieldTransformed(const G4ThreeVector& position,
69 const G4double t) const;
70
76 virtual void SetTransform(const G4Transform3D& transformIn);
77
79 virtual G4bool DoesFieldChangeEnergy() const {return true;}
80
82 void SetModulator(BDSModulator* modulatorIn) {modulator = modulatorIn;}
83
85 inline G4bool FiniteStrength() const {return finiteStrength;}
86
87protected:
89
91 G4Transform3D transform;
92
95
96private:
98 G4Transform3D inverseTransform;
99};
100
101#endif
Interface for BDSIM electro-magnetic fields that may or may not be local.
Definition: BDSFieldEM.hh:41
void SetModulator(BDSModulator *modulatorIn)
Set the optional modulator.
Definition: BDSFieldEM.hh:82
virtual G4bool TimeVarying() const
Definition: BDSFieldEM.hh:58
G4Transform3D inverseTransform
The complimentary transform used to initially rotate the point of query.
Definition: BDSFieldEM.hh:98
virtual void SetTransform(const G4Transform3D &transformIn)
Definition: BDSFieldEM.cc:88
G4bool transformIsNotIdentity
Cache of whether to use transform at all.
Definition: BDSFieldEM.hh:93
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:91
BDSModulator * modulator
Optional modulator;.
Definition: BDSFieldEM.hh:94
virtual void GetFieldValue(const G4double point[4], G4double *field) const
Definition: BDSFieldEM.cc:76
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:44
G4bool FiniteStrength() const
Accessor.
Definition: BDSFieldEM.hh:85
virtual G4bool DoesFieldChangeEnergy() const
Required overload by Geant4.
Definition: BDSFieldEM.hh:79
G4bool finiteStrength
Flag to cache whether finite nor not.
Definition: BDSFieldEM.hh:88
Base class for a modulator.
Definition: BDSModulator.hh:37