BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSFieldLoaderPoisson.cc
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#include "BDSArray2DCoords.hh"
20#include "BDSDebug.hh"
21#include "BDSException.hh"
22#include "BDSFieldLoaderPoisson.hh"
23#include "BDSFieldValue.hh"
24
25#include "globals.hh" // geant4 types / globals
26
27#include "CLHEP/Units/SystemOfUnits.h"
28
29#include <algorithm>
30#include <fstream>
31#include <iostream>
32#include <regex>
33#include <sstream>
34#include <string>
35#include <vector>
36
37#ifdef USE_GZSTREAM
38#include "src-external/gzstream/gzstream.h"
39#endif
40
41template <class T>
43{;}
44
45template <class T>
47{;}
48
49template <class T>
51{
52 file.open(fileName);
53
54 // test if file is valid
55#ifdef USE_GZSTREAM
56 bool validFile = file.rdbuf()->is_open();
57#else
58 bool validFile = file.is_open();
59#endif
60 if (!validFile)
61 {throw BDSException(__METHOD_NAME__, "Invalid file name or no such file named \"" + fileName + "\"");}
62 else
63 {G4cout << "Loading \"" << fileName << "\"" << G4endl;}
64
65 // Pointer to where result will be stored. Can't be constructed until we know
66 // the dimensions.
67 BDSArray2DCoords* result = nullptr;
68 G4int nX = 1;
69 G4int nY = 1;
70 G4double xMin = 0;
71 G4double yMin = 0;
72 G4double xMax = 0;
73 G4double yMax = 0;
74 G4int nColumns = 0;
75 std::vector<G4String> keys;
76 std::vector<G4String> units;
77
78 // temporary variables
79 G4bool intoData = false;
80 G4bool pastNStep = false;
81 G4bool dataFinished = false;
82 std::string line;
83 G4int indX = 0;
84 G4int indY = 0;
85
86 while (std::getline(file, line))
87 {// read a line only if it's not a blank one
88
89 // Skip a line if it's only whitespace
90 if (std::all_of(line.begin(), line.end(), isspace))
91 {continue;}
92
93 if (intoData)
94 {
95 if (dataFinished)
96 {continue;} // just keep skipping lines after the appropriate amount of data is loaded
97
98 // General data entry
99 std::istringstream liness(line);
100 G4double value = 0;
101 std::vector<G4float> lineData(nColumns); // reserve to avoid copying as it expands
102 for (G4int i = 0; i < nColumns; ++i)
103 {
104 liness >> value;
105 lineData[i] = value;
106 }
107 // Construct field value
108 BDSFieldValue fv = BDSFieldValue(lineData[2], lineData[3], 0);
109 // we could use the gradients here, but we don't
110
111 // Copy into array.
112 (*result)(indX,indY) = fv;
113
114 indX++;
115 if (indX == nX)
116 {// we've completed one set of x
117 indX = 0;
118 indY++;
119 }
120 if (indY == nY)
121 {// we've completed one set of y
122 indY = 0;
123 dataFinished = true;
124 }
125
126 continue;
127 }
128
129 // else try and parser what must be header bits.
130
131 // Min
132 std::regex min("^.*?\\(Xmin,Ymin\\)\\s*=\\s*", std::regex_constants::icase);
133 if (std::regex_search(line, min))
134 {// found line beginning with '(Xmin,YMin)'
135 std::regex minValues("\\(([0-9eE.+-]+),([0-9eE.+-]+)\\)");
136 std::smatch match;
137 std::regex_search(line, match, minValues);
138 xMin = std::stod(match[1])*CLHEP::cm;
139 yMin = std::stod(match[2])*CLHEP::cm;
140 continue;
141 }
142
143 // Max
144 std::regex max("^.*?\\(Xmax,Ymax\\)\\s*=\\s*", std::regex_constants::icase);
145 if (std::regex_search(line, max))
146 {// found line beginning with '(Xmax,YMax)'
147 std::regex maxValues("\\(([0-9eE.+-]+),([0-9eE.+-]+)\\)");
148 std::smatch match;
149 std::regex_search(line, match, maxValues);
150 xMax = std::stod(match[1])*CLHEP::cm;
151 yMax = std::stod(match[2])*CLHEP::cm;
152 continue;
153 }
154
155 // Steps
156 std::regex incr("\\s*X\\s*and\\s*Y", std::regex_constants::icase);
157 if (std::regex_search(line, incr))
158 {//found line for number of increments
159 std::regex nStepValues(":\\s*([0-9eE.+-]+)\\s*([0-9eE.+-]+)");
160 std::smatch match;
161 std::regex_search(line, match, nStepValues);
162 nX = std::stoi(match[1]) + 1;
163 nY = std::stoi(match[2]) + 1; // number of increments, so number of points is +1
164 pastNStep = true;
165 continue;
166 }
167
168 // Columns
169 std::regex noNumbers("^(\\D)*$");
170 G4bool noNumbersR = std::regex_search(line, noNumbers);
171 G4bool hasParenthesis = (line.find("(") != std::string::npos) || (line.find(")") != std::string::npos);
172 G4bool hasNoParenthesis = !hasParenthesis;
173 if (pastNStep && noNumbersR && hasNoParenthesis)
174 {
175 std::istringstream liness(line);
176 for (G4String key; liness >> key;)
177 {keys.push_back(key);}
178 nColumns = G4int(std::distance(keys.begin(), keys.end()));
179 continue;
180 }
181
182 // Units
183 if (pastNStep && noNumbersR && hasParenthesis)
184 {
185 // match text inside (XXX) parenthesis, ie strip them off
186 std::regex unitsRE("[\\(](\\w+\\-*\\/*\\w*)[\\)]");
187 std::smatch match;
188 std::regex_search(line, match, unitsRE);
189 for (const auto& u : match)
190 {units.push_back(G4String(u));}
191 result = new BDSArray2DCoords(nX, nY, xMin, xMax, yMin, yMax);
192 intoData = true;
193 }
194 //nothing heap allocated so safe to let all go out of scope
195 }
196
197 file.close();
198
199 if (!result)
200 {throw BDSException(__METHOD_NAME__, "Invalid file contents - no column header found");}
201
202 return result;
203}
204
205
207
208#ifdef USE_GZSTREAM
210#endif
2D array with spatial mapping derived from BDSArray4DCoords.
General exception with possible name of object and message.
Definition: BDSException.hh:35
Loader for 2D Poisson SuperFish SF7 files.
BDSArray2DCoords * LoadMag2D(G4String fileName)
Load the 2D array of 3-Vector field values.