BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSFieldEMInterpolated3D.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 BDSFIELDEMINTERPOLATED3D_H
20#define BDSFIELDEMINTERPOLATED3D_H
21
22#include "BDSFieldEMInterpolated.hh"
23
24#include "G4ThreeVector.hh"
25#include "G4Transform3D.hh"
26#include "G4Types.hh"
27
28#include <utility>
29
31
44{
45public:
46 BDSFieldEMInterpolated3D() = delete;
48 BDSInterpolator3D* bInterpolatorIn,
49 const G4Transform3D& offset = G4Transform3D::Identity,
50 G4double eScalingIn = 1.0,
51 G4double bScalingIn = 1.0);
52
54
56 virtual std::pair<G4ThreeVector,G4ThreeVector> GetField(const G4ThreeVector& position,
57 const G4double t = 0) const;
58
59 inline const BDSInterpolator3D* EInterpolator() const {return eInterpolator;}
60 inline const BDSInterpolator3D* BInterpolator() const {return bInterpolator;}
61
62private:
66 const G4bool eFirstTime;
68 const G4bool eSecondTime;
70 const G4bool eThirdTime;
72 const G4bool bFirstTime;
74 const G4bool bSecondTime;
76 const G4bool bThirdTime;
77};
78
79#endif
A 3D field from an interpolated array with any interpolation.
const G4int bFirstDimensionIndex
Integer index to dimension to use.
const G4bool eFirstTime
Cache of whether to use time coordinate.
const G4bool eSecondTime
Cache of whether to use time coordinate.
const G4int eFirstDimensionIndex
Integer index to dimension to use.
const G4bool eThirdTime
Cache of whether to use time coordinate.
const G4bool bThirdTime
Cache of whether to use time coordinate.
BDSInterpolator3D * eInterpolator
E Interplator the field is based on.
const G4bool bFirstTime
Cache of whether to use time coordinate.
const G4bool bSecondTime
Cache of whether to use time coordinate.
const G4int eThirdDimensionIndex
Integer index to dimension to use.
virtual std::pair< G4ThreeVector, G4ThreeVector > GetField(const G4ThreeVector &position, const G4double t=0) const
Return the interpolated field value at a given point.
const G4int eSecondDimensionIndex
Integer index to dimension to use.
BDSInterpolator3D * bInterpolator
B Interpolator the field is based on.
const G4int bThirdDimensionIndex
Integer index to dimension to use.
const G4int bSecondDimensionIndex
Integer index to dimension to use.
Class to provide scaling and a base class pointer for interpolator fields.
Interface for all 3D interpolators.