ASPECT
utilities.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2014 - 2024 by the authors of the ASPECT code.
3 
4  This file is part of ASPECT.
5 
6  ASPECT is free software; you can redistribute it and/or modify
7  it under the terms of the GNU General Public License as published by
8  the Free Software Foundation; either version 2, or (at your option)
9  any later version.
10 
11  ASPECT is distributed in the hope that it will be useful,
12  but WITHOUT ANY WARRANTY; without even the implied warranty of
13  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  GNU General Public License for more details.
15 
16  You should have received a copy of the GNU General Public License
17  along with ASPECT; see the file LICENSE. If not see
18  <http://www.gnu.org/licenses/>.
19 */
20 
21 
22 #ifndef _aspect_utilities_h
23 #define _aspect_utilities_h
24 
25 #include <aspect/global.h>
26 
27 #include <array>
28 #include <random>
29 #include <deal.II/base/point.h>
30 #include <deal.II/base/conditional_ostream.h>
31 #include <deal.II/base/table_indices.h>
32 #include <deal.II/base/function_lib.h>
33 #include <deal.II/dofs/dof_handler.h>
34 #include <deal.II/fe/component_mask.h>
35 #include <deal.II/lac/solver_control.h>
36 
38 #include <aspect/structured_data.h>
39 
40 #include <mpi.h>
41 
42 
43 namespace aspect
44 {
45  template <int dim> class SimulatorAccess;
46 
47  namespace GeometryModel
48  {
49  template <int dim> class Interface;
50  }
51 
56  namespace Utilities
57  {
72  using namespace ::Utilities;
73 
74 
88  template <typename T>
89  std::vector<T>
90  possibly_extend_from_1_to_N (const std::vector<T> &values,
91  const unsigned int N,
92  const std::string &id_text);
93 
101  namespace MapParsing
102  {
107  struct Options
108  {
118  std::vector<std::string> list_of_allowed_keys;
119 
127  std::vector<std::string> list_of_required_keys;
128 
134  std::string property_name;
135 
144 
145  /*
146  * Whether to allow for some keys in list_of_required_keys to be
147  * not set to any values, i.e. they do not appear at all.
148  * This also allows a completely empty map.
149  */
151 
161 
171 
180  std::vector<unsigned int> n_values_per_key;
181 
186  Options() = delete;
187 
198  Options(const std::vector<std::string> &list_of_required_keys,
199  const std::string &property_name)
200  :
201  list_of_allowed_keys(list_of_required_keys),
202  list_of_required_keys(list_of_required_keys),
203  property_name(property_name),
204  allow_multiple_values_per_key(false),
205  allow_missing_keys(false),
206  store_values_per_key(false),
207  check_values_per_key(false),
208  n_values_per_key()
209  {}
210  };
211 
237  std::vector<double>
238  parse_map_to_double_array(const std::string &input_string,
239  Options &options);
240  }
241 
315  DEAL_II_DEPRECATED
316  std::vector<double>
317  parse_map_to_double_array (const std::string &key_value_map,
318  const std::vector<std::string> &list_of_keys,
319  const bool expects_background_field,
320  const std::string &property_name,
321  const bool allow_multiple_values_per_key = false,
322  const std::unique_ptr<std::vector<unsigned int>> &n_values_per_key = nullptr,
323  const bool allow_missing_keys = false);
324 
338  template <typename T>
339  Table<2,T>
340  parse_input_table (const std::string &input_string,
341  const unsigned int n_rows,
342  const unsigned int n_columns,
343  const std::string &property_name);
344 
357  template <int dim>
358  std::vector<std::string>
359  expand_dimensional_variable_names (const std::vector<std::string> &var_declarations);
360 
367  template <int dim>
368  IndexSet extract_locally_active_dofs_with_component(const DoFHandler<dim> &dof_handler,
369  const ComponentMask &component_mask);
370 
380  template <int dim>
381  std::vector<Point<dim>> get_unit_support_points(const SimulatorAccess<dim> &simulator_access);
382 
383 
384 
390  template <int dim>
391  bool
392  point_is_in_triangulation(const Mapping<dim> &mapping,
393  const parallel::distributed::Triangulation<dim> &triangulation,
394  const Point<dim> &point,
395  const MPI_Comm mpi_communicator);
396 
397 
398  namespace Coordinates
399  {
400 
406  template <int dim>
407  std::array<double,dim>
408  WGS84_coordinates(const ::Point<dim> &position);
409 
420  template <int dim>
421  std::array<double,dim>
422  cartesian_to_spherical_coordinates(const ::Point<dim> &position);
423 
429  template <int dim>
430  ::Point<dim>
431  spherical_to_cartesian_coordinates(const std::array<double,dim> &scoord);
432 
438  template <int dim>
439  Tensor<1,dim>
440  spherical_to_cartesian_vector(const Tensor<1,dim> &spherical_vector,
441  const ::Point<dim> &position);
442 
443 
449  template <int dim>
450  std::array<double,3>
451  cartesian_to_ellipsoidal_coordinates(const ::Point<3> &position,
452  const double semi_major_axis_a,
453  const double eccentricity);
454 
459  template <int dim>
460  ::Point<3>
461  ellipsoidal_to_cartesian_coordinates(const std::array<double,3> &phi_theta_d,
462  const double semi_major_axis_a,
463  const double eccentricity);
464 
465 
472  string_to_coordinate_system (const std::string &);
473  }
474 
475 
480  template <int dim>
481  bool
482  polygon_contains_point(const std::vector<Point<2>> &point_list,
483  const ::Point<2> &point);
484 
490  template <int dim>
491  double
492  signed_distance_to_polygon(const std::vector<Point<2>> &point_list,
493  const ::Point<2> &point);
494 
495 
502  double
503  distance_to_line(const std::array<::Point<2>,2> &point_list,
504  const ::Point<2> &point);
505 
513  template <int dim>
514  std::array<Tensor<1,dim>,dim-1>
515  orthogonal_vectors (const Tensor<1,dim> &v);
516 
521  Tensor<2,3>
522  rotation_matrix_from_axis (const Tensor<1,3> &rotation_axis,
523  const double rotation_angle);
524 
531  Tensor<2,3>
532  compute_rotation_matrix_for_slice (const Tensor<1,3> &point_one,
533  const Tensor<1,3> &point_two);
534 
571  std::pair<double,double> real_spherical_harmonic( unsigned int l, // degree
572  unsigned int m, // order
573  double theta, // colatitude (radians)
574  double phi ); // longitude (radians)
575 
579  struct ThousandSep : std::numpunct<char>
580  {
581  protected:
582  char do_thousands_sep() const override
583  {
584  return ',';
585  }
586 
587  std::string do_grouping() const override
588  {
589  return "\003"; // groups of 3 digits (this string is in octal format)
590  }
591 
592  };
593 
606  bool fexists(const std::string &filename);
607 
621  bool fexists(const std::string &filename,
622  const MPI_Comm comm);
623 
629  bool filename_is_url(const std::string &filename);
630 
650  std::string
651  read_and_distribute_file_content(const std::string &filename,
652  const MPI_Comm comm);
653 
666  void
667  collect_and_write_file_content(const std::string &filename,
668  const std::string &file_content,
669  const MPI_Comm comm);
670 
684  int
685  mkdirp(std::string pathname, const mode_t mode = 0755);
686 
698  void create_directory(const std::string &pathname,
699  const MPI_Comm comm,
700  const bool silent);
701 
706  namespace tk
707  {
711  class spline
712  {
713  public:
723  void set_points(const std::vector<double> &x,
724  const std::vector<double> &y,
725  const bool cubic_spline = true,
726  const bool monotone_spline = false);
730  double operator() (double x) const;
731 
732  private:
736  std::vector<double> m_x;
737 
744  std::vector<double> m_a, m_b, m_c, m_y;
745  };
746  }
747 
755  inline
756  void
757  extract_composition_values_at_q_point (const std::vector<std::vector<double>> &composition_values,
758  const unsigned int q,
759  std::vector<double> &composition_values_at_q_point);
760 
765  std::string
766  expand_ASPECT_SOURCE_DIR (const std::string &location);
767 
772  std::string parenthesize_if_nonempty (const std::string &s);
773 
777  bool
778  string_to_bool(const std::string &s);
779 
783  std::vector<bool>
784  string_to_bool(const std::vector<std::string> &s);
785 
789  unsigned int
790  string_to_unsigned_int(const std::string &s);
791 
795  std::vector<unsigned int>
796  string_to_unsigned_int(const std::vector<std::string> &s);
797 
802  bool has_unique_entries (const std::vector<std::string> &strings);
803 
804 
805 
827  double weighted_p_norm_average (const std::vector<double> &weights,
828  const std::vector<double> &values,
829  const double p);
830 
831 
857  template <typename T>
858  T derivative_of_weighted_p_norm_average (const double averaged_parameter,
859  const std::vector<double> &weights,
860  const std::vector<double> &values,
861  const std::vector<T> &derivatives,
862  const double p);
894  template <int dim>
895  double compute_spd_factor(const double eta,
896  const SymmetricTensor<2,dim> &strain_rate,
897  const SymmetricTensor<2,dim> &dviscosities_dstrain_rate,
898  const double SPD_safety_factor);
899 
903  template <int dim>
904  Point<dim> convert_array_to_point(const std::array<double,dim> &array);
905 
909  template <int dim>
910  std::array<double,dim> convert_point_to_array(const Point<dim> &point);
911 
919  class Operator
920  {
921  public:
926  {
932  replace_if_valid
933  };
934 
939  Operator();
940 
944  Operator(const operation op);
945 
950  double operator() (const double x, const double y) const;
951 
956  bool operator== (const operation op) const;
957 
958  private:
963  };
964 
969  std::vector<Operator> create_model_operator_list(const std::vector<std::string> &operator_names);
970 
974  const std::string get_model_operator_options();
975 
981  template <int dim>
982  SymmetricTensor<2,dim> nth_basis_for_symmetric_tensors (const unsigned int k);
983 
987  template <int dim>
989  {
990  public:
994  NaturalCoordinate(Point<dim> &position,
995  const GeometryModel::Interface<dim> &geometry_model);
996 
1001  NaturalCoordinate(const std::array<double, dim> &coord,
1002  const Utilities::Coordinates::CoordinateSystem &coord_system);
1003 
1008  std::array<double,dim> &get_coordinates();
1009 
1014  const std::array<double,dim> &get_coordinates() const;
1015 
1020  std::array<double,dim-1> get_surface_coordinates() const;
1021 
1026  double get_depth_coordinate() const;
1027 
1028  private:
1034 
1038  std::array<double,dim> coordinates;
1039  };
1040 
1041 
1060  template <int dim, typename VectorType>
1061  void
1062  project_cellwise(const Mapping<dim> &mapping,
1063  const DoFHandler<dim> &dof_handler,
1064  const unsigned int component_index,
1065  const Quadrature<dim> &quadrature,
1066  const std::function<void(
1067  const typename DoFHandler<dim>::active_cell_iterator &,
1068  const std::vector<Point<dim>> &,
1069  std::vector<double> &)> &function,
1070  VectorType &vec_result);
1071 
1095  void throw_linear_solver_failure_exception(const std::string &solver_name,
1096  const std::string &function_name,
1097  const std::vector<SolverControl> &solver_controls,
1098  const std::exception &exc,
1099  const MPI_Comm mpi_communicator,
1100  const std::string &output_filename = "");
1101 
1112  template <int dim>
1114  {
1115  public:
1125  VectorFunctionFromVelocityFunctionObject (const unsigned int n_components,
1126  const std::function<Tensor<1,dim> (const Point<dim> &)> &function_object);
1127  };
1128 
1135  template <typename T>
1136  inline
1137  std::vector<std::size_t>
1138  compute_sorting_permutation(const std::vector<T> &vector);
1139 
1147  template <typename T>
1148  inline
1149  std::vector<T>
1151  const std::vector<T> &vector,
1152  const std::vector<std::size_t> &permutation_vector);
1153 
1166  std::vector<Tensor<2,3>>
1167  rotation_matrices_random_draw_volume_weighting(const std::vector<double> &volume_fractions,
1168  const std::vector<Tensor<2,3>> &rotation_matrices,
1169  const unsigned int n_output_matrices,
1170  std::mt19937 &random_number_generator);
1171 
1175  double wrap_angle(const double angle);
1176 
1181  std::array<double,3> zxz_euler_angles_from_rotation_matrix(const Tensor<2,3> &rotation_matrix);
1182 
1187  Tensor<2,3> zxz_euler_angles_to_rotation_matrix(const double phi1,
1188  const double theta,
1189  const double phi2);
1190 
1191  }
1192 
1193 
1194 // inline implementations:
1195 #ifndef DOXYGEN
1196  namespace Utilities
1197  {
1198 
1199  template <typename T>
1200  inline
1201  std::vector<T>
1202  possibly_extend_from_1_to_N (const std::vector<T> &values,
1203  const unsigned int N,
1204  const std::string &id_text)
1205  {
1206  if (values.size() == 1)
1207  {
1208  return std::vector<T> (N, values[0]);
1209  }
1210  else if (values.size() == N)
1211  {
1212  return values;
1213  }
1214  else
1215  {
1216  // Non-specified behavior
1217  AssertThrow(false,
1218  ExcMessage("Length of " + id_text + " list must be " +
1219  "either one or " + Utilities::to_string(N) +
1220  ". Currently it is " + Utilities::to_string(values.size()) + "."));
1221  }
1222 
1223  // This should never happen, but return an empty vector so the compiler
1224  // will be happy
1225  return std::vector<T> ();
1226  }
1227 
1228  inline
1229  void
1230  extract_composition_values_at_q_point (const std::vector<std::vector<double>> &composition_values,
1231  const unsigned int q,
1232  std::vector<double> &composition_values_at_q_point)
1233  {
1234  Assert(q<composition_values.size(), ExcInternalError());
1235  Assert(composition_values_at_q_point.size() > 0,
1236  ExcInternalError());
1237 
1238  for (unsigned int k=0; k < composition_values_at_q_point.size(); ++k)
1239  {
1240  Assert(composition_values[k].size() == composition_values_at_q_point.size(),
1241  ExcInternalError());
1242  composition_values_at_q_point[k] = composition_values[k][q];
1243  }
1244  }
1245 
1246  template <typename T>
1247  inline
1248  std::vector<std::size_t>
1249  compute_sorting_permutation(const std::vector<T> &vector)
1250  {
1251  std::vector<std::size_t> p(vector.size());
1252  std::iota(p.begin(), p.end(), 0);
1253  std::sort(p.begin(), p.end(),
1254  [&](std::size_t i, std::size_t j)
1255  {
1256  return vector[i] < vector[j];
1257  });
1258  return p;
1259  }
1260 
1261  template <typename T>
1262  inline
1263  std::vector<T>
1265  const std::vector<T> &vector,
1266  const std::vector<std::size_t> &permutation_vector)
1267  {
1268  std::vector<T> sorted_vec(vector.size());
1269  std::transform(permutation_vector.begin(), permutation_vector.end(), sorted_vec.begin(),
1270  [&](std::size_t i)
1271  {
1272  return vector[i];
1273  });
1274  return sorted_vec;
1275  }
1276 
1280  namespace Tensors
1281  {
1287  template <int dim, class Iterator>
1288  inline
1289  SymmetricTensor<2,dim>
1290  to_symmetric_tensor(const Iterator begin,
1291  const Iterator end)
1292  {
1293  AssertDimension(std::distance(begin, end), (SymmetricTensor<2,dim>::n_independent_components));
1294  (void) end;
1295 
1296  SymmetricTensor<2,dim> output;
1297 
1298  Iterator next = begin;
1299  for (unsigned int i=0; i < SymmetricTensor<2,dim>::n_independent_components; ++i, ++next)
1300  output[SymmetricTensor<2,dim>::unrolled_to_component_indices(i)] = *next;
1301 
1302  return output;
1303  }
1304 
1309  template <int dim, class Iterator>
1310  inline
1311  void
1312  unroll_symmetric_tensor_into_array(const SymmetricTensor<2,dim> &tensor,
1313  const Iterator begin,
1314  const Iterator end)
1315  {
1316  AssertDimension(std::distance(begin, end), (SymmetricTensor<2,dim>::n_independent_components));
1317  (void) end;
1318 
1319  Iterator next = begin;
1320  for (unsigned int i=0; i < SymmetricTensor<2,dim>::n_independent_components; ++i, ++next)
1321  *next = tensor[SymmetricTensor<2,dim>::unrolled_to_component_indices(i)];
1322  }
1323 
1327  SymmetricTensor<4,3>
1328  rotate_full_stiffness_tensor(const Tensor<2,3> &rotation_tensor, const SymmetricTensor<4,3> &input_tensor);
1329 
1334  SymmetricTensor<2,6>
1335  rotate_voigt_stiffness_matrix(const Tensor<2,3> &rotation_tensor, const SymmetricTensor<2,6> &input_tensor);
1336 
1341  SymmetricTensor<2,6>
1342  to_voigt_stiffness_matrix(const SymmetricTensor<4,3> &input_tensor);
1343 
1348  SymmetricTensor<4,3>
1349  to_full_stiffness_tensor(const SymmetricTensor<2,6> &input_tensor);
1350 
1355  Tensor<1,21>
1356  to_voigt_stiffness_vector(const SymmetricTensor<2,6> &input_tensor);
1357 
1362  SymmetricTensor<2,6>
1363  to_voigt_stiffness_matrix(const Tensor<1,21> &input_tensor);
1364 
1369  Tensor<1,21>
1370  to_voigt_stiffness_vector(const SymmetricTensor<4,3> &input);
1371 
1372 
1377  template <int dim>
1378  const Tensor<dim,dim> &levi_civita();
1379 
1380  // Declare the existence of a specialization:
1381  template <>
1382  const Tensor<3,3> &levi_civita<3>();
1383 
1395  template <int dim>
1396  SymmetricTensor<2,dim>
1397  consistent_deviator(const SymmetricTensor<2,dim> &input);
1398 
1412  template <int dim>
1413  double
1414  consistent_second_invariant_of_deviatoric_tensor(const SymmetricTensor<2,dim> &input);
1415  }
1416 
1417  }
1418 #endif
1419 }
1420 
1421 #endif
T derivative_of_weighted_p_norm_average(const double averaged_parameter, const std::vector< double > &weights, const std::vector< double > &values, const std::vector< T > &derivatives, const double p)
std::string read_and_distribute_file_content(const std::string &filename, const MPI_Comm comm)
Point< dim > convert_array_to_point(const std::array< double, dim > &array)
std::vector< double > m_y
Definition: utilities.h:744
Tensor< 1, dim > spherical_to_cartesian_vector(const Tensor< 1, dim > &spherical_vector, const ::Point< dim > &position)
Tensor< 2, 3 > compute_rotation_matrix_for_slice(const Tensor< 1, 3 > &point_one, const Tensor< 1, 3 > &point_two)
std::vector< std::string > list_of_required_keys
Definition: utilities.h:127
double weighted_p_norm_average(const std::vector< double > &weights, const std::vector< double > &values, const double p)
double distance_to_line(const std::array<::Point< 2 >, 2 > &point_list, const ::Point< 2 > &point)
std::vector< std::string > expand_dimensional_variable_names(const std::vector< std::string > &var_declarations)
std::array< double, dim > coordinates
Definition: utilities.h:1038
void collect_and_write_file_content(const std::string &filename, const std::string &file_content, const MPI_Comm comm)
double compute_spd_factor(const double eta, const SymmetricTensor< 2, dim > &strain_rate, const SymmetricTensor< 2, dim > &dviscosities_dstrain_rate, const double SPD_safety_factor)
Table< 2, T > parse_input_table(const std::string &input_string, const unsigned int n_rows, const unsigned int n_columns, const std::string &property_name)
std::pair< double, double > real_spherical_harmonic(unsigned int l, unsigned int m, double theta, double phi)
std::vector< double > m_x
Definition: utilities.h:736
std::string expand_ASPECT_SOURCE_DIR(const std::string &location)
Tensor< 2, 3 > zxz_euler_angles_to_rotation_matrix(const double phi1, const double theta, const double phi2)
std::array< Tensor< 1, dim >, dim-1 > orthogonal_vectors(const Tensor< 1, dim > &v)
std::vector< Tensor< 2, 3 > > rotation_matrices_random_draw_volume_weighting(const std::vector< double > &volume_fractions, const std::vector< Tensor< 2, 3 >> &rotation_matrices, const unsigned int n_output_matrices, std::mt19937 &random_number_generator)
std::string to_string(const Newton::Parameters::Stabilization preconditioner_stabilization)
std::array< double, dim > WGS84_coordinates(const ::Point< dim > &position)
::Point< 3 > ellipsoidal_to_cartesian_coordinates(const std::array< double, 3 > &phi_theta_d, const double semi_major_axis_a, const double eccentricity)
void project_cellwise(const Mapping< dim > &mapping, const DoFHandler< dim > &dof_handler, const unsigned int component_index, const Quadrature< dim > &quadrature, const std::function< void(const typename DoFHandler< dim >::active_cell_iterator &, const std::vector< Point< dim >> &, std::vector< double > &)> &function, VectorType &vec_result)
bool point_is_in_triangulation(const Mapping< dim > &mapping, const parallel::distributed::Triangulation< dim > &triangulation, const Point< dim > &point, const MPI_Comm mpi_communicator)
void throw_linear_solver_failure_exception(const std::string &solver_name, const std::string &function_name, const std::vector< SolverControl > &solver_controls, const std::exception &exc, const MPI_Comm mpi_communicator, const std::string &output_filename="")
bool polygon_contains_point(const std::vector< Point< 2 >> &point_list, const ::Point< 2 > &point)
int mkdirp(std::string pathname, const mode_t mode=0755)
Utilities::Coordinates::CoordinateSystem coordinate_system
Definition: utilities.h:1033
void create_directory(const std::string &pathname, const MPI_Comm comm, const bool silent)
::Point< dim > spherical_to_cartesian_coordinates(const std::array< double, dim > &scoord)
std::string do_grouping() const override
Definition: utilities.h:587
std::vector< T > possibly_extend_from_1_to_N(const std::vector< T > &values, const unsigned int N, const std::string &id_text)
std::array< double, 3 > cartesian_to_ellipsoidal_coordinates(const ::Point< 3 > &position, const double semi_major_axis_a, const double eccentricity)
Options(const std::vector< std::string > &list_of_required_keys, const std::string &property_name)
Definition: utilities.h:198
std::array< double, dim > convert_point_to_array(const Point< dim > &point)
std::vector< bool > string_to_bool(const std::vector< std::string > &s)
std::vector< Point< dim > > get_unit_support_points(const SimulatorAccess< dim > &simulator_access)
std::vector< T > apply_permutation(const std::vector< T > &vector, const std::vector< std::size_t > &permutation_vector)
IndexSet extract_locally_active_dofs_with_component(const DoFHandler< dim > &dof_handler, const ComponentMask &component_mask)
std::vector< unsigned int > n_values_per_key
Definition: utilities.h:180
std::vector< unsigned int > string_to_unsigned_int(const std::vector< std::string > &s)
DEAL_II_DEPRECATED std::vector< double > parse_map_to_double_array(const std::string &key_value_map, const std::vector< std::string > &list_of_keys, const bool expects_background_field, const std::string &property_name, const bool allow_multiple_values_per_key=false, const std::unique_ptr< std::vector< unsigned int >> &n_values_per_key=nullptr, const bool allow_missing_keys=false)
double wrap_angle(const double angle)
bool fexists(const std::string &filename, const MPI_Comm comm)
std::vector< Operator > create_model_operator_list(const std::vector< std::string > &operator_names)
double signed_distance_to_polygon(const std::vector< Point< 2 >> &point_list, const ::Point< 2 > &point)
std::string parenthesize_if_nonempty(const std::string &s)
std::array< double, 3 > zxz_euler_angles_from_rotation_matrix(const Tensor< 2, 3 > &rotation_matrix)
void extract_composition_values_at_q_point(const std::vector< std::vector< double >> &composition_values, const unsigned int q, std::vector< double > &composition_values_at_q_point)
const std::string get_model_operator_options()
char do_thousands_sep() const override
Definition: utilities.h:582
std::vector< std::string > list_of_allowed_keys
Definition: utilities.h:118
bool filename_is_url(const std::string &filename)
SymmetricTensor< 2, dim > nth_basis_for_symmetric_tensors(const unsigned int k)
std::vector< std::size_t > compute_sorting_permutation(const std::vector< T > &vector)
CoordinateSystem string_to_coordinate_system(const std::string &)
std::array< double, dim > cartesian_to_spherical_coordinates(const ::Point< dim > &position)
bool has_unique_entries(const std::vector< std::string > &strings)
Tensor< 2, 3 > rotation_matrix_from_axis(const Tensor< 1, 3 > &rotation_axis, const double rotation_angle)