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