BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSFieldEMInterpolated3D.cc
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#include "BDSFieldEMInterpolated3D.hh"
20#include "BDSInterpolator3D.hh"
21
22#include "G4ThreeVector.hh"
23#include "G4Types.hh"
24
25#include <utility>
26
27BDSFieldEMInterpolated3D::BDSFieldEMInterpolated3D(BDSInterpolator3D* eInterpolatorIn,
28 BDSInterpolator3D* bInterpolatorIn,
29 const G4Transform3D& offset,
30 G4double eScalingIn,
31 G4double bScalingIn):
32 BDSFieldEMInterpolated(eInterpolatorIn, bInterpolatorIn, offset, eScalingIn, bScalingIn),
33 eInterpolator(eInterpolatorIn),
34 bInterpolator(bInterpolatorIn),
35 eFirstDimensionIndex((eInterpolatorIn->FirstDimension()).underlying()),
36 eFirstTime((eInterpolatorIn->FirstDimension()).underlying() > 2),
37 eSecondDimensionIndex((eInterpolatorIn->SecondDimension()).underlying()),
38 eSecondTime((eInterpolatorIn->SecondDimension()).underlying() > 2),
39 eThirdDimensionIndex((eInterpolatorIn->ThirdDimension()).underlying()),
40 eThirdTime((eInterpolatorIn->ThirdDimension()).underlying() > 2),
41 bFirstDimensionIndex((bInterpolatorIn->FirstDimension()).underlying()),
42 bFirstTime((bInterpolatorIn->FirstDimension()).underlying() > 2),
43 bSecondDimensionIndex((bInterpolatorIn->SecondDimension()).underlying()),
44 bSecondTime((bInterpolatorIn->SecondDimension()).underlying() > 2),
45 bThirdDimensionIndex((bInterpolatorIn->ThirdDimension()).underlying()),
46 bThirdTime((bInterpolatorIn->ThirdDimension()).underlying() > 2)
47{;}
48
49BDSFieldEMInterpolated3D::~BDSFieldEMInterpolated3D()
50{
51 delete eInterpolator;
52 delete bInterpolator;
53}
54
55std::pair<G4ThreeVector,G4ThreeVector> BDSFieldEMInterpolated3D::GetField(const G4ThreeVector& position,
56 const G4double t) const
57{
58 G4double eFCoordinate = 0;
59 if (eFirstTime)
60 {eFCoordinate = t;}
61 else
62 {eFCoordinate = position[eFirstDimensionIndex];}
63 G4double eSCoordinate = 0;
64 if (eSecondTime)
65 {eSCoordinate = t;}
66 else
67 {eSCoordinate = position[eSecondDimensionIndex];}
68 G4double eTCoordinate = 0;
69 if (eThirdTime)
70 {eTCoordinate = t;}
71 else
72 {eTCoordinate = position[eThirdDimensionIndex];}
73 G4double bFCoordinate = 0;
74 if (bFirstTime)
75 {bFCoordinate = t;}
76 else
77 {bFCoordinate = position[bFirstDimensionIndex];}
78 G4double bSCoordinate = 0;
79 if (bSecondTime)
80 {bSCoordinate = t;}
81 else
82 {bSCoordinate = position[bSecondDimensionIndex];}
83 G4double bTCoordinate = 0;
84 if (bThirdTime)
85 {bTCoordinate = t;}
86 else
87 {bTCoordinate = position[bThirdDimensionIndex];}
88
89 G4ThreeVector e = eInterpolator->GetInterpolatedValue(eFCoordinate,
90 eSCoordinate,
91 eTCoordinate) * EScaling();
92 G4ThreeVector b = bInterpolator->GetInterpolatedValue(bFCoordinate,
93 bSCoordinate,
94 bTCoordinate) * BScaling();
95 return std::make_pair(b,e);
96}
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.
G4double EScaling() const
Accessor.
G4double BScaling() const
Accessor.
Interface for all 3D interpolators.
G4ThreeVector GetInterpolatedValue(G4double x, G4double y, G4double z) const
Public interface to a 3D interpolator. Returns Geant4 type as that's what will be needed.