BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSInterpolator4D.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 BDSINTERPOLATOR4D_H
20#define BDSINTERPOLATOR4D_H
21#include "BDSArray4DCoords.hh"
22#include "BDSExtent.hh"
23#include "BDSFieldValue.hh"
24#include "BDSInterpolator.hh"
25
26#include "G4Types.hh"
27#include "G4ThreeVector.hh"
28
30
40{
41public:
42 BDSInterpolator4D() = delete;
43 explicit BDSInterpolator4D(BDSArray4DCoords* arrayIn);
44 virtual ~BDSInterpolator4D(){;}
45
47 G4ThreeVector GetInterpolatedValue(G4double x, G4double y, G4double z, G4double t) const;
48
49 inline const BDSArray4DCoords* Array() const {return array;}
50
52 virtual BDSExtent Extent() const {return array->Extent();}
53
54protected:
59 G4double y,
60 G4double z,
61 G4double t) const = 0;
62
65};
66
67#endif
Overlay of 4D array that provides uniform only spatial coordinate mapping.
virtual BDSExtent Extent() const
Return the SPATIAL (only) extent of this field without any offset. Ignores time.
Holder for +- extents in 3 dimensions.
Definition: BDSExtent.hh:39
Interface for all 4D interpolators.
virtual BDSExtent Extent() const
Extent of field.
virtual BDSFieldValue GetInterpolatedValueT(G4double x, G4double y, G4double z, G4double t) const =0
BDSArray4DCoords * array
The field data.
G4ThreeVector GetInterpolatedValue(G4double x, G4double y, G4double z, G4double t) const
Public interface to a 4D interpolator. Returns Geant4 type as that's what will be needed.
Interface for all interpolators containing basic extent of validity.