22 #ifndef _aspect_boundary_velocity_interface_h 23 #define _aspect_boundary_velocity_interface_h 30 #include <deal.II/base/point.h> 31 #include <deal.II/base/parameter_handler.h> 33 #include <boost/core/demangle.hpp> 45 namespace BoundaryVelocity
55 class Interface :
public Plugins::InterfaceBase
72 boundary_velocity (
const types::boundary_id boundary_indicator,
73 const Point<dim> &position)
const = 0;
82 class Manager :
public SimulatorAccess<dim>
110 boundary_velocity (
const types::boundary_id boundary_indicator,
111 const Point<dim> &position)
const;
132 register_boundary_velocity (
const std::string &name,
133 const std::string &description,
134 void (*declare_parameters_function) (ParameterHandler &),
135 std::unique_ptr<Interface<dim>> (*factory_function) ());
150 const std::map<types::boundary_id, std::pair<std::string,std::vector<std::string>>> &
151 get_active_boundary_velocity_names ()
const;
162 const std::map<types::boundary_id,std::vector<std::unique_ptr<BoundaryVelocity::Interface<dim>>>> &
163 get_active_boundary_velocity_conditions ()
const;
169 const std::set<types::boundary_id> &
170 get_zero_boundary_velocity_indicators ()
const;
176 const std::set<types::boundary_id> &
177 get_tangential_boundary_velocity_indicators ()
const;
193 parse_parameters (ParameterHandler &prm);
204 template <
typename BoundaryVelocityType,
205 typename =
typename std::enable_if_t<std::is_base_of<Interface<dim>,BoundaryVelocityType>::value>>
207 has_matching_boundary_velocity_model ()
const;
220 template <
typename BoundaryVelocityType,
221 typename =
typename std::enable_if_t<std::is_base_of<Interface<dim>,BoundaryVelocityType>::value>>
222 const BoundaryVelocityType &
223 get_matching_boundary_velocity_model ()
const;
244 <<
"Could not find entry <" 246 <<
"> among the names of registered boundary velocity objects.");
280 template <
typename BoundaryVelocityType,
typename>
285 for (
const auto &boundary : boundary_velocity_objects)
286 for (
const auto &p : boundary.second)
287 if (Plugins::plugin_type_matches<BoundaryVelocityType>(*p))
294 template <
typename BoundaryVelocityType,
typename>
296 const BoundaryVelocityType &
299 AssertThrow(has_matching_boundary_velocity_model<BoundaryVelocityType> (),
300 ExcMessage(
"You asked BoundaryVelocity::Manager::get_boundary_velocity_model() for a " 301 "boundary velocity model of type <" + boost::core::demangle(
typeid(BoundaryVelocityType).name()) +
"> " 302 "that could not be found in the current model. Activate this " 303 "boundary velocity model in the input file."));
305 for (
const auto &boundary : boundary_velocity_objects)
306 for (
const auto &p : boundary)
307 if (Plugins::plugin_type_matches<BoundaryVelocityType>(*p))
308 return Plugins::get_plugin_as_type<BoundaryVelocityType>(*p);
311 return Plugins::get_plugin_as_type<BoundaryVelocityType>(**(boundary_velocity_objects.begin()));
333 #define ASPECT_REGISTER_BOUNDARY_VELOCITY_MODEL(classname, name, description) \ 334 template class classname<2>; \ 335 template class classname<3>; \ 336 namespace ASPECT_REGISTER_BOUNDARY_VELOCITY_MODEL_ ## classname \ 338 aspect::internal::Plugins::RegisterHelper<aspect::BoundaryVelocity::Interface<2>,classname<2>> \ 339 dummy_ ## classname ## _2d (&aspect::BoundaryVelocity::Manager<2>::register_boundary_velocity, \ 340 name, description); \ 341 aspect::internal::Plugins::RegisterHelper<aspect::BoundaryVelocity::Interface<3>,classname<3>> \ 342 dummy_ ## classname ## _3d (&aspect::BoundaryVelocity::Manager<3>::register_boundary_velocity, \ 343 name, description); \
void write_plugin_graph(std::ostream &output_stream)
std::set< types::boundary_id > tangential_velocity_boundary_indicators
std::map< types::boundary_id, std::pair< std::string, std::vector< std::string > > > boundary_velocity_indicators
std::string get_valid_model_names_pattern()
std::map< types::boundary_id, std::vector< std::unique_ptr< BoundaryVelocity::Interface< dim > > > > boundary_velocity_objects
bool has_matching_boundary_velocity_model() const
const BoundaryVelocityType & get_matching_boundary_velocity_model() const
void declare_parameters(ParameterHandler &prm)
std::set< types::boundary_id > zero_velocity_boundary_indicators
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.")