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 
37 #include <aspect/structured_data.h>
38 
39 
40 
41 namespace aspect
42 {
43  template <int dim> class SimulatorAccess;
44 
45  namespace GeometryModel
46  {
47  template <int dim> class Interface;
48  }
49 
54  namespace Utilities
55  {
70  using namespace ::Utilities;
71 
72 
86  template <typename T>
87  std::vector<T>
88  possibly_extend_from_1_to_N (const std::vector<T> &values,
89  const unsigned int N,
90  const std::string &id_text);
91 
99  namespace MapParsing
100  {
105  struct Options
106  {
116  std::vector<std::string> list_of_allowed_keys;
117 
125  std::vector<std::string> list_of_required_keys;
126 
132  std::string property_name;
133 
142 
143  /*
144  * Whether to allow for some keys in list_of_required_keys to be
145  * not set to any values, i.e. they do not appear at all.
146  * This also allows a completely empty map.
147  */
149 
159 
169 
178  std::vector<unsigned int> n_values_per_key;
179 
184  Options() = delete;
185 
196  Options(const std::vector<std::string> &list_of_required_keys,
197  const std::string &property_name)
198  :
199  list_of_allowed_keys(list_of_required_keys),
200  list_of_required_keys(list_of_required_keys),
201  property_name(property_name),
202  allow_multiple_values_per_key(false),
203  allow_missing_keys(false),
204  store_values_per_key(false),
205  check_values_per_key(false),
206  n_values_per_key()
207  {}
208  };
209 
235  std::vector<double>
236  parse_map_to_double_array(const std::string &input_string,
237  Options &options);
238  }
239 
313  DEAL_II_DEPRECATED
314  std::vector<double>
315  parse_map_to_double_array (const std::string &key_value_map,
316  const std::vector<std::string> &list_of_keys,
317  const bool expects_background_field,
318  const std::string &property_name,
319  const bool allow_multiple_values_per_key = false,
320  const std::unique_ptr<std::vector<unsigned int>> &n_values_per_key = nullptr,
321  const bool allow_missing_keys = false);
322 
336  template <typename T>
337  Table<2,T>
338  parse_input_table (const std::string &input_string,
339  const unsigned int n_rows,
340  const unsigned int n_columns,
341  const std::string &property_name);
342 
355  template <int dim>
356  std::vector<std::string>
357  expand_dimensional_variable_names (const std::vector<std::string> &var_declarations);
358 
365  template <int dim>
366  IndexSet extract_locally_active_dofs_with_component(const DoFHandler<dim> &dof_handler,
367  const ComponentMask &component_mask);
368 
378  template <int dim>
379  std::vector<Point<dim>> get_unit_support_points(const SimulatorAccess<dim> &simulator_access);
380 
381 
382 
388  template <int dim>
389  bool
390  point_is_in_triangulation(const Mapping<dim> &mapping,
391  const parallel::distributed::Triangulation<dim> &triangulation,
392  const Point<dim> &point,
393  const MPI_Comm mpi_communicator);
394 
395 
396  namespace Coordinates
397  {
398 
404  template <int dim>
405  std::array<double,dim>
406  WGS84_coordinates(const ::Point<dim> &position);
407 
418  template <int dim>
419  std::array<double,dim>
420  cartesian_to_spherical_coordinates(const ::Point<dim> &position);
421 
427  template <int dim>
428  ::Point<dim>
429  spherical_to_cartesian_coordinates(const std::array<double,dim> &scoord);
430 
436  template <int dim>
437  Tensor<1,dim>
438  spherical_to_cartesian_vector(const Tensor<1,dim> &spherical_vector,
439  const ::Point<dim> &position);
440 
441 
447  template <int dim>
448  std::array<double,3>
449  cartesian_to_ellipsoidal_coordinates(const ::Point<3> &position,
450  const double semi_major_axis_a,
451  const double eccentricity);
452 
457  template <int dim>
458  ::Point<3>
459  ellipsoidal_to_cartesian_coordinates(const std::array<double,3> &phi_theta_d,
460  const double semi_major_axis_a,
461  const double eccentricity);
462 
463 
470  string_to_coordinate_system (const std::string &);
471  }
472 
473 
478  template <int dim>
479  bool
480  polygon_contains_point(const std::vector<Point<2>> &point_list,
481  const ::Point<2> &point);
482 
488  template <int dim>
489  double
490  signed_distance_to_polygon(const std::vector<Point<2>> &point_list,
491  const ::Point<2> &point);
492 
493 
500  double
501  distance_to_line(const std::array<::Point<2>,2> &point_list,
502  const ::Point<2> &point);
503 
511  template <int dim>
512  std::array<Tensor<1,dim>,dim-1>
513  orthogonal_vectors (const Tensor<1,dim> &v);
514 
519  Tensor<2,3>
520  rotation_matrix_from_axis (const Tensor<1,3> &rotation_axis,
521  const double rotation_angle);
522 
529  Tensor<2,3>
530  compute_rotation_matrix_for_slice (const Tensor<1,3> &point_one,
531  const Tensor<1,3> &point_two);
532 
569  std::pair<double,double> real_spherical_harmonic( unsigned int l, // degree
570  unsigned int m, // order
571  double theta, // colatitude (radians)
572  double phi ); // longitude (radians)
573 
577  struct ThousandSep : std::numpunct<char>
578  {
579  protected:
580  char do_thousands_sep() const override
581  {
582  return ',';
583  }
584 
585  std::string do_grouping() const override
586  {
587  return "\003"; // groups of 3 digits (this string is in octal format)
588  }
589 
590  };
591 
604  bool fexists(const std::string &filename);
605 
619  bool fexists(const std::string &filename,
620  const MPI_Comm comm);
621 
627  bool filename_is_url(const std::string &filename);
628 
648  std::string
649  read_and_distribute_file_content(const std::string &filename,
650  const MPI_Comm comm);
651 
664  void
665  collect_and_write_file_content(const std::string &filename,
666  const std::string &file_content,
667  const MPI_Comm comm);
668 
682  int
683  mkdirp(std::string pathname, const mode_t mode = 0755);
684 
695  void create_directory(const std::string &pathname,
696  const MPI_Comm comm,
697  const bool silent);
698 
703  namespace tk
704  {
708  class spline
709  {
710  public:
720  void set_points(const std::vector<double> &x,
721  const std::vector<double> &y,
722  const bool cubic_spline = true,
723  const bool monotone_spline = false);
727  double operator() (double x) const;
728 
729  private:
733  std::vector<double> m_x;
734 
741  std::vector<double> m_a, m_b, m_c, m_y;
742  };
743  }
744 
752  inline
753  void
754  extract_composition_values_at_q_point (const std::vector<std::vector<double>> &composition_values,
755  const unsigned int q,
756  std::vector<double> &composition_values_at_q_point);
757 
762  std::string
763  expand_ASPECT_SOURCE_DIR (const std::string &location);
764 
769  std::string parenthesize_if_nonempty (const std::string &s);
770 
774  bool
775  string_to_bool(const std::string &s);
776 
780  std::vector<bool>
781  string_to_bool(const std::vector<std::string> &s);
782 
786  unsigned int
787  string_to_unsigned_int(const std::string &s);
788 
792  std::vector<unsigned int>
793  string_to_unsigned_int(const std::vector<std::string> &s);
794 
799  bool has_unique_entries (const std::vector<std::string> &strings);
800 
801 
802 
824  double weighted_p_norm_average (const std::vector<double> &weights,
825  const std::vector<double> &values,
826  const double p);
827 
828 
854  template <typename T>
855  T derivative_of_weighted_p_norm_average (const double averaged_parameter,
856  const std::vector<double> &weights,
857  const std::vector<double> &values,
858  const std::vector<T> &derivatives,
859  const double p);
891  template <int dim>
892  double compute_spd_factor(const double eta,
893  const SymmetricTensor<2,dim> &strain_rate,
894  const SymmetricTensor<2,dim> &dviscosities_dstrain_rate,
895  const double SPD_safety_factor);
896 
900  template <int dim>
901  Point<dim> convert_array_to_point(const std::array<double,dim> &array);
902 
906  template <int dim>
907  std::array<double,dim> convert_point_to_array(const Point<dim> &point);
908 
916  class Operator
917  {
918  public:
923  {
929  replace_if_valid
930  };
931 
936  Operator();
937 
941  Operator(const operation op);
942 
947  double operator() (const double x, const double y) const;
948 
953  bool operator== (const operation op) const;
954 
955  private:
960  };
961 
966  std::vector<Operator> create_model_operator_list(const std::vector<std::string> &operator_names);
967 
971  const std::string get_model_operator_options();
972 
978  template <int dim>
979  SymmetricTensor<2,dim> nth_basis_for_symmetric_tensors (const unsigned int k);
980 
984  template <int dim>
986  {
987  public:
991  NaturalCoordinate(Point<dim> &position,
992  const GeometryModel::Interface<dim> &geometry_model);
993 
998  NaturalCoordinate(const std::array<double, dim> &coord,
999  const Utilities::Coordinates::CoordinateSystem &coord_system);
1000 
1005  std::array<double,dim> &get_coordinates();
1006 
1011  const std::array<double,dim> &get_coordinates() const;
1012 
1017  std::array<double,dim-1> get_surface_coordinates() const;
1018 
1023  double get_depth_coordinate() const;
1024 
1025  private:
1031 
1035  std::array<double,dim> coordinates;
1036  };
1037 
1038 
1057  template <int dim, typename VectorType>
1058  void
1059  project_cellwise(const Mapping<dim> &mapping,
1060  const DoFHandler<dim> &dof_handler,
1061  const unsigned int component_index,
1062  const Quadrature<dim> &quadrature,
1063  const std::function<void(
1064  const typename DoFHandler<dim>::active_cell_iterator &,
1065  const std::vector<Point<dim>> &,
1066  std::vector<double> &)> &function,
1067  VectorType &vec_result);
1068 
1092  void throw_linear_solver_failure_exception(const std::string &solver_name,
1093  const std::string &function_name,
1094  const std::vector<SolverControl> &solver_controls,
1095  const std::exception &exc,
1096  const MPI_Comm mpi_communicator,
1097  const std::string &output_filename = "");
1098 
1106  template <int dim>
1108  {
1109  public:
1119  VectorFunctionFromVelocityFunctionObject (const unsigned int n_components,
1120  const std::function<Tensor<1,dim> (const Point<dim> &)> &function_object);
1121 
1129  double value (const Point<dim> &p,
1130  const unsigned int component = 0) const override;
1131 
1141  void vector_value (const Point<dim> &p,
1142  Vector<double> &values) const override;
1143 
1144  private:
1149  const std::function<Tensor<1,dim> (const Point<dim> &)> function_object;
1150  };
1151 
1158  template <typename T>
1159  inline
1160  std::vector<std::size_t>
1161  compute_sorting_permutation(const std::vector<T> &vector);
1162 
1170  template <typename T>
1171  inline
1172  std::vector<T>
1174  const std::vector<T> &vector,
1175  const std::vector<std::size_t> &permutation_vector);
1176 
1189  std::vector<Tensor<2,3>>
1190  rotation_matrices_random_draw_volume_weighting(const std::vector<double> &volume_fractions,
1191  const std::vector<Tensor<2,3>> &rotation_matrices,
1192  const unsigned int n_output_matrices,
1193  std::mt19937 &random_number_generator);
1194 
1198  double wrap_angle(const double angle);
1199 
1204  std::array<double,3> zxz_euler_angles_from_rotation_matrix(const Tensor<2,3> &rotation_matrix);
1205 
1210  Tensor<2,3> zxz_euler_angles_to_rotation_matrix(const double phi1,
1211  const double theta,
1212  const double phi2);
1213 
1214  }
1215 }
1216 
1217 
1218 // inline implementations:
1219 #ifndef DOXYGEN
1220 namespace aspect
1221 {
1222  namespace Utilities
1223  {
1224 
1225  template <typename T>
1226  inline
1227  std::vector<T>
1228  possibly_extend_from_1_to_N (const std::vector<T> &values,
1229  const unsigned int N,
1230  const std::string &id_text)
1231  {
1232  if (values.size() == 1)
1233  {
1234  return std::vector<T> (N, values[0]);
1235  }
1236  else if (values.size() == N)
1237  {
1238  return values;
1239  }
1240  else
1241  {
1242  // Non-specified behavior
1243  AssertThrow(false,
1244  ExcMessage("Length of " + id_text + " list must be " +
1245  "either one or " + Utilities::to_string(N) +
1246  ". Currently it is " + Utilities::to_string(values.size()) + "."));
1247  }
1248 
1249  // This should never happen, but return an empty vector so the compiler
1250  // will be happy
1251  return std::vector<T> ();
1252  }
1253 
1254  inline
1255  void
1256  extract_composition_values_at_q_point (const std::vector<std::vector<double>> &composition_values,
1257  const unsigned int q,
1258  std::vector<double> &composition_values_at_q_point)
1259  {
1260  Assert(q<composition_values.size(), ExcInternalError());
1261  Assert(composition_values_at_q_point.size() > 0,
1262  ExcInternalError());
1263 
1264  for (unsigned int k=0; k < composition_values_at_q_point.size(); ++k)
1265  {
1266  Assert(composition_values[k].size() == composition_values_at_q_point.size(),
1267  ExcInternalError());
1268  composition_values_at_q_point[k] = composition_values[k][q];
1269  }
1270  }
1271 
1272  template <typename T>
1273  inline
1274  std::vector<std::size_t>
1275  compute_sorting_permutation(const std::vector<T> &vector)
1276  {
1277  std::vector<std::size_t> p(vector.size());
1278  std::iota(p.begin(), p.end(), 0);
1279  std::sort(p.begin(), p.end(),
1280  [&](std::size_t i, std::size_t j)
1281  {
1282  return vector[i] < vector[j];
1283  });
1284  return p;
1285  }
1286 
1287  template <typename T>
1288  inline
1289  std::vector<T>
1291  const std::vector<T> &vector,
1292  const std::vector<std::size_t> &permutation_vector)
1293  {
1294  std::vector<T> sorted_vec(vector.size());
1295  std::transform(permutation_vector.begin(), permutation_vector.end(), sorted_vec.begin(),
1296  [&](std::size_t i)
1297  {
1298  return vector[i];
1299  });
1300  return sorted_vec;
1301  }
1302 
1306  namespace Tensors
1307  {
1313  template <int dim, class Iterator>
1314  inline
1315  SymmetricTensor<2,dim>
1316  to_symmetric_tensor(const Iterator begin,
1317  const Iterator end)
1318  {
1319  AssertDimension(std::distance(begin, end), (SymmetricTensor<2,dim>::n_independent_components));
1320  (void) end;
1321 
1322  SymmetricTensor<2,dim> output;
1323 
1324  Iterator next = begin;
1325  for (unsigned int i=0; i < SymmetricTensor<2,dim>::n_independent_components; ++i, ++next)
1326  output[SymmetricTensor<2,dim>::unrolled_to_component_indices(i)] = *next;
1327 
1328  return output;
1329  }
1330 
1335  template <int dim, class Iterator>
1336  inline
1337  void
1338  unroll_symmetric_tensor_into_array(const SymmetricTensor<2,dim> &tensor,
1339  const Iterator begin,
1340  const Iterator end)
1341  {
1342  AssertDimension(std::distance(begin, end), (SymmetricTensor<2,dim>::n_independent_components));
1343  (void) end;
1344 
1345  Iterator next = begin;
1346  for (unsigned int i=0; i < SymmetricTensor<2,dim>::n_independent_components; ++i, ++next)
1347  *next = tensor[SymmetricTensor<2,dim>::unrolled_to_component_indices(i)];
1348  }
1349 
1353  SymmetricTensor<4,3>
1354  rotate_full_stiffness_tensor(const Tensor<2,3> &rotation_tensor, const SymmetricTensor<4,3> &input_tensor);
1355 
1360  SymmetricTensor<2,6>
1361  rotate_voigt_stiffness_matrix(const Tensor<2,3> &rotation_tensor, const SymmetricTensor<2,6> &input_tensor);
1362 
1367  SymmetricTensor<2,6>
1368  to_voigt_stiffness_matrix(const SymmetricTensor<4,3> &input_tensor);
1369 
1374  SymmetricTensor<4,3>
1375  to_full_stiffness_tensor(const SymmetricTensor<2,6> &input_tensor);
1376 
1381  Tensor<1,21>
1382  to_voigt_stiffness_vector(const SymmetricTensor<2,6> &input_tensor);
1383 
1388  SymmetricTensor<2,6>
1389  to_voigt_stiffness_matrix(const Tensor<1,21> &input_tensor);
1390 
1395  Tensor<1,21>
1396  to_voigt_stiffness_vector(const SymmetricTensor<4,3> &input);
1397 
1398 
1403  template <int dim>
1404  const Tensor<dim,dim> &levi_civita();
1405 
1406  // Declare the existence of a specialization:
1407  template <>
1408  const Tensor<3,3> &levi_civita<3>();
1409  }
1410 
1411  }
1412 }
1413 #endif
1414 
1415 #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:741
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:125
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:1035
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:733
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:1030
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:585
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:196
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:178
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)
Definition: compat.h:59
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:580
std::vector< std::string > list_of_allowed_keys
Definition: utilities.h:116
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)