OPAL (Object Oriented Parallel Accelerator Library)  2024.1
OPAL
Util.cpp
Go to the documentation of this file.
1 //
2 // Namespace Util
3 // This namespace contains useful global methods.
4 //
5 // Copyright (c) 200x - 2022, Paul Scherrer Institut, Villigen PSI, Switzerland
6 // All rights reserved
7 //
8 // This file is part of OPAL.
9 //
10 // OPAL is free software: you can redistribute it and/or modify
11 // it under the terms of the GNU General Public License as published by
12 // the Free Software Foundation, either version 3 of the License, or
13 // (at your option) any later version.
14 //
15 // You should have received a copy of the GNU General Public License
16 // along with OPAL. If not, see <https://www.gnu.org/licenses/>.
17 //
18 #include "Utilities/Util.h"
19 
20 #include "OPALrevision.h"
22 
23 #include <boost/regex.hpp>
24 
25 #include <cctype>
26 #include <fstream>
27 #include <filesystem>
28 #include <iostream>
29 #include <iterator>
30 #include <queue>
31 
32 namespace Util {
33  std::string getGitRevision() {
34  return std::string(GIT_VERSION);
35  }
36 
37 #define erfinv_a3 -0.140543331
38 #define erfinv_a2 0.914624893
39 #define erfinv_a1 -1.645349621
40 #define erfinv_a0 0.886226899
41 
42 #define erfinv_b4 0.012229801
43 #define erfinv_b3 -0.329097515
44 #define erfinv_b2 1.442710462
45 #define erfinv_b1 -2.118377725
46 #define erfinv_b0 1
47 
48 #define erfinv_c3 1.641345311
49 #define erfinv_c2 3.429567803
50 #define erfinv_c1 -1.62490649
51 #define erfinv_c0 -1.970840454
52 
53 #define erfinv_d2 1.637067800
54 #define erfinv_d1 3.543889200
55 #define erfinv_d0 1
56 
57  double erfinv (double x) // inverse error function
58  {
59  double r;
60  int sign_x;
61 
62  if (x < -1 || x > 1)
63  return NAN;
64 
65  if (x == 0)
66  return 0;
67 
68  if (x > 0)
69  sign_x = 1;
70  else {
71  sign_x = -1;
72  x = -x;
73  }
74 
75  if (x <= 0.7) {
76  double x2 = x * x;
77  r =
78  x * (((erfinv_a3 * x2 + erfinv_a2) * x2 + erfinv_a1) * x2 + erfinv_a0);
79  r /= (((erfinv_b4 * x2 + erfinv_b3) * x2 + erfinv_b2) * x2 +
80  erfinv_b1) * x2 + erfinv_b0;
81  }
82  else {
83  double y = std::sqrt (-std::log ((1 - x) / 2));
84  r = (((erfinv_c3 * y + erfinv_c2) * y + erfinv_c1) * y + erfinv_c0);
85  r /= ((erfinv_d2 * y + erfinv_d1) * y + erfinv_d0);
86  }
87 
88  r = r * sign_x;
89  x = x * sign_x;
90 
91  r -= (std::erf (r) - x) / (2 / std::sqrt (M_PI) * std::exp (-r * r));
92  r -= (std::erf (r) - x) / (2 / std::sqrt (M_PI) * std::exp (-r * r));
93 
94  return r;
95  }
96 
97 #undef erfinv_a3
98 #undef erfinv_a2
99 #undef erfinv_a1
100 #undef erfinv_a0
101 
102 #undef erfinv_b4
103 #undef erfinv_b3
104 #undef erfinv_b2
105 #undef erfinv_b1
106 #undef erfinv_b0
107 
108 #undef erfinv_c3
109 #undef erfinv_c2
110 #undef erfinv_c1
111 #undef erfinv_c0
112 
113 #undef erfinv_d2
114 #undef erfinv_d1
115 #undef erfinv_d0
116 
117  Vector_t getTaitBryantAngles(Quaternion rotation, const std::string& /*elementName*/) {
118  Quaternion rotationBAK = rotation;
119 
120  // y axis
121  Vector_t tmp = rotation.rotate(Vector_t(0, 0, 1));
122  tmp(1) = 0.0;
123  // tmp /= euclidean_norm(tmp);
124  double theta = std::fmod(std::atan2(tmp(0), tmp(2)) + Physics::two_pi, Physics::two_pi);
125 
126  Quaternion rotTheta(std::cos(0.5 * theta), 0, std::sin(0.5 * theta), 0);
127  rotation = rotTheta.conjugate() * rotation;
128 
129  // x axis
130  tmp = rotation.rotate(Vector_t(0, 0, 1));
131  tmp(0) = 0.0;
132  tmp /= euclidean_norm(tmp);
133  double phi = std::fmod(std::atan2(-tmp(1), tmp(2)) + Physics::two_pi, Physics::two_pi);
134 
135  Quaternion rotPhi(std::cos(0.5 * phi), std::sin(0.5 * phi), 0, 0);
136  rotation = rotPhi.conjugate() * rotation;
137 
138  // z axis
139  tmp = rotation.rotate(Vector_t(1, 0, 0));
140  tmp(2) = 0.0;
141  tmp /= euclidean_norm(tmp);
142  double psi = std::fmod(std::atan2(tmp(1), tmp(0)) + Physics::two_pi, Physics::two_pi);
143 
144  return Vector_t(theta, phi, psi);
145  }
146 
147  std::string toUpper(const std::string& str) {
148  std::string output = str;
149  std::transform(output.begin(), output.end(), output.begin(), [](const char c) { return std::toupper(c);});
150 
151  return output;
152  }
153 
154  std::string boolToUpperString(const bool& b) {
155  std::ostringstream valueStream;
156  valueStream << std::boolalpha << b;
157  std::string output = Util::toUpper(valueStream.str());
158  return output;
159  }
160 
161  std::string boolVectorToUpperString(const std::vector<bool>& b) {
162  std::ostringstream output;
163  if (b.size() > 1) {
164  output << "(";
165  }
166  for (size_t i = 0; i < b.size(); ++i) {
167  output << std::boolalpha << boolToUpperString(b[i]);
168  if (b.size() > 1) {
169  (i < (b.size()-1)) ? (output << ", ") : (output << ")");
170  }
171  }
172 
173  return output.str();
174  }
175 
176  std::string doubleVectorToString(const std::vector<double>& v) {
177  std::vector<std::string> stringVec;
178  stringVec.reserve(v.size());
179  Util::toString(std::begin(v), std::end(v), std::back_inserter(stringVec));
180 
181  std::ostringstream output;
182  if (v.size() > 1) {
183  output << "(";
184  }
185  unsigned int i = 0;
186  for (auto& s: stringVec) {
187  ++i;
188  output << s;
189  if (v.size() > 1) {
190  (i < stringVec.size()) ? (output << ", ") : (output << ")");
191  }
192  }
193 
194  return output.str();
195  }
196 
197  std::string combineFilePath(std::initializer_list<std::string> ilist) {
198  std::filesystem::path path;
199  for (auto entry : ilist) {
200  path /= entry;
201  }
202  return path.string();
203  }
204 
205  void checkInt(double real, std::string name, double tolerance) {
206  real += tolerance; // prevent rounding error
207  if (std::abs(std::floor(real) - real) > 2*tolerance) {
208  throw OpalException("Util::checkInt",
209  "Value for " + name +
210  " should be an integer but a real value was found");
211  }
212  if (std::floor(real) < 0.5) {
213  throw OpalException("Util::checkInt",
214  "Value for " + name + " should be 1 or more");
215  }
216  }
217 
218  bool isAllDigits(const std::string& str) {
219  return std::all_of(str.begin(),
220  str.end(),
221  [](char c) { return std::isdigit(c); });
222  }
223 
225  sum(0.0),
226  correction(0.0)
227  { }
228 
230  long double y = value - this->correction;
231  long double t = this->sum + y;
232  this->correction = (t - this->sum) - y;
233  this->sum = t;
234  return *this;
235  }
236 
240  unsigned int rewindLinesSDDS(const std::string& fileName, double maxSPos, bool checkForTime) {
241  if (Ippl::myNode() > 0) return 0;
242 
243  std::fstream fs(fileName.c_str(), std::fstream::in);
244  if (!fs.is_open()) return 0;
245 
246  std::string line;
247  std::queue<std::string> allLines;
248  unsigned int numParameters = 0;
249  unsigned int numColumns = 0;
250  unsigned int sposColumnNr = 0;
251  unsigned int timeColumnNr = 0;
252  double spos, time = 0.0;
253  double lastTime = -1.0;
254 
255  boost::regex parameters("&parameter");
256  boost::regex column("&column");
257  boost::regex data("&data");
258  boost::regex end("&end");
259  boost::regex name("name=([a-zA-Z0-9\\$_]+)");
260  boost::smatch match;
261 
262  std::istringstream linestream;
263 
264  while (std::getline(fs, line)) {
265  allLines.push(line);
266  }
267  fs.close();
268 
269  fs.open (fileName.c_str(), std::fstream::out);
270 
271  if (!fs.is_open()) return 0;
272 
273  do {
274  line = allLines.front();
275  allLines.pop();
276  fs << line << "\n";
277  if (boost::regex_search(line, match, parameters)) {
278  ++numParameters;
279  while (!boost::regex_search(line, match, end)) {
280  line = allLines.front();
281  allLines.pop();
282  fs << line << "\n";
283  }
284  } else if (boost::regex_search(line, match, column)) {
285  ++numColumns;
286  while (!boost::regex_search(line, match, name)) {
287  line = allLines.front();
288  allLines.pop();
289  fs << line << "\n";
290  }
291  if (match[1] == "s") {
292  sposColumnNr = numColumns;
293  }
294  if (match[1] == "t") {
295  timeColumnNr = numColumns;
296  }
297  while (!boost::regex_search(line, match, end)) {
298  line = allLines.front();
299  allLines.pop();
300  fs << line << "\n";
301  }
302  }
303  } while (!boost::regex_search(line, match, data));
304 
305  while (!boost::regex_search(line, match, end)) {
306  line = allLines.front();
307  allLines.pop();
308  fs << line << "\n";
309  }
310 
311  for (unsigned int i = 0; i < numParameters; ++ i) {
312  fs << allLines.front() << "\n";
313  allLines.pop();
314  }
315 
316  while (!allLines.empty()) {
317  line = allLines.front();
318 
319  linestream.str(line);
320  if (checkForTime) {
321  for (unsigned int i = 0; i < timeColumnNr; ++ i) {
322  linestream >> time;
323  }
324  }
325 
326  linestream.str(line);
327  for (unsigned int i = 0; i < sposColumnNr; ++ i) {
328  linestream >> spos;
329  }
330 
331  if ((spos - maxSPos) > 1e-20 * Physics::c) break;
332 
333  allLines.pop();
334 
335  if (!checkForTime || (time - lastTime) > 1e-20)
336  fs << line << "\n";
337 
338  lastTime = time;
339  }
340 
341  fs.close();
342 
343  if (!allLines.empty())
344  INFOMSG(level2 << "rewind " + fileName + " to " + std::to_string(maxSPos) << " m" << endl);
345 
346  return allLines.size();
347  }
348 
349  /*
350  base64.cpp and base64.h
351 
352  Copyright (C) 2004-2008 RenĂ© Nyffenegger
353 
354  This source code is provided 'as-is', without any express or implied
355  warranty. In no event will the author be held liable for any damages
356  arising from the use of this software.
357 
358  Permission is granted to anyone to use this software for any purpose,
359  including commercial applications, and to alter it and redistribute it
360  freely, subject to the following restrictions:
361 
362  1. The origin of this source code must not be misrepresented; you must not
363  claim that you wrote the original source code. If you use this source code
364  in a product, an acknowledgment in the product documentation would be
365  appreciated but is not required.
366 
367  2. Altered source versions must be plainly marked as such, and must not be
368  misrepresented as being the original source code.
369 
370  3. This notice may not be removed or altered from any source distribution.
371 
372  RenĂ© Nyffenegger rene.nyffenegger@adp-gmbh.ch
373 
374  */
375 
376  static const std::string base64_chars = "ABCDEFGHIJKLMNOPQRSTUVWXYZ"
377  "abcdefghijklmnopqrstuvwxyz"
378  "0123456789+/";
379 
380  static inline bool is_base64(unsigned char c) {
381  return (std::isalnum(c) || (c == '+') || (c == '/'));
382  }
383 
384  std::string base64_encode(const std::string& string_to_encode) {
385  const char* bytes_to_encode = string_to_encode.c_str();
386  unsigned int in_len = string_to_encode.size();
387  std::string ret;
388  int i = 0;
389  int j = 0;
390  unsigned char char_array_3[3];
391  unsigned char char_array_4[4];
392 
393  while (in_len--) {
394  char_array_3[i++] = *(bytes_to_encode++);
395  if (i == 3) {
396  char_array_4[0] = (char_array_3[0] & 0xfc) >> 2;
397  char_array_4[1] = ((char_array_3[0] & 0x03) << 4) + ((char_array_3[1] & 0xf0) >> 4);
398  char_array_4[2] = ((char_array_3[1] & 0x0f) << 2) + ((char_array_3[2] & 0xc0) >> 6);
399  char_array_4[3] = char_array_3[2] & 0x3f;
400 
401  for (i = 0; (i <4) ; i++)
402  ret += base64_chars[char_array_4[i]];
403  i = 0;
404  }
405  }
406 
407  if (i)
408  {
409  for (j = i; j < 3; j++)
410  char_array_3[j] = '\0';
411 
412  char_array_4[0] = (char_array_3[0] & 0xfc) >> 2;
413  char_array_4[1] = ((char_array_3[0] & 0x03) << 4) + ((char_array_3[1] & 0xf0) >> 4);
414  char_array_4[2] = ((char_array_3[1] & 0x0f) << 2) + ((char_array_3[2] & 0xc0) >> 6);
415  char_array_4[3] = char_array_3[2] & 0x3f;
416 
417  for (j = 0; (j < i + 1); j++)
418  ret += base64_chars[char_array_4[j]];
419 
420  while((i++ < 3))
421  ret += '=';
422 
423  }
424 
425  return ret;
426  }
427 
428  std::string base64_decode(std::string const& encoded_string) {
429  int in_len = encoded_string.size();
430  int i = 0;
431  int j = 0;
432  int in_ = 0;
433  unsigned char char_array_4[4], char_array_3[3];
434  std::string ret;
435 
436  while (in_len-- && ( encoded_string[in_] != '=') && is_base64(encoded_string[in_])) {
437  char_array_4[i++] = encoded_string[in_]; in_++;
438  if (i ==4) {
439  for (i = 0; i <4; i++)
440  char_array_4[i] = base64_chars.find(char_array_4[i]);
441 
442  char_array_3[0] = (char_array_4[0] << 2) + ((char_array_4[1] & 0x30) >> 4);
443  char_array_3[1] = ((char_array_4[1] & 0xf) << 4) + ((char_array_4[2] & 0x3c) >> 2);
444  char_array_3[2] = ((char_array_4[2] & 0x3) << 6) + char_array_4[3];
445 
446  for (i = 0; (i < 3); i++)
447  ret += char_array_3[i];
448  i = 0;
449  }
450  }
451 
452  if (i) {
453  for (j = i; j <4; j++)
454  char_array_4[j] = 0;
455 
456  for (j = 0; j <4; j++)
457  char_array_4[j] = base64_chars.find(char_array_4[j]);
458 
459  char_array_3[0] = (char_array_4[0] << 2) + ((char_array_4[1] & 0x30) >> 4);
460  char_array_3[1] = ((char_array_4[1] & 0xf) << 4) + ((char_array_4[2] & 0x3c) >> 2);
461  char_array_3[2] = ((char_array_4[2] & 0x3) << 6) + char_array_4[3];
462 
463  for (j = 0; (j < i - 1); j++) ret += char_array_3[j];
464  }
465 
466  return ret;
467  }
468 }
#define erfinv_c3
Definition: Util.cpp:48
Tps< T > sqrt(const Tps< T > &x)
Square root.
Definition: TpsMath.h:91
constexpr double c
The velocity of light in m/s.
Definition: Physics.h:45
std::string boolVectorToUpperString(const std::vector< bool > &b)
Definition: Util.cpp:161
Quaternion conjugate() const
Definition: Quaternion.h:103
Vector_t rotate(const Vector_t &) const
Definition: Quaternion.cpp:122
Vector_t getTaitBryantAngles(Quaternion rotation, const std::string &)
Definition: Util.cpp:117
#define erfinv_a0
Definition: Util.cpp:40
constexpr double two_pi
The value of .
Definition: Physics.h:33
PETE_TUTree< FnAbs, typename T::PETE_Expr_t > abs(const PETE_Expr< T > &l)
FTps< T, N > erf(const FTps< T, N > &x, int trunc=(FTps< T, N >::EXACT))
Error function.
Definition: FTpsMath.h:385
#define erfinv_a1
Definition: Util.cpp:39
long double correction
Definition: Util.h:235
static int myNode()
Definition: IpplInfo.cpp:691
std::string boolToUpperString(const bool &b)
Definition: Util.cpp:154
std::string getGitRevision()
Definition: Util.cpp:33
Vektor< double, 3 > Vector_t
Definition: Vektor.h:6
Definition: TSVMeta.h:24
#define erfinv_d1
Definition: Util.cpp:54
#define erfinv_d0
Definition: Util.cpp:55
unsigned int rewindLinesSDDS(const std::string &fileName, double maxSPos, bool checkForTime)
rewind the SDDS file such that the spos of the last step is less or equal to maxSPos ...
Definition: Util.cpp:240
Tps< T > exp(const Tps< T > &x)
Exponential.
Definition: TpsMath.h:165
std::string base64_encode(const std::string &string_to_encode)
Definition: Util.cpp:384
clearpage the user may choose between constant or variable radius This model includes fringe fields begin
Definition: multipole_t.tex:6
std::string toUpper(const std::string &str)
Definition: Util.cpp:147
#define INFOMSG(msg)
Definition: IpplInfo.h:348
Inform & endl(Inform &inf)
Definition: Inform.cpp:42
#define erfinv_c2
Definition: Util.cpp:49
PETE_TBTree< FnArcTan2, PETE_Scalar< Vektor< T1, Dim > >, typename T2::PETE_Expr_t > atan2(const Vektor< T1, Dim > &l, const PETE_Expr< T2 > &r)
#define erfinv_a3
Definition: Util.cpp:37
T euclidean_norm(const Vector< T > &)
Euclidean norm.
Definition: Vector.h:243
long double sum
Definition: Util.h:234
T::PETE_Expr_t::PETE_Return_t sum(const PETE_Expr< T > &expr)
Definition: PETE.h:1111
FRONT * fs
Definition: hypervolume.cpp:59
#define erfinv_b2
Definition: Util.cpp:44
The base class for all OPAL exceptions.
Definition: OpalException.h:28
#define erfinv_a2
Definition: Util.cpp:38
bool isAllDigits(const std::string &str)
Definition: Util.cpp:218
#define erfinv_b0
Definition: Util.cpp:46
#define erfinv_b3
Definition: Util.cpp:43
#define erfinv_b1
Definition: Util.cpp:45
having only three parameters(the centre length $s_0 $and the fringe field lengths $\lambda_{left}$,$\lambda_{right}$)
Definition: multipole_t.tex:58
Tps< T > cos(const Tps< T > &x)
Cosine.
Definition: TpsMath.h:129
double erfinv(double x)
Definition: Util.cpp:57
PETE_TBTree< FnFmod, PETE_Scalar< Vektor< T1, Dim > >, typename T2::PETE_Expr_t > fmod(const Vektor< T1, Dim > &l, const PETE_Expr< T2 > &r)
Tps< T > log(const Tps< T > &x)
Natural logarithm.
Definition: TpsMath.h:182
const std::string name
std::string combineFilePath(std::initializer_list< std::string > ilist)
Definition: Util.cpp:197
#define erfinv_b4
Definition: Util.cpp:42
std::string doubleVectorToString(const std::vector< double > &v)
Definition: Util.cpp:176
constexpr double e
The value of .
Definition: Physics.h:39
Inform & level2(Inform &inf)
Definition: Inform.cpp:46
#define erfinv_d2
Definition: Util.cpp:53
#define erfinv_c0
Definition: Util.cpp:51
#define erfinv_c1
Definition: Util.cpp:50
Tps< T > sin(const Tps< T > &x)
Sine.
Definition: TpsMath.h:111
void checkInt(double real, std::string name, double tolerance)
Definition: Util.cpp:205
PETE_TUTree< FnFloor, typename T::PETE_Expr_t > floor(const PETE_Expr< T > &l)
Definition: PETE.h:733
std::string base64_decode(std::string const &encoded_string)
Definition: Util.cpp:428
FLieGenerator< T, N > real(const FLieGenerator< std::complex< T >, N > &)
Take real part of a complex generator.
end
Definition: multipole_t.tex:9
void toString(IteratorIn first, IteratorIn last, IteratorOut out)
Definition: Util.h:288
KahanAccumulation & operator+=(double value)
Definition: Util.cpp:229