MARLEY (Model of Argon Reaction Low Energy Yields)  v1.2.0
A Monte Carlo event generator for tens-of-MeV neutrino interactions
marley_utils.hh
1 //
5 // This file is part of MARLEY (Model of Argon Reaction Low Energy Yields)
6 //
7 // MARLEY is free software: you can redistribute it and/or modify it under the
8 // terms of version 3 of the GNU General Public License as published by the
9 // Free Software Foundation.
10 //
11 // For the full text of the license please see COPYING or
12 // visit http://opensource.org/licenses/GPL-3.0
13 //
14 // Please respect the MCnet academic usage guidelines. See GUIDELINES
15 // or visit https://www.montecarlonet.org/GUIDELINES for details.
16 
17 // Nonstandard but widely-supported (see
18 // http://en.wikipedia.org/wiki/Pragma_once) preprocessor directive that
19 // prevents this file from being included multiple times. Another option is an
20 // include guard (http://en.wikipedia.org/wiki/Include_guard).
21 #pragma once
22 
23 #include <algorithm>
24 #include <chrono>
25 #include <cmath>
26 #include <complex>
27 #include <functional>
28 #include <limits>
29 #include <random>
30 #include <regex>
31 #include <sstream>
32 #include <string>
33 #include <unordered_map>
34 
35 namespace marley_utils {
36 
37  // Frequently used particle IDs
38  constexpr int PHOTON = 22;
39  constexpr int ELECTRON = 11;
40  constexpr int POSITRON = -11;
41  constexpr int ELECTRON_NEUTRINO = 12;
42  constexpr int ELECTRON_ANTINEUTRINO = -12;
43  constexpr int MUON = 13;
44  constexpr int MUON_NEUTRINO = 14;
45  constexpr int MUON_ANTINEUTRINO = -14;
46  constexpr int TAU = 15;
47  constexpr int TAU_NEUTRINO = 16;
48  constexpr int TAU_ANTINEUTRINO = -16;
49  constexpr int NEUTRON = 2112;
50  constexpr int PROTON = 2212;
51  constexpr int DEUTERON = 1000010020;
52  constexpr int TRITON = 1000010030;
53  constexpr int HELION = 1000020030;
54  constexpr int ALPHA = 1000020040;
55 
56  // Dummy double value representing an unknown maximum PDF value. Signals to
57  // marley::Generator::rejection_sample that it needs to search for the
58  // maximum before doing the sampling.
59  constexpr double UNKNOWN_MAX = std::numeric_limits<double>::infinity();
60 
61  // Fermi coupling constant (MeV^(-2))
62  constexpr double GF = 1.16637e-11;
63  // Square of the Fermi coupling constant (MeV^(-4))
64  constexpr double GF2 = GF * GF;
65 
66  // Absolute value of the CKM matrix element for mixing between the up and
67  // down quarks abs(V_ud)
68  constexpr double Vud = 0.97427;
69  // Square of abs(V_ud)
70  constexpr double Vud2 = Vud * Vud;
71 
72  // sin^2(theta_W) (sine squared of the weak mixing angle)
73  // Effective value taken from 2014 PDG Review of Particle Physics,
74  // Table 1.1, "Physical Constants"
75  constexpr double sin2thetaw = 0.23155;
76 
77  // Conversion factor to use when expressing ENSDF energies (keV) in
78  // standard MARLEY energy units (MeV)
79  constexpr double MeV = 1e-3;
80 
81  // Conversion factor to use when expressing atomic masses (micro-amu)
82  // in standard MARLEY energy units (MeV)
83  constexpr double micro_amu = 0.000931494061;
84 
85  // Infinities
86  constexpr double infinity = std::numeric_limits<double>::max();
87  constexpr double minus_infinity = -infinity;
88 
89  // Muon mass
90  constexpr double m_mu = 113428.9267; // micro-amu
91 
92  // Consistent value of pi to use throughout all of MARLEY
93  constexpr double pi = M_PI;
94  constexpr double two_pi = 2.*pi;
95  const double sqrt_two_pi = std::sqrt(pi);
96  constexpr double half_pi = pi/2.0;
97 
98  // Imaginary unit
99  constexpr std::complex<double> i(0, 1);
100 
101  // Natural logarithm of 2
102  const double log_2 = std::log(2);
103 
104  // The physical constants given here were all taken from
105  // the 2014 edition of the Review of Particle Physics
106  // published by the Particle Data Group.
107 
108  // Fine structure constant
109  constexpr double alpha = 7.2973525698e-3;
110  // Conversion factor used to switch to natural units (hbar = c = 1)
111  constexpr double hbar_c = 197.3269718; // MeV*fm
112  constexpr double hbar_c2 = hbar_c * hbar_c; // MeV^2 * fm^2
113  // Electron mass
114  constexpr double m_e = 0.510998928; // MeV
115  // Constant to use when converting from mb to MeV^(-2)
116  constexpr double mb = 1/3.89379338e5; // MeV^(-2) mb^(-1)
117  // Constant to use to convert from fm^2 to 10^(-40) cm^2
118  constexpr double fm2_to_minus40_cm2 = 1e14;
119  // Square of the elementary charge
120  constexpr double e2 = hbar_c * alpha; // MeV*fm
121  // Constant to use when approximating nuclear radii via
122  // r = r0 * A^(1/3), where A is the nucleus's mass number.
123  // See, for example, Introductory Nuclear Physics by Kenneth S. Krane.
124  constexpr double r0 = 1.2; // fm
125  // Handy constants for the fractions 1/2 and 1/3
126  constexpr double ONE_HALF = 1.0/2.0;
127  constexpr double ONE_THIRD = 1.0/3.0;
128 
129  // Strings to use for latex table output of DecayScheme objects
130  extern std::string latex_table_1, latex_table_2, latex_table_3, latex_table_4;
131 
132  // Create an ENSDF nucid string given a nuclide's atomic number Z
133  // and mass number A
134  std::string nuc_id(int Z, int A);
135 
136  // Return the PDG particle ID that corresponds to a ground-state
137  // nucleus with atomic number Z and mass number A
138  inline int get_nucleus_pid(int Z, int A) {
139  if (Z == 0 && A == 1) return NEUTRON;
140  else if (Z == 1 && A == 1) return PROTON;
141  else return 10000*Z + 10*A + 1000000000;
142  }
143 
144  inline int get_particle_Z(int pid) {
145  if (pid == marley_utils::PROTON) return 1;
146  else if (pid == marley_utils::NEUTRON) return 0;
147  // nuclear fragment
148  else if (pid > 1000000000) return (pid % 10000000)/10000;
149  // other particle
150  else return 0;
151  }
152 
153  inline int get_particle_A(int pid) {
154  if (pid == marley_utils::PROTON) return 1;
155  else if (pid == marley_utils::NEUTRON) return 1;
156  // nuclear fragment
157  else if (pid > 1000000000) return (pid % 10000)/10;
158  // other particle
159  else return 0;
160  }
161 
166  bool string_to_neutrino_pdg(const std::string& str, int& pdg);
167 
171  std::string neutrino_pdg_to_string(int pdg);
172 
176  inline bool is_lepton( int pdg ) {
177  int abs_pdg = std::abs( pdg );
178  bool is_a_lepton = ( abs_pdg >= ELECTRON && abs_pdg <= TAU_NEUTRINO );
179  return is_a_lepton;
180  }
181 
187  inline bool is_ion( int pdg ) {
188  bool is_an_ion = ( pdg > 1000000000 && pdg < 2000000000 );
189  return is_an_ion;
190  }
191 
192  // Take the square root of a number. Assume that a negative argument is
193  // due to roundoff error and return zero in such cases rather than NaN.
194  double real_sqrt(double num);
195 
196  // A function template that will raise a number to an integer power.
197  // We can usually use std::pow for this sort of thing (and, unlike this
198  // approach, fractional powers are also supported by that one). However, since
199  // std::pow is not constexpr, this is a workaround.
200  // This function was taken from from https://tinyurl.com/constexpr-pow.
201  template <typename T> constexpr T ipow(T num, unsigned int pow)
202  {
203  return ( pow >= sizeof(unsigned int)*8 ) ? 0 :
204  pow == 0 ? 1 : num * ipow(num, pow - 1);
205  }
206 
207  // Compute the complex gamma function using the Lanczos approximation
208  std::complex<double> gamma(std::complex<double> z);
209 
210  // Numerically integrate a 1D function using Clenshaw-Curtis quadrature
211  double num_integrate(const std::function<double(double)> &f,
212  double a, double b);
213 
214  // Numerically minimize or maximize a function of one variable using
215  // Brent's method (see http://en.wikipedia.org/wiki/Brent%27s_method)
216  double minimize(const std::function<double(double)> f, double leftEnd,
217  double rightEnd, double epsilon, double& minLoc);
218 
219  double maximize(const std::function<double(double)> f, double leftEnd,
220  double rightEnd, double epsilon, double& maxLoc);
221 
222  // Find both solutions of a quadratic equation while attempting
223  // to avoid floating-point arithmetic issues
224  void solve_quadratic_equation(double A, double B,
225  double C, double &solPlus, double &solMinus);
226 
227  // Efficiently read in an entire file as a std::string
228  std::string get_file_contents(std::string filename);
229 
230  // Advance to the next line of an ifstream that either matches (match == true)
231  // or does not match (match == false) a given regular expression
232  std::string get_next_line(std::ifstream &file_in, const std::regex &rx,
233  bool match);
234  // Do the same, but store the number of lines used in num_lines
235  std::string get_next_line(std::ifstream &file_in, const std::regex &rx,
236  bool match, int& num_lines);
237 
238  // String containing all of the characters that will be
239  // considered whitespace by default in the string
240  // manipulation functions below
241  const std::string whitespace = " \f\n\r\t\v";
242 
243  // This version of std::stod will return 0 if it encounters
244  // an empty string or an all-whitespace string.
245  inline double str_to_double(const std::string& s) {
246  size_t endpos = s.find_last_not_of(whitespace);
247  if (endpos == std::string::npos) {
248  return 0.0; // string was all whitespace
249  }
250  else {
251  return std::stod(s);
252  }
253  }
254 
255  // Function that creates a copy of a std::string object
256  // that has been converted to all lowercase
257  inline std::string to_lowercase(const std::string& s) {
258  std::string new_s = s;
259  std::transform(new_s.begin(), new_s.end(), new_s.begin(), ::tolower);
260  return new_s;
261  }
262 
263  // Function that converts a std::string object to
264  // all lowercase in place and returns a reference to
265  // it afterwards
266  inline std::string& to_lowercase_inplace(std::string& s) {
267  std::transform(s.begin(), s.end(), s.begin(), ::tolower);
268  return s;
269  }
270 
271  // Function that converts a std::string object to
272  // all uppercase in place and returns a reference to
273  // it afterwards
274  inline std::string& to_uppercase_inplace(std::string& s) {
275  std::transform(s.begin(), s.end(), s.begin(), ::toupper);
276  return s;
277  }
278 
279  // Functions for padding std::string objects in place. They all return
280  // references to the string afterwards. These functions are based on
281  // http://stackoverflow.com/a/667219/4081973
282  inline std::string& pad_left_inplace(std::string &str,
283  const size_t len, const char pad_char = ' ')
284  {
285  if(len > str.size())
286  str.insert(0, len - str.size(), pad_char);
287  return str;
288  }
289 
290  inline std::string& pad_right_inplace(std::string &str,
291  const size_t len, const char pad_char = ' ')
292  {
293  if(len > str.size())
294  str.append(len - str.size(), pad_char);
295  return str;
296  }
297 
298  // These std::string trimming functions were taken from code
299  // presented here: http://www.cplusplus.com/faq/sequences/strings/trim/
300  // and here: http://stackoverflow.com/questions/216823/whats-the-best-way-to-trim-stdstring
301  // The first three (with the suffix copy) return a trimmed copy of
302  // the string without modifying the original.
303  inline std::string trim_right_copy(const std::string& s,
304  const std::string& delimiters = whitespace)
305  {
306  size_t endpos = s.find_last_not_of(delimiters);
307  return (endpos == std::string::npos) ? "" : s.substr(0, endpos + 1);
308  }
309 
310  inline std::string trim_left_copy(const std::string& s,
311  const std::string& delimiters = whitespace)
312  {
313  size_t startpos = s.find_first_not_of(delimiters);
314  return (startpos == std::string::npos) ? "" : s.substr(startpos);
315  }
316 
317  inline std::string trim_copy(const std::string& s,
318  const std::string& delimiters = whitespace)
319  {
320  return trim_left_copy(trim_right_copy(s, delimiters), delimiters);
321  }
322 
323  // The second three alter the original string, returning a
324  // reference to it after it has been trimmed.
325  inline std::string& trim_right_inplace(std::string& s,
326  const std::string& delimiters = whitespace)
327  {
328  size_t endpos = s.find_last_not_of(delimiters);
329  if (endpos == std::string::npos) {
330  s.clear();
331  }
332  else {
333  s.erase(endpos + 1);
334  }
335  return s;
336  }
337 
338  inline std::string& trim_left_inplace(std::string& s,
339  const std::string& delimiters = whitespace)
340  {
341  size_t startpos = s.find_first_not_of(delimiters);
342  if (startpos == std::string::npos) {
343  s.clear();
344  }
345  else {
346  s.erase(0, startpos);
347  }
348  return s;
349  }
350 
351  inline std::string& trim_inplace(std::string& s,
352  const std::string& delimiters = whitespace)
353  {
354  return trim_left_inplace(trim_right_inplace(s,delimiters), delimiters);
355  }
356 
357  // Split a string into substrings separated by a single-character
358  // delimiter. Return a vector loaded with the resulting array of strings.
359  // Based on http://www.martinbroadhurst.com/how-to-split-a-string-in-c.html
360  inline std::vector<std::string> split_string(const std::string& str,
361  char delim = ' ')
362  {
363  std::vector<std::string> vec;
364  std::stringstream ss( str );
365  std::string token;
366  while ( std::getline(ss, token, delim) ) {
367  vec.push_back( token );
368  }
369  return vec;
370  }
371 
372  // Function that takes a number of bytes and returns a string
373  // representing the amount of memory in more readable units
374  std::string num_bytes_to_string(double bytes, unsigned precision = 3);
375 
376  // Trim an ENSDF nucid string and make two-letter element symbols have a
377  // lowercase last letter. Currently, no checking is done to see if the
378  // string is a valid nucid.
379  std::string nucid_to_symbol(std::string nucid);
380 
381  // Similar to nucid_to_symbol, but returns the atomic number as an integer
382  // instead
383  int nucid_to_Z(std::string nucid);
384 
385  // Generalized std::chrono::duration helper types
386  template <typename repType> using
387  seconds = std::chrono::duration< repType >;
388  template <typename repType> using
389  minutes = std::chrono::duration< repType, std::ratio<60> >;
390  template <typename repType> using
391  hours = std::chrono::duration< repType, std::ratio<3600> >;
392  template <typename repType> using
393  days = std::chrono::duration< repType, std::ratio<86400> >;
394 
395  // This function is a generalized version of code taken from the accepted answer at
396  // http://stackoverflow.com/questions/15957805/extract-year-month-day-etc-from-stdchronotime-point-in-c
397  //
398  // Function that returns a string representation (in the
399  // format days hours:minutes:seconds) of a std::chrono::duration object
400  template <typename repType, typename periodType = std::ratio<1>> std::string duration_to_string(
401  std::chrono::duration<repType, periodType> duration)
402  {
403  int day_count = static_cast<int>(std::chrono::duration_cast
404  <marley_utils::days<repType>>(duration) / (marley_utils::days<repType>(1)));
405  duration -= marley_utils::days<repType>(day_count);
406 
407  int hour_count = static_cast<int>(std::chrono::duration_cast
408  <marley_utils::hours<repType>>(duration) / (marley_utils::hours<repType>(1)));
409  duration -= marley_utils::hours<repType>(hour_count);
410 
411  int minute_count = static_cast<int>(std::chrono::duration_cast
412  <marley_utils::minutes<repType>>(duration) / (marley_utils::minutes<repType>(1)));
413  duration -= marley_utils::minutes<repType>(minute_count);
414 
415  int second_count = static_cast<int>(std::chrono::duration_cast
416  <marley_utils::seconds<repType>>(duration) / (marley_utils::seconds<repType>(1)));
417  duration -= marley_utils::seconds<repType>(second_count);
418 
419  std::ostringstream out;
420  if (day_count > 1) {
421  out << day_count << " days ";
422  }
423  if (day_count == 1) {
424  out << day_count << " day ";
425  }
426  if (hour_count < 10) out << "0";
427  out << hour_count << ":";
428  if (minute_count < 10) out << "0";
429  out << minute_count << ":";
430  if (second_count < 10) out << "0";
431  out << second_count;
432 
433  return out.str();
434  }
435 
436  template <typename durationType> std::string duration_to_string(
437  durationType duration)
438  {
439  return duration_to_string<typename durationType::rep,
440  typename durationType::period>(duration);
441  }
442 
443  // Function that takes two std::system_clock::time_point objects and returns
444  // a string (in the same format as marley_utils::duration_to_string)
445  // representing the time between them
446  std::string elapsed_time_string(
447  std::chrono::system_clock::time_point &start_time,
448  std::chrono::system_clock::time_point &end_time);
449 
450  // Lookup table for particle symbols (keys are PDG particle IDs,
451  // values are symbols).
452  const std::unordered_map<int, std::string> particle_symbols = {
453  { 11, "e" },
454  { 12, "\u03BDe" },
455  { 13, "\u03BC" },
456  { 14, "\u03BD\u03BC" },
457  { 15, "\u03C4" },
458  { 16, "\u03BD\u03C4" },
459  { 22, "\u03B3" },
460  { 2112, "n" },
461  { 2212, "p" },
462  { 1000010020, "d" },
463  { 1000010030, "t" },
464  { 1000020030, "h" },
465  { 1000020040, "\u03B1" },
466  };
467 
468  // Lookup table for particle electric charges (keys are PDG particle IDs,
469  // values are charges expressed as integer multiples of the proton charge).
470  const std::unordered_map<int, int> particle_electric_charges = {
471  { 11, -1 },
472  { 12, 0 },
473  { 13, -1 },
474  { 14, 0 },
475  { 15, -1 },
476  { 16, 0 },
477  { 22, 0 },
478  { 2112, 0 },
479  { 2212, 1 }
480  };
481 
482  // Looks up the electric charge of a particle based on its PDG particle ID
483  inline int get_particle_charge(int pid) {
484  // If a nuclear particle ID is supplied to this function, assume
485  // that it is a bare nucleus, and return its atomic number Z.
486  if (pid > 1000000000) return (pid % 10000000)/10000;
487  // Otherwise, use the lookup table to determine the charge
488  int charge = particle_electric_charges.at( std::abs(pid) );
489  // The lookup table only contains particles (as opposed to antiparticles).
490  // If an antiparticle was requested, return the opposite electric charge.
491  if ( pid < 0 ) charge *= -1;
492  return charge;
493  }
494 
495  inline std::string get_particle_symbol(int pid) {
496  int charge = get_particle_charge( pid );
497  std::string result = particle_symbols.at( std::abs(pid) );
498  if ( charge < 0 ) result += "\u207B";
499  else if ( charge > 0 ) result += "\u207A";
500  else if ( pid < 0 ) result = "anti-" + result;
501  return result;
502  }
503 
504  // Prompt the user with a yes/no question and retrieve the result
505  bool prompt_yes_no(const std::string& message);
506 
507  // Lookup table for element symbols (keys are atomic numbers Z,
508  // values are symbols on the periodic table). The symbol "Nn" is
509  // used for a neutron to match the ENSDF convention.
510  extern const std::unordered_map<int, std::string> element_symbols;
511 
512  // Lookup table for atomic numbers (keys are symbols on the periodic table,
513  // values are atomic numbers Z). The symbol "Nn" is used for a neutron to
514  // match the ENSDF convention.
515  extern const std::unordered_map<std::string, int> atomic_numbers;
516 
517  extern const std::string marley_logo;
518 
519  extern const std::string marley_pic;
520 }