BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSFieldLoaderBDSIM.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 "BDSArray1DCoords.hh"
20#include "BDSArray2DCoords.hh"
21#include "BDSArray3DCoords.hh"
22#include "BDSArray4DCoords.hh"
23#include "BDSDebug.hh"
24#include "BDSException.hh"
25#include "BDSFieldLoaderBDSIM.hh"
26#include "BDSFieldValue.hh"
27#include "BDSUtilities.hh"
28
29#include "globals.hh"
30#include "G4String.hh"
31
32#include "CLHEP/Units/SystemOfUnits.h"
33
34#include <algorithm>
35#include <array>
36#include <exception>
37#include <fstream>
38#include <iostream>
39#include <map>
40#include <regex>
41#include <sstream>
42#include <string>
43#include <vector>
44
45#ifdef USE_GZSTREAM
46#include "src-external/gzstream/gzstream.h"
47#endif
48
49template <class T>
51 nColumns(0),
52 fv(BDSFieldValue()),
53 result(nullptr),
54 headerMustBePositiveKeys({"nx", "ny", "nz", "nt"}),
55 indexOfFirstFieldValue(0)
56{
57 dimKeyMap = {
58 {BDSDimensionType::x, {"nx", "xmin", "xmax"}},
59 {BDSDimensionType::y, {"ny", "ymin", "ymax"}},
60 {BDSDimensionType::z, {"nz", "zmin", "zmax"}},
61 {BDSDimensionType::t, {"nt", "tmin", "tmax"}}
62 };
63}
64
65template <class T>
67{;}
68
69template <class T>
71{
72 nColumns = 0;
73 lineData.clear();
74 fv = BDSFieldValue();
75 std::vector<G4String> allKeys = {"nx", "ny", "nz", "nt",
76 "xmin", "xmax", "ymin", "ymax",
77 "zmin", "zmax", "tmin", "tmax"};
78 for (auto& key : allKeys)
79 {header[key] = 0;}
80 for (auto& key : headerMustBePositiveKeys)
81 {header[key] = 1;}
82 result = nullptr;
83 loopOrder = "xyzt";
84 indexOfFirstFieldValue = 0;
85}
86
87template <class T>
88void BDSFieldLoaderBDSIM<T>::Terminate(const G4String& message)
89{
90 file.close();
91 throw BDSException("BDSFieldLoaderBDSIM", message);
92}
93
94template <class T>
96{
97 Load(fileName,1);
98 return static_cast<BDSArray1DCoords*>(result);
99}
100
101template <class T>
103{
104 Load(fileName,2);
105 return static_cast<BDSArray2DCoords*>(result);
106}
107
108template <class T>
110{
111 Load(fileName,3);
112 return static_cast<BDSArray3DCoords*>(result);
113}
114
115template <class T>
117{
118 Load(fileName,4);
119 return result;
120}
121
122template <class T>
123void BDSFieldLoaderBDSIM<T>::Load(const G4String& fileName,
124 const unsigned int nDim)
125{
126 G4String functionName = "BDSIM Field Format> ";
127 CleanUp();
128
129 file.open(fileName);
130 long long unsigned int currentLineNumber = 0;
131 // test if file is valid
132#ifdef USE_GZSTREAM
133 bool validFile = file.rdbuf()->is_open();
134#else
135 bool validFile = file.is_open();
136#endif
137 if (!validFile)
138 {throw BDSException(__METHOD_NAME__, "Invalid file name or no such file named \"" + fileName + "\"");}
139 else
140 {G4cout << functionName << "Loading \"" << fileName << "\"" << G4endl;}
141
142 // temporary variables
143 unsigned long xIndex = 0;
144 unsigned long yIndex = 0;
145 unsigned long zIndex = 0;
146 G4bool intoData = false;
147 std::string line;
148
149 // Number of points in each dimension - not necessarily x,y,z,t.
150 // For example, could be xz
151 // These are always in xyzt order
152 G4int n1 = 1, n2 = 1, n3 = 1, n4 = 1;
153
154 // wrap in vectors for easy assignment
155 G4String nominalOrder = "xyzt";
156
157 float maximumFieldValue = 0;
158 float minimumFieldValue = 0;
159 // read top of the file
160 while (!intoData)
161 {
162 std::getline(file, line);
163
164 currentLineNumber++;
165 // Skip a line if it's only whitespace
166 if (std::all_of(line.begin(), line.end(), isspace))
167 {continue;}
168
169 if (currentLineNumber > 100) // !intoData is always true as we're inside this while loop
170 {
171 G4String msg = functionName+"did not encounter column specification line starting with '!' in the first 100 lines\n";
172 msg += "This does not appear to be a BDSIM-format field map.";
173 Terminate(msg);
174 }
175
176 // now several checks
177 // 1) key > number
178 // 2) key > word
179 // 3) row starting '!' (column row ->then constructs data)
180
181 // key definition
182 // if (line.find(">") != std::string::npos)
183 // NOTE the use of a regex expression here allows us to safely parse the first line
184 // in a tar.gz file which normally has a load of bumf upfront for the tar file. Each
185 // subsequent line comes out normally with the gzstream reader.
186 std::smatch matchHeaderNumber;
187 std::regex keyValue(R"((\w*)\s*>\s*([0-9eE.+-]+))");
188 if (std::regex_search(line, matchHeaderNumber, keyValue))
189 {// must be key definition
190 if (matchHeaderNumber.size() < 2)
191 {Terminate(functionName+"Invalid key definition in field format: " + line);}
192 else
193 {
194 G4String key = G4String(matchHeaderNumber[1]);
195 key = BDS::LowerCase(key);
196
197 // check it's a valid key - header preloaded with valid keys
198 if (header.find(key) == header.end())
199 {Terminate(functionName+"Invalid key \"" + key +"\" in header");}
200
201 G4double value = 0;
202 try
203 {value = std::stod(matchHeaderNumber[2]);}
204 catch (const std::invalid_argument&)
205 {Terminate(functionName+"Invalid argument " + std::string(matchHeaderNumber[2]));}
206 catch (const std::out_of_range&)
207 {Terminate(functionName+"Number out of range " + std::string(matchHeaderNumber[2]));}
208
209 if (std::find(headerMustBePositiveKeys.begin(), headerMustBePositiveKeys.end(), key) != headerMustBePositiveKeys.end())
210 {
211 if (value < 1)
212 {Terminate(functionName+"Number of points in dimension must be greater than 0 -> see \"" + key + "\"");}
213 }
214
215 header[key] = value;
216 continue;
217 }
218 }
219
220 std::smatch matchHeaderString;
221 // matches "key > string" where string is 1-4 characters (not numbers)
222 // can be padded between each part with whitespace \s*
223 // not more than four characters (via \b for word boundary)
224 std::regex keyWord(R"((\w+)\s*>\s*([a-zA-Z]{1,4})\b\s*)");
225 if (std::regex_search(line, matchHeaderString, keyWord))
226 {
227 loopOrder = G4String(matchHeaderString[2]); // member variable
228
229 G4bool test1 = loopOrder == "xyzt";
230 G4bool test2 = loopOrder == "tzyx";
231 if (! (test1 || test2) )
232 {Terminate(functionName+"loopOrder> header variable is not one of (\"xyzt\" or \"tzyx\")");}
233 continue;
234 }
235
236 // if starts with '!' - columns
237 // check all required keys have been built up ok
238 // set nColumns
239 std::regex columnRow("^\\s*!"); // ignore any initial white space and look for '!'
240 if (std::regex_search(line, columnRow))
241 {
242 // we only need to record the number of columns and which ones are
243 // the x,y,z field component ones.
244 std::regex afterExclamation("\\s*!\\s*(.+)");
245 std::smatch match;
246 std::regex_search(line, match, afterExclamation);
247 std::string restOfLine = match[1];
248 std::string columnName;
249 std::istringstream restOfLineSS(restOfLine);
250 std::vector<G4String> columnNames;
251 while (restOfLineSS >> columnName)
252 {
253 nColumns++;
254 if (columnName.find("Fx") != std::string::npos)
255 {xIndex = nColumns; continue;}
256 else if (columnName.find("Fy") != std::string::npos)
257 {yIndex = nColumns; continue;}
258 else if (columnName.find("Fz") != std::string::npos)
259 {zIndex = nColumns; continue;}
260 else
261 {columnNames.emplace_back(columnName);}
262 }
263 lineData.resize(nColumns + 1); // +1 for default value
264 indexOfFirstFieldValue = std::min({xIndex, yIndex, zIndex});
265
266 for (const auto& key : headerMustBePositiveKeys)
267 {
268 if (header[key] < 1)
269 {Terminate(functionName+"Number of points in dimension must be greater than 0 -> see \"" + key + "\"");}
270 }
271
272 if (nColumns < (nDim + 3)) // 3 for field components
273 {Terminate(functionName+"Too few columns for " + std::to_string(nDim) + "D field loading");}
274
275 // we have all the information now, so initialise the container
276 switch (nDim)
277 {
278 case 1:
279 {
280 BDSDimensionType firstDim = BDS::DetermineDimensionType(columnNames[0]);
281 auto keys = dimKeyMap[firstDim];
282 n1 = G4int(header[keys.number]);
283 result = new BDSArray1DCoords(n1,
284 header[keys.min] * CLHEP::cm,
285 header[keys.max] * CLHEP::cm,
286 firstDim);
287 break;
288 }
289 case 2:
290 {
291 BDSDimensionType firstDim = BDS::DetermineDimensionType(columnNames[0]);
292 BDSDimensionType secondDim = BDS::DetermineDimensionType(columnNames[1]);
293 auto fKeys = dimKeyMap[firstDim];
294 auto sKeys = dimKeyMap[secondDim];
295 n1 = G4int(header[fKeys.number]);
296 n2 = G4int(header[sKeys.number]);
297 result = new BDSArray2DCoords(n1, n2,
298 header[fKeys.min] * CLHEP::cm,
299 header[fKeys.max] * CLHEP::cm,
300 header[sKeys.min] * CLHEP::cm,
301 header[sKeys.max] * CLHEP::cm,
302 firstDim,
303 secondDim);
304 break;
305 }
306 case 3:
307 {
308 BDSDimensionType firstDim = BDS::DetermineDimensionType(columnNames[0]);
309 BDSDimensionType secondDim = BDS::DetermineDimensionType(columnNames[1]);
310 BDSDimensionType thirdDim = BDS::DetermineDimensionType(columnNames[2]);
311 auto fKeys = dimKeyMap[firstDim];
312 auto sKeys = dimKeyMap[secondDim];
313 auto tKeys = dimKeyMap[thirdDim];
314 n1 = G4int(header[fKeys.number]);
315 n2 = G4int(header[sKeys.number]);
316 n3 = G4int(header[tKeys.number]);
317 result = new BDSArray3DCoords(n1, n2, n3,
318 header[fKeys.min] * CLHEP::cm,
319 header[fKeys.max] * CLHEP::cm,
320 header[sKeys.min] * CLHEP::cm,
321 header[sKeys.max] * CLHEP::cm,
322 header[tKeys.min] * CLHEP::cm,
323 header[tKeys.max] * CLHEP::cm,
324 firstDim,
325 secondDim,
326 thirdDim);
327 break;
328 }
329 case 4:
330 {
331 n1 = G4int(header["nx"]);
332 n2 = G4int(header["ny"]);
333 n3 = G4int(header["nz"]);
334 n4 = G4int(header["nt"]);
335 result = new BDSArray4DCoords(n1, n2, n3, n4,
336 header["xmin"] * CLHEP::cm, header["xmax"] * CLHEP::cm,
337 header["ymin"] * CLHEP::cm, header["ymax"] * CLHEP::cm,
338 header["zmin"] * CLHEP::cm, header["zmax"] * CLHEP::cm,
339 header["tmin"] * CLHEP::s, header["tmax"] * CLHEP::s);
340 break;
341 }
342 default:
343 {break;}
344 }
345 intoData = true;
346 continue;
347 }
348 }
349
350 // now only read data - two loops, one for each way of looping
351 if (loopOrder == "tzyx")
352 {
353 for (G4int i = 0; i < n1; i++)
354 {
355 for (G4int j = 0; j < n2; j++)
356 {
357 for (G4int k = 0; k < n3; k++)
358 {
359 for (G4int l = 0; l < n4; l++)
360 {
361 if (!std::getline(file, line))
362 {Terminate(functionName + "unexpected end to file on line number " + std::to_string(currentLineNumber));}
363 currentLineNumber++;
364 if (std::all_of(line.begin(), line.end(), isspace))
365 {continue;}
366 ProcessData(line, xIndex, yIndex, zIndex); // changes member fv
367 (*result)(i, j, k, l) = fv;
368 float mag = fv.mag();
369 maximumFieldValue = std::max(maximumFieldValue, mag);
370 minimumFieldValue = std::min(minimumFieldValue, mag);
371 }
372 }
373 }
374 }
375 }
376 else
377 {// xyzt
378 for (G4int l = 0; l < n4; l++)
379 {
380 for (G4int k = 0; k < n3; k++)
381 {
382 for (G4int j = 0; j < n2; j++)
383 {
384 for (G4int i = 0; i < n1; i++)
385 {
386 if (!std::getline(file, line))
387 {Terminate(functionName + "unexpected end to file on line number " + std::to_string(currentLineNumber));}
388 currentLineNumber++;
389 if (std::all_of(line.begin(), line.end(), isspace))
390 {continue;}
391 ProcessData(line, xIndex, yIndex, zIndex); // changes member fv
392 (*result)(i, j, k, l) = fv;
393 float mag = fv.mag();
394 maximumFieldValue = std::max(maximumFieldValue, mag);
395 minimumFieldValue = std::min(minimumFieldValue, mag);
396 }
397 }
398 }
399 }
400 }
401
402 file.close();
403 G4cout << functionName << "Loaded " << currentLineNumber << " lines from file" << G4endl;
404 G4cout << functionName << "(Min | Max) field magnitudes in loaded file (before scaling): (" << minimumFieldValue << " | " << maximumFieldValue << ")" << G4endl;
405}
406
407template <class T>
408void BDSFieldLoaderBDSIM<T>::ProcessData(const std::string& line,
409 const unsigned long xIndex,
410 const unsigned long yIndex,
411 const unsigned long zIndex)
412{
413 std::istringstream liness(line);
414 G4float value = 0;
415 std::fill(lineData.begin(), lineData.end(), 0); // reset data - technically unnecessary
416
417 // read all columns - indices shifted +1 for default value offset
418 for (unsigned long i = 1; i < nColumns+1; ++i)
419 {
420 liness >> value;
421 if (i < indexOfFirstFieldValue)// coordinates before this index
422 {value *= CLHEP::cm;}
423 lineData[i] = value;
424 }
425
426 // Construct field value
427 fv = BDSFieldValue(lineData[xIndex],
428 lineData[yIndex],
429 lineData[zIndex]);
430}
431
432
434
435#ifdef USE_GZSTREAM
436template class BDSFieldLoaderBDSIM<igzstream>;
437#endif
1D array with spatial mapping derived from BDSArray4DCoords.
2D array with spatial mapping derived from BDSArray4DCoords.
3D array with spatial mapping derived from BDSArray4DCoords.
Overlay of 4D array that provides uniform only spatial coordinate mapping.
General exception with possible name of object and message.
Definition: BDSException.hh:35
Loader for BDSIM format fields.
BDSArray3DCoords * Load3D(const G4String &fileName)
Load a 3D array.
void CleanUp()
Ensure any member variables are reset between usages.
BDSArray4DCoords * Load4D(const G4String &fileName)
Load a 4D array.
BDSArray1DCoords * Load1D(const G4String &fileName)
Load a 1D array.
void Terminate(const G4String &message="")
Close file and exit program in case of an error.
void Load(const G4String &fileName, const unsigned int nDim)
General loader for any number of dimensions.
void ProcessData(const std::string &line, const unsigned long xIndex, const unsigned long yIndex=0, const unsigned long zIndex=0)
BDSArray2DCoords * Load2D(const G4String &fileName)
Load a 2D array.
BDSDimensionType DetermineDimensionType(G4String dimensionType)
Determine the output format to be used from the input string.
G4String LowerCase(const G4String &str)
Utility function to simplify lots of syntax changes for pedantic g4 changes.