22 #ifndef _aspect_utilities_h 23 #define _aspect_utilities_h 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> 45 template <
int dim>
class SimulatorAccess;
47 namespace GeometryModel
49 template <
int dim>
class Interface;
72 using namespace ::Utilities;
92 const std::string &id_text);
198 Options(
const std::vector<std::string> &list_of_required_keys,
199 const std::string &property_name)
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),
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);
338 template <
typename T>
341 const unsigned int n_rows,
342 const unsigned int n_columns,
343 const std::string &property_name);
358 std::vector<std::string>
369 const ComponentMask &component_mask);
393 const parallel::distributed::Triangulation<dim> &triangulation,
394 const Point<dim> &point,
395 const MPI_Comm mpi_communicator);
398 namespace Coordinates
407 std::array<double,dim>
421 std::array<double,dim>
441 const ::Point<dim> &position);
452 const double semi_major_axis_a,
453 const double eccentricity);
462 const double semi_major_axis_a,
463 const double eccentricity);
483 const ::Point<2> &point);
493 const ::Point<2> &point);
504 const ::Point<2> &point);
514 std::array<Tensor<1,dim>,dim-1>
523 const double rotation_angle);
533 const Tensor<1,3> &point_two);
606 bool fexists(
const std::string &filename);
621 bool fexists(
const std::string &filename,
622 const MPI_Comm comm);
652 const MPI_Comm comm);
668 const std::string &file_content,
669 const MPI_Comm comm);
685 mkdirp(std::string pathname,
const mode_t mode = 0755);
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;
744 std::vector<double> m_a, m_b, m_c,
m_y;
758 const unsigned int q,
759 std::vector<double> &composition_values_at_q_point);
795 std::vector<unsigned int>
828 const std::vector<double> &values,
857 template <
typename T>
859 const std::vector<double> &weights,
860 const std::vector<double> &values,
861 const std::vector<T> &derivatives,
897 const SymmetricTensor<2,dim> &dviscosities_dstrain_rate,
898 const double SPD_safety_factor);
950 double operator() (
const double x,
const double y)
const;
956 bool operator== (
const operation op)
const;
1008 std::array<double,dim> &get_coordinates();
1014 const std::array<double,dim> &get_coordinates()
const;
1020 std::array<double,dim-1> get_surface_coordinates()
const;
1026 double get_depth_coordinate()
const;
1060 template <
int dim,
typename VectorType>
1063 const DoFHandler<dim> &dof_handler,
1064 const unsigned int component_index,
1065 const Quadrature<dim> &quadrature,
1067 const typename DoFHandler<dim>::active_cell_iterator &,
1068 const std::vector<Point<dim>> &,
1069 std::vector<double> &)> &
function,
1070 VectorType &vec_result);
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 =
"");
1126 const std::function<Tensor<1,dim> (
const Point<dim> &)> &function_object);
1135 template <
typename T>
1137 std::vector<std::size_t>
1147 template <
typename T>
1151 const std::vector<T> &vector,
1152 const std::vector<std::size_t> &permutation_vector);
1166 std::vector<Tensor<2,3>>
1168 const std::vector<Tensor<2,3>> &rotation_matrices,
1169 const unsigned int n_output_matrices,
1170 std::mt19937 &random_number_generator);
1199 template <
typename T>
1203 const unsigned int N,
1204 const std::string &id_text)
1206 if (values.size() == 1)
1208 return std::vector<T> (N, values[0]);
1210 else if (values.size() == N)
1218 ExcMessage(
"Length of " + id_text +
" list must be " +
1225 return std::vector<T> ();
1231 const unsigned int q,
1232 std::vector<double> &composition_values_at_q_point)
1234 Assert(q<composition_values.size(), ExcInternalError());
1235 Assert(composition_values_at_q_point.size() > 0,
1236 ExcInternalError());
1238 for (
unsigned int k=0; k < composition_values_at_q_point.size(); ++k)
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];
1246 template <
typename T>
1248 std::vector<std::size_t>
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)
1256 return vector[i] < vector[j];
1261 template <
typename T>
1265 const std::vector<T> &vector,
1266 const std::vector<std::size_t> &permutation_vector)
1268 std::vector<T> sorted_vec(vector.size());
1269 std::transform(permutation_vector.begin(), permutation_vector.end(), sorted_vec.begin(),
1287 template <
int dim,
class Iterator>
1289 SymmetricTensor<2,dim>
1290 to_symmetric_tensor(
const Iterator begin,
1293 AssertDimension(std::distance(begin, end), (SymmetricTensor<2,dim>::n_independent_components));
1296 SymmetricTensor<2,dim> output;
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;
1309 template <
int dim,
class Iterator>
1312 unroll_symmetric_tensor_into_array(
const SymmetricTensor<2,dim> &tensor,
1313 const Iterator begin,
1316 AssertDimension(std::distance(begin, end), (SymmetricTensor<2,dim>::n_independent_components));
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)];
1327 SymmetricTensor<4,3>
1328 rotate_full_stiffness_tensor(
const Tensor<2,3> &rotation_tensor,
const SymmetricTensor<4,3> &input_tensor);
1334 SymmetricTensor<2,6>
1335 rotate_voigt_stiffness_matrix(
const Tensor<2,3> &rotation_tensor,
const SymmetricTensor<2,6> &input_tensor);
1341 SymmetricTensor<2,6>
1342 to_voigt_stiffness_matrix(
const SymmetricTensor<4,3> &input_tensor);
1348 SymmetricTensor<4,3>
1349 to_full_stiffness_tensor(
const SymmetricTensor<2,6> &input_tensor);
1356 to_voigt_stiffness_vector(
const SymmetricTensor<2,6> &input_tensor);
1362 SymmetricTensor<2,6>
1363 to_voigt_stiffness_matrix(
const Tensor<1,21> &input_tensor);
1370 to_voigt_stiffness_vector(
const SymmetricTensor<4,3> &input);
1378 const Tensor<dim,dim> &levi_civita();
1382 const Tensor<3,3> &levi_civita<3>();
1396 SymmetricTensor<2,dim>
1397 consistent_deviator(
const SymmetricTensor<2,dim> &input);
1414 consistent_second_invariant_of_deviatoric_tensor(
const SymmetricTensor<2,dim> &input);
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
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
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
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
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)
bool check_values_per_key
int mkdirp(std::string pathname, const mode_t mode=0755)
Utilities::Coordinates::CoordinateSystem coordinate_system
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
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)
std::string property_name
Options(const std::vector< std::string > &list_of_required_keys, const std::string &property_name)
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
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)
bool store_values_per_key
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
std::vector< std::string > list_of_allowed_keys
bool filename_is_url(const std::string &filename)
bool allow_multiple_values_per_key
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)