BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSFieldQueryInfo.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 BDSFIELDQUERYINFO_H
20#define BDSFIELDQUERYINFO_H
21
22#include "BDSFourVector.hh"
23
24#include "G4AffineTransform.hh"
25#include "G4String.hh"
26#include "G4Types.hh"
27
28#include <vector>
29
37{
38public:
40 {
43 QueryDimensionInfo(G4int nIn, G4double minIn, G4double maxIn):
44 n(nIn), min(minIn), max(maxIn)
45 {;}
46 QueryDimensionInfo(): n(1), min(0), max(0) {;}
47 G4int n = 1;
48 G4double min = 0;
49 G4double max = 0;
50 G4double StepSize() const {return (max - min)/n;}
51 };
52
54 BDSFieldQueryInfo(const G4String& nameIn,
55 const G4String& outfileMagneticIn,
56 const G4String& outfileElectricIn,
57 G4bool queryMagneticIn,
58 G4bool queryElectricIn,
59 QueryDimensionInfo xInfoIn,
60 QueryDimensionInfo yInfoIn,
61 QueryDimensionInfo zInfoIn,
62 QueryDimensionInfo tInfoIn,
63 const G4AffineTransform& globalTransformIn = G4AffineTransform(),
64 G4bool overwriteExistingFilesIn = false,
65 const G4String& fieldObjectIn = "",
66 G4bool printTransformIn = false,
67 G4bool checkParametersIn = true,
68 G4bool drawArrowsIn = true,
69 G4bool drawZeroValuePointsIn = true,
70 G4bool drawBoxesIn = true,
71 G4double boxAlphaIn = 0.2);
72
74 BDSFieldQueryInfo(const G4String& nameIn,
75 const G4String& outfileMagneticIn,
76 const G4String& outfileElectricIn,
77 G4bool queryMagneticIn,
78 G4bool queryElectricIn,
79 const std::vector<BDSFourVector<G4double>>& pointsToQueryIn,
80 const std::vector<G4String>& pointsColumnNamesIn,
81 G4bool overwriteExistingFilesIn = false,
82 const G4String& fieldObjectIn = "",
83 G4bool checkParametersIn = true,
84 G4bool drawArrowsIn = true,
85 G4bool drawZeroValuePointsIn = true,
86 G4bool drawBoxesIn = true,
87 G4double boxAlphaIn = 0.2);
89
90 G4String name;
91 G4String outfileMagnetic;
92 G4String outfileElectric;
93 G4bool queryMagnetic;
94 G4bool queryElectric;
95 QueryDimensionInfo xInfo;
96 QueryDimensionInfo yInfo;
97 QueryDimensionInfo zInfo;
98 QueryDimensionInfo tInfo;
99
100 std::vector<BDSFourVector<G4double>> pointsToQuery;
101 std::vector<G4String> pointsColumnNames;
102
103 G4AffineTransform globalTransform;
104
105 G4bool overwriteExistingFiles;
106 G4bool printTransform;
107
108 G4String fieldObject;
109
111
112 G4bool drawArrows;
113 G4bool drawZeroValuePoints;
114 G4bool drawBoxes;
115 G4double boxAlpha;
116
118 G4bool SpecificPoints() const {return !pointsToQuery.empty();}
119};
120
121#endif
Holder class for all information required for a field query.
G4bool checkParameters
For internal testing use only.
BDSFieldQueryInfo(const G4String &nameIn, const G4String &outfileMagneticIn, const G4String &outfileElectricIn, G4bool queryMagneticIn, G4bool queryElectricIn, QueryDimensionInfo xInfoIn, QueryDimensionInfo yInfoIn, QueryDimensionInfo zInfoIn, QueryDimensionInfo tInfoIn, const G4AffineTransform &globalTransformIn=G4AffineTransform(), G4bool overwriteExistingFilesIn=false, const G4String &fieldObjectIn="", G4bool printTransformIn=false, G4bool checkParametersIn=true, G4bool drawArrowsIn=true, G4bool drawZeroValuePointsIn=true, G4bool drawBoxesIn=true, G4double boxAlphaIn=0.2)
Usual constructor with number of points to query in each dimension.
G4bool SpecificPoints() const
Whether to query a specific set of points.
G4String fieldObject
Optional for use in interpolator.
QueryDimensionInfo(G4int nIn, G4double minIn, G4double maxIn)