OPAL (Object Oriented Parallel Accelerator Library)  2021.1.99
OPAL
Util.cpp
Go to the documentation of this file.
1 #include "Utilities/Util.h"
2 #include "Physics/Physics.h"
3 #include "OPALrevision.h"
4 
5 #include <boost/regex.hpp>
6 
7 #include <boost/filesystem.hpp>
8 
9 #include <queue>
10 #include <string>
11 #include <iostream>
12 #include <fstream>
13 #include <algorithm>
14 #include <cctype>
15 #include <cmath>
16 
17 namespace Util {
18  std::string getGitRevision() {
19  return std::string(GIT_VERSION);
20  }
21 
22 #define erfinv_a3 -0.140543331
23 #define erfinv_a2 0.914624893
24 #define erfinv_a1 -1.645349621
25 #define erfinv_a0 0.886226899
26 
27 #define erfinv_b4 0.012229801
28 #define erfinv_b3 -0.329097515
29 #define erfinv_b2 1.442710462
30 #define erfinv_b1 -2.118377725
31 #define erfinv_b0 1
32 
33 #define erfinv_c3 1.641345311
34 #define erfinv_c2 3.429567803
35 #define erfinv_c1 -1.62490649
36 #define erfinv_c0 -1.970840454
37 
38 #define erfinv_d2 1.637067800
39 #define erfinv_d1 3.543889200
40 #define erfinv_d0 1
41 
42  double erfinv (double x) // inverse error function
43  {
44  double r;
45  int sign_x;
46 
47  if (x < -1 || x > 1)
48  return NAN;
49 
50  if (x == 0)
51  return 0;
52 
53  if (x > 0)
54  sign_x = 1;
55  else {
56  sign_x = -1;
57  x = -x;
58  }
59 
60  if (x <= 0.7) {
61  double x2 = x * x;
62  r =
63  x * (((erfinv_a3 * x2 + erfinv_a2) * x2 + erfinv_a1) * x2 + erfinv_a0);
64  r /= (((erfinv_b4 * x2 + erfinv_b3) * x2 + erfinv_b2) * x2 +
65  erfinv_b1) * x2 + erfinv_b0;
66  }
67  else {
68  double y = std::sqrt (-std::log ((1 - x) / 2));
69  r = (((erfinv_c3 * y + erfinv_c2) * y + erfinv_c1) * y + erfinv_c0);
70  r /= ((erfinv_d2 * y + erfinv_d1) * y + erfinv_d0);
71  }
72 
73  r = r * sign_x;
74  x = x * sign_x;
75 
76  r -= (std::erf (r) - x) / (2 / std::sqrt (M_PI) * std::exp (-r * r));
77  r -= (std::erf (r) - x) / (2 / std::sqrt (M_PI) * std::exp (-r * r));
78 
79  return r;
80  }
81 
82 #undef erfinv_a3
83 #undef erfinv_a2
84 #undef erfinv_a1
85 #undef erfinv_a0
86 
87 #undef erfinv_b4
88 #undef erfinv_b3
89 #undef erfinv_b2
90 #undef erfinv_b1
91 #undef erfinv_b0
92 
93 #undef erfinv_c3
94 #undef erfinv_c2
95 #undef erfinv_c1
96 #undef erfinv_c0
97 
98 #undef erfinv_d2
99 #undef erfinv_d1
100 #undef erfinv_d0
101 
102  Vector_t getTaitBryantAngles(Quaternion rotation, const std::string &/*elementName*/) {
103  Quaternion rotationBAK = rotation;
104 
105  // y axis
106  Vector_t tmp = rotation.rotate(Vector_t(0, 0, 1));
107  tmp(1) = 0.0;
108  // tmp /= euclidean_norm(tmp);
109  double theta = std::fmod(std::atan2(tmp(0), tmp(2)) + Physics::two_pi, Physics::two_pi);
110 
111  Quaternion rotTheta(std::cos(0.5 * theta), 0, std::sin(0.5 * theta), 0);
112  rotation = rotTheta.conjugate() * rotation;
113 
114  // x axis
115  tmp = rotation.rotate(Vector_t(0, 0, 1));
116  tmp(0) = 0.0;
117  tmp /= euclidean_norm(tmp);
118  double phi = std::fmod(std::atan2(-tmp(1), tmp(2)) + Physics::two_pi, Physics::two_pi);
119 
120  Quaternion rotPhi(std::cos(0.5 * phi), std::sin(0.5 * phi), 0, 0);
121  rotation = rotPhi.conjugate() * rotation;
122 
123  // z axis
124  tmp = rotation.rotate(Vector_t(1, 0, 0));
125  tmp(2) = 0.0;
126  tmp /= euclidean_norm(tmp);
127  double psi = std::fmod(std::atan2(tmp(1), tmp(0)) + Physics::two_pi, Physics::two_pi);
128 
129  return Vector_t(theta, phi, psi);
130  }
131 
132  std::string toUpper(const std::string &str) {
133  std::string output = str;
134  std::transform(output.begin(), output.end(), output.begin(), [](const char c) { return toupper(c);});
135 
136  return output;
137  }
138 
139  std::string combineFilePath(std::initializer_list<std::string> ilist) {
140  boost::filesystem::path path;
141  for (auto entry : ilist) {
142  path /= entry;
143  }
144  return path.string();
145  }
146 
148  sum(0.0),
149  correction(0.0)
150  { }
151 
152 
154  long double y = value - this->correction;
155  long double t = this->sum + y;
156  this->correction = (t - this->sum) - y;
157  this->sum = t;
158  return *this;
159  }
160 
164  unsigned int rewindLinesSDDS(const std::string &fileName, double maxSPos, bool checkForTime) {
165  if (Ippl::myNode() > 0) return 0;
166 
167  std::fstream fs(fileName.c_str(), std::fstream::in);
168  if (!fs.is_open()) return 0;
169 
170  std::string line;
171  std::queue<std::string> allLines;
172  unsigned int numParameters = 0;
173  unsigned int numColumns = 0;
174  unsigned int sposColumnNr = 0;
175  unsigned int timeColumnNr = 0;
176  double spos, time = 0.0;
177  double lastTime = -1.0;
178 
179  boost::regex parameters("&parameter");
180  boost::regex column("&column");
181  boost::regex data("&data");
182  boost::regex end("&end");
183  boost::regex name("name=([a-zA-Z0-9\\$_]+)");
184  boost::smatch match;
185 
186  std::istringstream linestream;
187 
188  while (getline(fs, line)) {
189  allLines.push(line);
190  }
191  fs.close();
192 
193 
194  fs.open (fileName.c_str(), std::fstream::out);
195 
196  if (!fs.is_open()) return 0;
197 
198  do {
199  line = allLines.front();
200  allLines.pop();
201  fs << line << "\n";
202  if (boost::regex_search(line, match, parameters)) {
203  ++numParameters;
204  while (!boost::regex_search(line, match, end)) {
205  line = allLines.front();
206  allLines.pop();
207  fs << line << "\n";
208  }
209  } else if (boost::regex_search(line, match, column)) {
210  ++numColumns;
211  while (!boost::regex_search(line, match, name)) {
212  line = allLines.front();
213  allLines.pop();
214  fs << line << "\n";
215  }
216  if (match[1] == "s") {
217  sposColumnNr = numColumns;
218  }
219  if (match[1] == "t") {
220  timeColumnNr = numColumns;
221  }
222  while (!boost::regex_search(line, match, end)) {
223  line = allLines.front();
224  allLines.pop();
225  fs << line << "\n";
226  }
227  }
228  } while (!boost::regex_search(line, match, data));
229 
230  while (!boost::regex_search(line, match, end)) {
231  line = allLines.front();
232  allLines.pop();
233  fs << line << "\n";
234  }
235 
236  for (unsigned int i = 0; i < numParameters; ++ i) {
237  fs << allLines.front() << "\n";
238  allLines.pop();
239  }
240 
241  while (allLines.size() > 0) {
242  line = allLines.front();
243 
244  linestream.str(line);
245  if (checkForTime) {
246  for (unsigned int i = 0; i < timeColumnNr; ++ i) {
247  linestream >> time;
248  }
249  }
250 
251  linestream.str(line);
252  for (unsigned int i = 0; i < sposColumnNr; ++ i) {
253  linestream >> spos;
254  }
255 
256  if ((spos - maxSPos) > 1e-20 * Physics::c) break;
257 
258  allLines.pop();
259 
260  if (!checkForTime || (time - lastTime) > 1e-20)
261  fs << line << "\n";
262 
263  lastTime = time;
264  }
265 
266  fs.close();
267 
268  if (allLines.size() > 0)
269  INFOMSG(level2 << "rewind " + fileName + " to " + std::to_string(maxSPos) << " m" << endl);
270 
271  return allLines.size();
272  }
273 
274  /*
275  base64.cpp and base64.h
276 
277  Copyright (C) 2004-2008 RenĂ© Nyffenegger
278 
279  This source code is provided 'as-is', without any express or implied
280  warranty. In no event will the author be held liable for any damages
281  arising from the use of this software.
282 
283  Permission is granted to anyone to use this software for any purpose,
284  including commercial applications, and to alter it and redistribute it
285  freely, subject to the following restrictions:
286 
287  1. The origin of this source code must not be misrepresented; you must not
288  claim that you wrote the original source code. If you use this source code
289  in a product, an acknowledgment in the product documentation would be
290  appreciated but is not required.
291 
292  2. Altered source versions must be plainly marked as such, and must not be
293  misrepresented as being the original source code.
294 
295  3. This notice may not be removed or altered from any source distribution.
296 
297  RenĂ© Nyffenegger rene.nyffenegger@adp-gmbh.ch
298 
299  */
300 
301  static const std::string base64_chars = "ABCDEFGHIJKLMNOPQRSTUVWXYZ"
302  "abcdefghijklmnopqrstuvwxyz"
303  "0123456789+/";
304 
305 
306  static inline bool is_base64(unsigned char c) {
307  return (isalnum(c) || (c == '+') || (c == '/'));
308  }
309 
310  std::string base64_encode(const std::string &string_to_encode) {
311  const char* bytes_to_encode = string_to_encode.c_str();
312  unsigned int in_len = string_to_encode.size();
313  std::string ret;
314  int i = 0;
315  int j = 0;
316  unsigned char char_array_3[3];
317  unsigned char char_array_4[4];
318 
319  while (in_len--) {
320  char_array_3[i++] = *(bytes_to_encode++);
321  if (i == 3) {
322  char_array_4[0] = (char_array_3[0] & 0xfc) >> 2;
323  char_array_4[1] = ((char_array_3[0] & 0x03) << 4) + ((char_array_3[1] & 0xf0) >> 4);
324  char_array_4[2] = ((char_array_3[1] & 0x0f) << 2) + ((char_array_3[2] & 0xc0) >> 6);
325  char_array_4[3] = char_array_3[2] & 0x3f;
326 
327  for (i = 0; (i <4) ; i++)
328  ret += base64_chars[char_array_4[i]];
329  i = 0;
330  }
331  }
332 
333  if (i)
334  {
335  for (j = i; j < 3; j++)
336  char_array_3[j] = '\0';
337 
338  char_array_4[0] = (char_array_3[0] & 0xfc) >> 2;
339  char_array_4[1] = ((char_array_3[0] & 0x03) << 4) + ((char_array_3[1] & 0xf0) >> 4);
340  char_array_4[2] = ((char_array_3[1] & 0x0f) << 2) + ((char_array_3[2] & 0xc0) >> 6);
341  char_array_4[3] = char_array_3[2] & 0x3f;
342 
343  for (j = 0; (j < i + 1); j++)
344  ret += base64_chars[char_array_4[j]];
345 
346  while((i++ < 3))
347  ret += '=';
348 
349  }
350 
351  return ret;
352  }
353 
354  std::string base64_decode(std::string const& encoded_string) {
355  int in_len = encoded_string.size();
356  int i = 0;
357  int j = 0;
358  int in_ = 0;
359  unsigned char char_array_4[4], char_array_3[3];
360  std::string ret;
361 
362  while (in_len-- && ( encoded_string[in_] != '=') && is_base64(encoded_string[in_])) {
363  char_array_4[i++] = encoded_string[in_]; in_++;
364  if (i ==4) {
365  for (i = 0; i <4; i++)
366  char_array_4[i] = base64_chars.find(char_array_4[i]);
367 
368  char_array_3[0] = (char_array_4[0] << 2) + ((char_array_4[1] & 0x30) >> 4);
369  char_array_3[1] = ((char_array_4[1] & 0xf) << 4) + ((char_array_4[2] & 0x3c) >> 2);
370  char_array_3[2] = ((char_array_4[2] & 0x3) << 6) + char_array_4[3];
371 
372  for (i = 0; (i < 3); i++)
373  ret += char_array_3[i];
374  i = 0;
375  }
376  }
377 
378  if (i) {
379  for (j = i; j <4; j++)
380  char_array_4[j] = 0;
381 
382  for (j = 0; j <4; j++)
383  char_array_4[j] = base64_chars.find(char_array_4[j]);
384 
385  char_array_3[0] = (char_array_4[0] << 2) + ((char_array_4[1] & 0x30) >> 4);
386  char_array_3[1] = ((char_array_4[1] & 0xf) << 4) + ((char_array_4[2] & 0x3c) >> 2);
387  char_array_3[2] = ((char_array_4[2] & 0x3) << 6) + char_array_4[3];
388 
389  for (j = 0; (j < i - 1); j++) ret += char_array_3[j];
390  }
391 
392  return ret;
393  }
394 }
Tps< T > log(const Tps< T > &x)
Natural logarithm.
Definition: TpsMath.h:182
Tps< T > cos(const Tps< T > &x)
Cosine.
Definition: TpsMath.h:129
Tps< T > exp(const Tps< T > &x)
Exponential.
Definition: TpsMath.h:165
Tps< T > sin(const Tps< T > &x)
Sine.
Definition: TpsMath.h:111
Tps< T > sqrt(const Tps< T > &x)
Square root.
Definition: TpsMath.h:91
T euclidean_norm(const Vector< T > &)
Euclidean norm.
Definition: Vector.h:243
PartBunchBase< T, Dim >::ConstIterator end(PartBunchBase< T, Dim > const &bunch)
FTps< T, N > erf(const FTps< T, N > &x, int trunc=(FTps< T, N >::EXACT))
Error function.
Definition: FTpsMath.h:385
#define erfinv_c1
Definition: Util.cpp:35
#define erfinv_d0
Definition: Util.cpp:40
#define erfinv_a0
Definition: Util.cpp:25
#define erfinv_b2
Definition: Util.cpp:29
#define erfinv_b1
Definition: Util.cpp:30
#define erfinv_a2
Definition: Util.cpp:23
#define erfinv_c2
Definition: Util.cpp:34
#define erfinv_c3
Definition: Util.cpp:33
#define erfinv_d1
Definition: Util.cpp:39
#define erfinv_a1
Definition: Util.cpp:24
#define erfinv_a3
Definition: Util.cpp:22
#define erfinv_b4
Definition: Util.cpp:27
#define erfinv_d2
Definition: Util.cpp:38
#define erfinv_b0
Definition: Util.cpp:31
#define erfinv_b3
Definition: Util.cpp:28
#define erfinv_c0
Definition: Util.cpp:36
#define GIT_VERSION
Definition: OPALrevision.h:1
PETE_TBTree< FnFmod, PETE_Scalar< Vektor< T1, Dim > >, typename T2::PETE_Expr_t > fmod(const Vektor< T1, Dim > &l, const PETE_Expr< T2 > &r)
PETE_TBTree< FnArcTan2, PETE_Scalar< Vektor< T1, Dim > >, typename T2::PETE_Expr_t > atan2(const Vektor< T1, Dim > &l, const PETE_Expr< T2 > &r)
T::PETE_Expr_t::PETE_Return_t sum(const PETE_Expr< T > &expr)
Definition: PETE.h:1111
Inform & endl(Inform &inf)
Definition: Inform.cpp:42
Inform & level2(Inform &inf)
Definition: Inform.cpp:46
#define INFOMSG(msg)
Definition: IpplInfo.h:348
const std::string name
constexpr double two_pi
The value of.
Definition: Physics.h:33
constexpr double e
The value of.
Definition: Physics.h:39
constexpr double c
The velocity of light in m/s.
Definition: Physics.h:51
Definition: Util.cpp:17
std::string combineFilePath(std::initializer_list< std::string > ilist)
Definition: Util.cpp:139
Vector_t getTaitBryantAngles(Quaternion rotation, const std::string &)
Definition: Util.cpp:102
std::string toUpper(const std::string &str)
Definition: Util.cpp:132
double erfinv(double x)
Definition: Util.cpp:42
std::string base64_decode(std::string const &encoded_string)
Definition: Util.cpp:354
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:164
std::string getGitRevision()
Definition: Util.cpp:18
std::string base64_encode(const std::string &string_to_encode)
Definition: Util.cpp:310
FRONT * fs
Definition: hypervolume.cpp:59
Vector_t rotate(const Vector_t &) const
Definition: Quaternion.cpp:122
Quaternion conjugate() const
Definition: Quaternion.h:105
long double sum
Definition: Util.h:188
long double correction
Definition: Util.h:189
KahanAccumulation & operator+=(double value)
Definition: Util.cpp:153
static int myNode()
Definition: IpplInfo.cpp:691
Vektor< double, 3 > Vector_t
Definition: Vektor.h:6