BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSArray4DCoords.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 BDSARRAY4DCOORDS_H
20#define BDSARRAY4DCOORDS_H
21
22#include "BDSArray4D.hh"
23#include "BDSDimensionType.hh"
24#include "BDSFourVector.hh"
25
26#include "globals.hh"
27
28#include <array>
29#include <ostream>
30
31class BDSExtent;
32
49{
50public:
53 BDSArray4DCoords() = delete;
54
58 BDSArray4DCoords(G4int nXIn, G4int nYIn, G4int nZIn, G4int nTIn,
59 G4double xMinIn, G4double xMaxIn,
60 G4double yMinIn, G4double yMaxIn,
61 G4double zMinIn, G4double zMaxIn,
62 G4double tMinIn, G4double tMaxIn,
63 BDSDimensionType xDimensionIn = BDSDimensionType::x,
64 BDSDimensionType yDimensionIn = BDSDimensionType::y,
65 BDSDimensionType zDimensionIn = BDSDimensionType::z,
66 BDSDimensionType tDimensionIn = BDSDimensionType::t);
67
68 virtual ~BDSArray4DCoords(){;}
69
71 inline G4double XStep() const {return xStep;}
72 inline G4double YStep() const {return yStep;}
73 inline G4double ZStep() const {return zStep;}
74 inline G4double TStep() const {return tStep;}
76
77 inline G4double SmallestSpatialStep() const {return smallestSpatialStep;}
78
80 virtual G4bool OutsideCoords(G4double x,
81 G4double y,
82 G4double z,
83 G4double t) const;
84
88 virtual void OutsideCoordsWarn(G4double x,
89 G4double y,
90 G4double z,
91 G4double t) const;
92
94 virtual G4double ArrayCoordsFromX(G4double x) const {return (x - xMin) / xStep;}
95 virtual G4double ArrayCoordsFromY(G4double y) const {return (y - yMin) / yStep;}
96 virtual G4double ArrayCoordsFromZ(G4double z) const {return (z - zMin) / zStep;}
97 virtual G4double ArrayCoordsFromT(G4double t) const {return (t - tMin) / tStep;}
99
101 void ArrayCoordsFromXY(G4double& x, G4double& xArr,
102 G4double& y, G4double& yArr) const
103 {xArr = ArrayCoordsFromX(x); yArr = ArrayCoordsFromY(y);}
104 void ArrayCoordsFromXYZ(G4double& x, G4double& xArr,
105 G4double& y, G4double& yArr,
106 G4double& z, G4double& zArr) const
107 {xArr = ArrayCoordsFromX(x); yArr = ArrayCoordsFromY(y); zArr = ArrayCoordsFromZ(z);}
108 void ArrayCoordsFromXYZT(G4double& x, G4double& xArr,
109 G4double& y, G4double& yArr,
110 G4double& z, G4double& zArr,
111 G4double& t, G4double& tArr) const
112 {xArr = ArrayCoordsFromX(x); yArr = ArrayCoordsFromY(y); zArr = ArrayCoordsFromZ(z); tArr = ArrayCoordsFromT(t);}
114
117 G4double y,
118 G4double z,
119 G4double t) const
120 {
125 }
126
128 inline G4double XFromArrayCoords(G4double xCoord) const {return xMin + xCoord*xStep;}
129 inline G4double YFromArrayCoords(G4double yCoord) const {return yMin + yCoord*yStep;}
130 inline G4double ZFromArrayCoords(G4double zCoord) const {return zMin + zCoord*zStep;}
131 inline G4double TFromArrayCoords(G4double tCoord) const {return tMin + tCoord*tStep;}
133
137 G4double y,
138 G4double z,
139 G4double t) const
140 {
145 }
146
148 virtual G4int NearestX(G4double x) const {return (G4int)round((x - xMin) / xStep);}
149 virtual G4int NearestY(G4double y) const {return (G4int)round((y - yMin) / yStep);}
150 virtual G4int NearestZ(G4double z) const {return (G4int)round((z - zMin) / zStep);}
151 virtual G4int NearestT(G4double t) const {return (G4int)round((t - tMin) / tStep);}
153
156 G4double y,
157 G4double z,
158 G4double t) const
159 {
161 NearestY(y),
162 NearestZ(z),
163 NearestT(t));
164 }
165
168 {
169 return BDSFourVector<G4int>(NearestX(pos[0]),
170 NearestY(pos[1]),
171 NearestZ(pos[2]),
172 NearestT(pos[3]));
173 }
174
176 virtual void ExtractSection2x2x2x2(G4double x,
177 G4double y,
178 G4double z,
179 G4double t,
180 BDSFieldValue (&localData)[2][2][2][2],
181 G4double& xFrac,
182 G4double& yFrac,
183 G4double& zFrac,
184 G4double& tFrac) const;
185
187 virtual void ExtractSection4x4x4x4(G4double x,
188 G4double y,
189 G4double z,
190 G4double t,
191 BDSFieldValue (&localData)[4][4][4][4],
192 G4double& xFrac,
193 G4double& yFrac,
194 G4double& zFrac,
195 G4double& tFrac) const;
196
198 virtual BDSFieldValue ExtractNearest(G4double x,
199 G4double y = 0,
200 G4double z = 0,
201 G4double t = 0) const;
202
205 virtual std::ostream& Print(std::ostream& out) const;
206
208 friend std::ostream& operator<< (std::ostream& out, BDSArray4DCoords const &a);
209
211 virtual BDSExtent Extent() const;
212
213 inline G4double XMin() const {return xMin;}
214 inline G4double YMin() const {return yMin;}
215 inline G4double ZMin() const {return zMin;}
216 inline G4double TMin() const {return tMin;}
217 inline G4double XMax() const {return xMax;}
218 inline G4double YMax() const {return yMax;}
219 inline G4double ZMax() const {return zMax;}
220 inline G4double TMax() const {return tMax;}
221
224 inline BDSDimensionType SecondDimension() const {return yDimension;}
225 inline BDSDimensionType ThirdDimension() const {return zDimension;}
226 inline BDSDimensionType FourthDimension() const {return tDimension;}
228
231 G4int DimensionIndex(BDSDimensionType spatialDimension) const;
232
234 inline std::array<G4int, 4> ArrayToSpatialDimensionIndices() const {return dimensions;}
235
236protected:
238 void BuildDimensionIndex();
239
241 G4double xMin;
242 G4double xMax;
243 G4double yMin;
244 G4double yMax;
245 G4double zMin;
246 G4double zMax;
247 G4double tMin;
248 G4double tMax;
249
250 G4double xStep;
251 G4double yStep;
252 G4double zStep;
253 G4double tStep;
255
256 G4double smallestSpatialStep;
257
261 BDSDimensionType yDimension;
262 BDSDimensionType zDimension;
263 BDSDimensionType tDimension;
264 std::array<G4int, 4> dimensions;
265
266private:
267 static void CheckStep(G4double step, const G4String& name);
268};
269
270#endif
Overlay of 4D array that provides uniform only spatial coordinate mapping.
G4double xMin
Dimension parameter - protected for derived class access.
void ArrayCoordsFromXY(G4double &x, G4double &xArr, G4double &y, G4double &yArr) const
Utility version to forward to individual function.
virtual G4bool OutsideCoords(G4double x, G4double y, G4double z, G4double t) const
Whether the spatial coordinates lie outside the range of the array or not.
G4double yMin
Dimension parameter - protected for derived class access.
G4double XStep() const
The distance in spatial coordinates between any two points in the array.
G4double YFromArrayCoords(G4double yCoord) const
Return spatial value from a continuous array coordinate in one dimension.
virtual G4double ArrayCoordsFromX(G4double x) const
Not much point in being both virtual and inline (in our use case) but has to be virtual.
BDSFourVector< G4int > NearestXYZT(const BDSFourVector< G4double > &pos) const
Return the index of the nearest field value in space via four-vector.
G4double yStep
Dimension parameter - protected for derived class access.
BDSDimensionType ThirdDimension() const
Accessor for each dimension label. e.g. array 'x' = spatial z.
virtual std::ostream & Print(std::ostream &out) const
virtual G4int NearestY(G4double y) const
Not much point in being both virtual and inline (in our use case) but has to be virtual.
virtual G4double ArrayCoordsFromZ(G4double z) const
Not much point in being both virtual and inline (in our use case) but has to be virtual.
virtual G4int NearestX(G4double x) const
Not much point in being both virtual and inline (in our use case) but has to be virtual.
BDSDimensionType SecondDimension() const
Accessor for each dimension label. e.g. array 'x' = spatial z.
BDSFourVector< G4double > XYZTFromArrayCoords(G4double x, G4double y, G4double z, G4double t) const
virtual BDSExtent Extent() const
Return the SPATIAL (only) extent of this field without any offset. Ignores time.
virtual G4double ArrayCoordsFromY(G4double y) const
Not much point in being both virtual and inline (in our use case) but has to be virtual.
std::array< G4int, 4 > ArrayToSpatialDimensionIndices() const
Access all indices at once.
virtual G4double ArrayCoordsFromT(G4double t) const
Not much point in being both virtual and inline (in our use case) but has to be virtual.
G4double zMax
Dimension parameter - protected for derived class access.
virtual G4int NearestT(G4double t) const
Not much point in being both virtual and inline (in our use case) but has to be virtual.
G4double tMin
Dimension parameter - protected for derived class access.
friend std::ostream & operator<<(std::ostream &out, BDSArray4DCoords const &a)
Output stream.
BDSDimensionType FourthDimension() const
Accessor for each dimension label. e.g. array 'x' = spatial z.
G4double zMin
Dimension parameter - protected for derived class access.
BDSDimensionType xDimension
virtual G4int NearestZ(G4double z) const
Not much point in being both virtual and inline (in our use case) but has to be virtual.
G4double xStep
Dimension parameter - protected for derived class access.
virtual void OutsideCoordsWarn(G4double x, G4double y, G4double z, G4double t) const
virtual void ExtractSection2x2x2x2(G4double x, G4double y, G4double z, G4double t, BDSFieldValue(&localData)[2][2][2][2], G4double &xFrac, G4double &yFrac, G4double &zFrac, G4double &tFrac) const
Extract 2x2x2x2 points lying around coordinate x.
G4double YStep() const
The distance in spatial coordinates between any two points in the array.
G4int DimensionIndex(BDSDimensionType spatialDimension) const
G4double yMax
Dimension parameter - protected for derived class access.
G4double xMax
Dimension parameter - protected for derived class access.
virtual void ExtractSection4x4x4x4(G4double x, G4double y, G4double z, G4double t, BDSFieldValue(&localData)[4][4][4][4], G4double &xFrac, G4double &yFrac, G4double &zFrac, G4double &tFrac) const
Extract 4x4x4x4 points lying around coordinate x.
G4double TFromArrayCoords(G4double tCoord) const
Return spatial value from a continuous array coordinate in one dimension.
G4double tMax
Dimension parameter - protected for derived class access.
G4double zStep
Dimension parameter - protected for derived class access.
G4double TStep() const
The distance in spatial coordinates between any two points in the array.
void ArrayCoordsFromXYZ(G4double &x, G4double &xArr, G4double &y, G4double &yArr, G4double &z, G4double &zArr) const
Utility version to forward to individual function.
void ArrayCoordsFromXYZT(G4double &x, G4double &xArr, G4double &y, G4double &yArr, G4double &z, G4double &zArr, G4double &t, G4double &tArr) const
Utility version to forward to individual function.
BDSArray4DCoords()=delete
G4double ZFromArrayCoords(G4double zCoord) const
Return spatial value from a continuous array coordinate in one dimension.
G4double XFromArrayCoords(G4double xCoord) const
Return spatial value from a continuous array coordinate in one dimension.
G4double ZStep() const
The distance in spatial coordinates between any two points in the array.
BDSFourVector< G4int > NearestXYZT(G4double x, G4double y, G4double z, G4double t) const
Return the index of the nearest field value in space.
virtual BDSFieldValue ExtractNearest(G4double x, G4double y=0, G4double z=0, G4double t=0) const
Extract nearest field value from array.
BDSFourVector< G4double > ArrayCoordsFromXYZT(G4double x, G4double y, G4double z, G4double t) const
Convenience function to easily get array coords in all dimensions at once.
G4double tStep
Dimension parameter - protected for derived class access.
BDSDimensionType FirstDimension() const
Accessor for each dimension label. e.g. array 'x' = spatial z.
void BuildDimensionIndex()
Build up an array of ints based on the order of dimensions stored in the array.
4D array and base class for 3,2 & 1D arrays.
Definition: BDSArray4D.hh:51
Holder for +- extents in 3 dimensions.
Definition: BDSExtent.hh:39