BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSFieldInfo.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 BDSFIELDINFO_H
20#define BDSFIELDINFO_H
21
22#include "BDSArrayReflectionType.hh"
23#include "BDSFieldFormat.hh"
24#include "BDSFieldType.hh"
25#include "BDSIntegratorType.hh"
26#include "BDSInterpolatorType.hh"
27
28#include "globals.hh" // geant4 types / globals
29#include "G4Transform3D.hh"
30#include "G4ThreeVector.hh"
31
32#include <ostream>
33
36class G4UserLimits;
37
66{
67public:
70 BDSFieldInfo(BDSFieldType fieldTypeIn,
71 G4double brhoIn,
72 BDSIntegratorType integratorTypeIn,
73 BDSMagnetStrength* magnetStrengthIn = nullptr,
74 G4bool provideGlobalTransformIn = true,
75 const G4Transform3D& transformIn = G4Transform3D(),
76 const G4String& magneticFieldFilePathIn = "",
77 BDSFieldFormat magneticFieldFormatIn = BDSFieldFormat::bdsim1d,
78 BDSInterpolatorType magneticInterpolatorTypeIn = BDSInterpolatorType::nearest3d,
79 const G4String& electricFieldFilePathIn = "",
80 BDSFieldFormat electricFieldFormatIn = BDSFieldFormat::bdsim1d,
81 BDSInterpolatorType electricInterpolatorTypeIn = BDSInterpolatorType::nearest3d,
82 G4bool cacheTransformsIn = false,
83 G4double eScalingIn = 1.0,
84 G4double bScalingIn = 1.0,
85 G4double timeOffsetIn = 0,
86 G4bool autoScaleIn = false,
87 G4UserLimits* stepLimitIn = nullptr,
88 G4double poleTipRadiusIn = 1,
89 G4double beamPipeRadiusIn = 0,
90 G4bool left = true,
91 const G4String& magneticSubFieldNameIn = "",
92 const G4String& electricSubFieldNameIn = "");
94
96 BDSFieldInfo(const BDSFieldInfo& other);
99
101 inline BDSFieldType FieldType() const {return fieldType;}
102 inline G4double BRho() const {return brho;}
103 inline BDSIntegratorType IntegratorType() const {return integratorType;}
104 inline BDSMagnetStrength* MagnetStrength() const {return magnetStrength;}
105 inline G4bool ProvideGlobal() const {return provideGlobalTransform;}
106 inline G4String MagneticFile() const {return magneticFieldFilePath;}
107 inline BDSFieldFormat MagneticFormat() const {return magneticFieldFormat;}
108 inline BDSInterpolatorType MagneticInterpolatorType() const {return magneticInterpolatorType;}
109 inline G4String ElectricFile() const {return electricFieldFilePath;}
110 inline BDSFieldFormat ElectricFormat() const {return electricFieldFormat;}
111 inline BDSInterpolatorType ElectricInterpolatorType() const {return electricInterpolatorType;}
112 inline G4bool CacheTransforms() const {return cacheTransforms;}
113 inline G4double EScaling() const {return eScaling;}
114 inline G4double BScaling() const {return bScaling;}
115 inline G4double TimeOffset() const {return timeOffset;}
116 inline G4bool AutoScale() const {return autoScale;}
117 inline G4UserLimits* UserLimits() const {return stepLimit;}
118 inline G4double PoleTipRadius() const {return poleTipRadius;}
119 inline G4double BeamPipeRadius() const {return beamPipeRadius;}
120 inline G4double ChordStepMinimum() const {return chordStepMinimum;}
121 inline G4double Tilt() const {return tilt;}
122 inline G4bool SecondFieldOnLeft() const {return secondFieldOnLeft;}
123 inline G4String MagneticSubFieldName() const {return magneticSubFieldName;}
124 inline G4String ElectricSubFieldName() const {return electricSubFieldName;}
125 inline G4String NameOfParserDefinition() const {return nameOfParserDefinition;}
126 inline G4bool UsePlacementWorldTransform() const {return usePlacementWorldTransform;}
127 inline const BDSArrayReflectionTypeSet& MagneticArrayReflectionType() const {return magneticArrayReflectionTypeSet;}
128 inline const BDSArrayReflectionTypeSet& ElectricArrayReflectionType() const {return electricArrayReflectionTypeSet;}
129 inline BDSModulatorInfo* ModulatorInfo() const {return modulatorInfo;}
131
132 G4Transform3D Transform() const;
133 G4Transform3D TransformBeamline() const;
134 G4Transform3D TransformComplete() const;
135
137 inline void SetFieldType(BDSFieldType fieldTypeIn) {fieldType = fieldTypeIn;}
138 inline void SetIntegratorType(BDSIntegratorType typeIn) {integratorType = typeIn;}
139 inline void SetProvideGlobalTransform(G4bool provideGlobalTransformIn) {provideGlobalTransform = provideGlobalTransformIn;}
140 inline void SetMagneticInterpolatorType(BDSInterpolatorType typeIn) {magneticInterpolatorType = typeIn;}
141 inline void SetMagneticArrayReflectionType(const BDSArrayReflectionTypeSet& typeIn) {magneticArrayReflectionTypeSet = typeIn;}
142 inline void SetElectricArrayReflectionType(const BDSArrayReflectionTypeSet& typeIn) {electricArrayReflectionTypeSet = typeIn;}
143 inline void SetBScaling(G4double bScalingIn) {bScaling = bScalingIn;}
144 inline void SetAutoScale(G4bool autoScaleIn) {autoScale = autoScaleIn;}
145 inline void SetScalingRadius(G4double poleTipRadiusIn) {poleTipRadius = poleTipRadiusIn;}
146 inline void SetBeamPipeRadius(G4double beamPipeRadiusIn) {beamPipeRadius = beamPipeRadiusIn;}
147 inline void SetChordStepMinimum(G4double chordStepMinimumIn) {chordStepMinimum = chordStepMinimumIn;}
148 inline void SetSecondFieldOnLeft(G4bool leftIn) {secondFieldOnLeft = leftIn;}
149 inline void SetMagneticSubField(const G4String& mfnIn) {magneticSubFieldName = mfnIn;}
150 inline void SetElectricSubField(const G4String& efnIn) {electricSubFieldName = efnIn;}
151 inline void SetUsePlacementWorldTransform(G4bool use) {usePlacementWorldTransform = use;}
152 inline void SetModulatorInfo(BDSModulatorInfo* modulatorInfoIn) {modulatorInfo = modulatorInfoIn;}
153
155 inline void CompoundBScaling(G4double extraBScalingIn) {bScaling *= extraBScalingIn;}
157 inline void CompoundEScaling(G4double extraEScalingIn) {eScaling *= extraEScalingIn;}
158
159 void SetTransform(const G4Transform3D& transformIn);
160 void SetTransformBeamline(const G4Transform3D& transformIn);
161
163 void SetUserLimits(G4UserLimits* userLimitsIn);
164
165 void SetNameOfParserDefinition(const G4String& nameIn) {nameOfParserDefinition = nameIn;}
166
169 void UpdateUserLimitsLengthMaximumStepSize(G4double maximumStepSize,
170 G4bool warn = false) const;
171
174 void Translate(const G4ThreeVector& translationIn);
175
177 inline void CacheTransforms(G4bool cacheTransformsIn) {cacheTransforms = cacheTransformsIn;}
178
180 friend std::ostream& operator<< (std::ostream &out, BDSFieldInfo const &info);
181
182 static G4UserLimits* defaultUL;
183
184private:
185 BDSFieldType fieldType;
186 G4double brho;
187 BDSIntegratorType integratorType;
188 BDSMagnetStrength* magnetStrength;
189 G4bool provideGlobalTransform;
190 G4Transform3D* transform;
191 G4String magneticFieldFilePath;
192 BDSFieldFormat magneticFieldFormat;
193 BDSInterpolatorType magneticInterpolatorType;
194 BDSArrayReflectionTypeSet magneticArrayReflectionTypeSet;
195 G4String electricFieldFilePath;
196 BDSFieldFormat electricFieldFormat;
197 BDSInterpolatorType electricInterpolatorType;
198 BDSArrayReflectionTypeSet electricArrayReflectionTypeSet;
199 G4bool cacheTransforms;
200 G4double eScaling;
201 G4double bScaling;
202 G4double timeOffset;
203 G4bool autoScale;
204 mutable G4UserLimits* stepLimit;
205 G4double poleTipRadius;
206 G4double beamPipeRadius;
207 G4double chordStepMinimum;
208 G4double tilt;
210 G4String magneticSubFieldName;
211 G4String electricSubFieldName;
212 G4bool usePlacementWorldTransform;
213 BDSModulatorInfo* modulatorInfo;
214
216 G4Transform3D* transformBeamline;
217
218 G4String nameOfParserDefinition;
219
220 // We need a default to pass back if none is specified.
221 const static G4ThreeVector defaultUnitDirection;
222};
223
224#endif
All info required to build complete field of any type.
Definition: BDSFieldInfo.hh:66
void CacheTransforms(G4bool cacheTransformsIn)
Turn on or off transform caching.
BDSFieldInfo & operator=(const BDSFieldInfo &)=delete
Assignment operator not used.
G4double beamPipeRadius
Optional radius of beam pipe.
G4String MagneticSubFieldName() const
Accessor.
void CompoundEScaling(G4double extraEScalingIn)
*= for EScaling.
void UpdateUserLimitsLengthMaximumStepSize(G4double maximumStepSize, G4bool warn=false) const
static G4UserLimits * defaultUL
Cache of default user limits.
G4Transform3D * transformBeamline
Transform from curvilinear frame to this field - ie beam line bit only.
G4double BScaling() const
Accessor.
const BDSArrayReflectionTypeSet & ElectricArrayReflectionType() const
Accessor.
G4double TimeOffset() const
Accessor.
BDSFieldFormat MagneticFormat() const
Accessor.
G4double ChordStepMinimum() const
Accessor.
G4Transform3D TransformComplete() const
Compound transform of field + beam line transform.
G4String NameOfParserDefinition() const
Accessor.
G4String ElectricFile() const
Accessor.
void SetTransformBeamline(const G4Transform3D &transformIn)
Set the beam line transform.
BDSInterpolatorType ElectricInterpolatorType() const
Accessor.
void CompoundBScaling(G4double extraBScalingIn)
*= for BScaling.
G4bool SecondFieldOnLeft() const
Accessor.
BDSInterpolatorType MagneticInterpolatorType() const
Accessor.
G4double poleTipRadius
Radius at which point the field will be scaled to.
G4double EScaling() const
Accessor.
G4double PoleTipRadius() const
Accessor.
G4double tilt
Cache of tilt of field.
G4bool ProvideGlobal() const
Accessor.
void SetFieldType(BDSFieldType fieldTypeIn)
Set Transform - could be done afterwards once instance of this class is passed around.
G4String MagneticFile() const
Accessor.
G4bool CacheTransforms() const
Accessor.
friend std::ostream & operator<<(std::ostream &out, BDSFieldInfo const &info)
output stream
BDSIntegratorType IntegratorType() const
Accessor.
void SetTransform(const G4Transform3D &transformIn)
Set the field definition transform.
G4double BeamPipeRadius() const
Accessor.
G4String ElectricSubFieldName() const
Accessor.
G4Transform3D Transform() const
Transform for the field definition only.
G4bool UsePlacementWorldTransform() const
Accessor.
BDSFieldFormat ElectricFormat() const
Accessor.
const BDSArrayReflectionTypeSet & MagneticArrayReflectionType() const
Accessor.
G4UserLimits * UserLimits() const
Accessor.
G4bool secondFieldOnLeft
Flag for case of two-beam field - if not left, it's right.
void Translate(const G4ThreeVector &translationIn)
BDSModulatorInfo * ModulatorInfo() const
Accessor.
G4double Tilt() const
Accessor.
BDSFieldType FieldType() const
Accessor.
BDSMagnetStrength * MagnetStrength() const
Accessor.
G4double BRho() const
Accessor.
G4bool AutoScale() const
Accessor.
G4Transform3D TransformBeamline() const
Transform from the curvilinear coordinates to the beam line component.
void SetUserLimits(G4UserLimits *userLimitsIn)
Delete and replace the user limits which this class owns (only if not default ul).
G4Transform3D * transform
Transform w.r.t. solid field will be attached to.
BDSFieldInfo()
Default constructor for zero field effectively.
Definition: BDSFieldInfo.cc:40
Efficient storage of magnet strengths.
Holder class for all information required to describe a modulator.