22 #ifndef _aspect_mesh_deformation_interface_h 23 #define _aspect_mesh_deformation_interface_h 29 #include <deal.II/fe/fe_system.h> 30 #include <deal.II/dofs/dof_handler.h> 31 #include <deal.II/lac/affine_constraints.h> 32 #include <deal.II/base/index_set.h> 33 #include <deal.II/base/mg_level_object.h> 34 #include <deal.II/lac/la_parallel_vector.h> 35 #include <deal.II/multigrid/mg_constrained_dofs.h> 36 #include <deal.II/multigrid/mg_transfer_matrix_free.h> 37 #include <deal.II/multigrid/mg_transfer_global_coarsening.templates.h> 77 namespace MeshDeformation
95 virtual bool needs_surface_stabilization()
const;
107 compute_initial_deformation_on_boundary(
const types::boundary_id boundary_indicator,
108 const Point<dim> &position)
const;
123 compute_initial_deformation_as_constraints(
const Mapping<dim> &mapping,
124 const DoFHandler<dim> &mesh_deformation_dof_handler,
125 const types::boundary_id boundary_indicator,
126 AffineConstraints<double> &constraints)
const;
137 compute_velocity_constraints_on_boundary(
const DoFHandler<dim> &mesh_deformation_dof_handler,
138 AffineConstraints<double> &mesh_velocity_constraints,
139 const std::set<types::boundary_id> &boundary_ids)
const;
216 template <
class Archive>
217 void save (Archive &ar,
218 const unsigned int version)
const;
224 template <
class Archive>
225 void load (Archive &ar,
226 const unsigned int version);
228 BOOST_SERIALIZATION_SPLIT_MEMBER()
249 register_mesh_deformation
250 (
const std::string &name,
251 const std::string &description,
252 void (*declare_parameters_function) (ParameterHandler &),
259 const std::map<types::boundary_id, std::vector<std::string>> &
260 get_active_mesh_deformation_names ()
const;
266 const std::map<types::boundary_id,std::vector<std::unique_ptr<Interface<dim>>>> &
267 get_active_mesh_deformation_models ()
const;
273 const std::set<types::boundary_id> &
274 get_active_mesh_deformation_boundary_indicators ()
const;
281 const std::set<types::boundary_id> &
282 get_tangential_velocity_with_active_mesh_deformation_boundary_indicators ()
const;
289 const std::set<types::boundary_id> &
290 get_tangential_velocity_without_active_mesh_deformation_boundary_indicators ()
const;
296 const std::set<types::boundary_id> &
297 get_boundary_indicators_requiring_stabilization ()
const;
304 const std::set<types::boundary_id> &
305 get_free_surface_boundary_indicators ()
const;
310 double get_free_surface_theta ()
const;
327 get_initial_topography ()
const;
334 get_mesh_displacements ()
const;
339 const DoFHandler<dim> &
340 get_mesh_deformation_dof_handler ()
const;
351 template <
typename MeshDeformationType,
352 typename =
typename std::enable_if_t<std::is_base_of<Interface<dim>,MeshDeformationType>::value>>
354 has_matching_mesh_deformation_object ()
const;
367 template <
typename MeshDeformationType,
368 typename =
typename std::enable_if_t<std::is_base_of<Interface<dim>,MeshDeformationType>::value>>
369 const MeshDeformationType &
370 get_matching_mesh_deformation_object ()
const;
379 get_level_mapping(
const unsigned int level)
const;
399 <<
"Could not find entry <" 401 <<
"> among the names of registered mesh deformation objects.");
411 void make_initial_constraints ();
425 void make_constraints ();
431 void compute_mesh_displacements ();
437 void compute_mesh_displacements_gmg ();
454 void set_initial_topography ();
459 void interpolate_mesh_velocity ();
464 void update_multilevel_deformation ();
640 template <
class Archive>
642 const unsigned int)
const 647 std::map<std::string,std::string> saved_text;
648 for (
const auto &boundary_id_and_mesh_deformation_objects : mesh_deformation_objects)
649 for (
const auto &p : boundary_id_and_mesh_deformation_objects.second)
650 p->
save (saved_text);
657 template <
class Archive>
665 std::map<std::string,std::string> saved_text;
668 for (
const auto &boundary_id_and_mesh_deformation_objects : mesh_deformation_objects)
669 for (
const auto &p : boundary_id_and_mesh_deformation_objects.second)
670 p->
load (saved_text);
677 template <
typename MeshDeformationType,
typename>
682 for (
const auto &object_iterator : mesh_deformation_objects)
683 for (
const auto &p : object_iterator.second)
684 if (Plugins::plugin_type_matches<MeshDeformationType>(*p))
693 template <
typename MeshDeformationType,
typename>
695 const MeshDeformationType &
698 AssertThrow(has_matching_mesh_deformation_object<MeshDeformationType> (),
699 ExcMessage(
"You asked MeshDeformation::MeshDeformationHandler::get_matching_mesh_deformation_object() for a " 700 "mesh deformation object of type <" + boost::core::demangle(
typeid(MeshDeformationType).name()) +
"> " 701 "that could not be found in the current model. Activate this " 702 "mesh deformation in the input file."));
704 for (
const auto &object_iterator : mesh_deformation_objects)
705 for (
const auto &p : object_iterator.second)
706 if (Plugins::plugin_type_matches<MeshDeformationType>(*p))
707 return Plugins::get_plugin_as_type<MeshDeformationType>(*p);
710 typename std::vector<std::unique_ptr<Interface<dim>>>::const_iterator mesh_def;
711 return Plugins::get_plugin_as_type<MeshDeformationType>(*(*mesh_def));
735 #define ASPECT_REGISTER_MESH_DEFORMATION_MODEL(classname,name,description) \ 736 template class classname<2>; \ 737 template class classname<3>; \ 738 namespace ASPECT_REGISTER_MESH_DEFORMATION_MODEL_ ## classname \ 740 aspect::internal::Plugins::RegisterHelper<aspect::MeshDeformation::Interface<2>,classname<2>> \ 741 dummy_ ## classname ## _2d (&aspect::MeshDeformation::MeshDeformationHandler<2>::register_mesh_deformation, \ 742 name, description); \ 743 aspect::internal::Plugins::RegisterHelper<aspect::MeshDeformation::Interface<3>,classname<3>> \ 744 dummy_ ## classname ## _3d (&aspect::MeshDeformation::MeshDeformationHandler<3>::register_mesh_deformation, \ 745 name, description); \
virtual void parse_parameters(ParameterHandler &prm)
const Simulator< dim > * simulator
::TrilinosWrappers::MPI::BlockVector BlockVector
ApplyStabilization(const double stabilization_theta)
void write_plugin_graph(std::ostream &output_stream)
::TrilinosWrappers::MPI::Vector Vector
virtual void load(const std::map< std::string, std::string > &status_strings)
const double free_surface_theta
std::string get_valid_model_names_pattern()
::MGTransferMatrixFree< dim, NumberType > MGTransferMF
void execute(internal::Assembly::Scratch::ScratchBase< dim > &scratch, internal::Assembly::CopyData::CopyDataBase< dim > &data) const override
virtual void save(std::map< std::string, std::string > &status_strings) const
virtual void initialize()
static void declare_parameters(ParameterHandler &prm)
DeclException1(ProbabilityFunctionNegative, Point< dim >,<< "Your probability density function in the particle generator " "returned a negative probability density for the following position: "<< arg1<< ". Please check your function expression.")