22 #ifndef _aspect_gravity_model_interface_h 23 #define _aspect_gravity_model_interface_h 26 #include <deal.II/base/point.h> 27 #include <deal.II/base/parameter_handler.h> 37 namespace GravityModel
45 class Interface :
public Plugins::InterfaceBase
51 virtual Tensor<1,dim>
gravity_vector (
const Point<dim> &position)
const = 0;
75 const std::string &description,
76 void (*declare_parameters_function) (ParameterHandler &),
77 std::unique_ptr<Interface<dim>> (*factory_function) ());
90 std::unique_ptr<Interface<dim>>
125 #define ASPECT_REGISTER_GRAVITY_MODEL(classname,name,description) \ 126 template class classname<2>; \ 127 template class classname<3>; \ 128 namespace ASPECT_REGISTER_GRAVITY_MODEL_ ## classname \ 130 aspect::internal::Plugins::RegisterHelper<aspect::GravityModel::Interface<2>,classname<2>> \ 131 dummy_ ## classname ## _2d (&aspect::GravityModel::register_gravity_model<2>, \ 132 name, description); \ 133 aspect::internal::Plugins::RegisterHelper<aspect::GravityModel::Interface<3>,classname<3>> \ 134 dummy_ ## classname ## _3d (&aspect::GravityModel::register_gravity_model<3>, \ 135 name, description); \ virtual Tensor< 1, dim > gravity_vector(const Point< dim > &position) const =0
void register_gravity_model(const std::string &name, const std::string &description, void(*declare_parameters_function)(ParameterHandler &), std::unique_ptr< Interface< dim >>(*factory_function)())
std::unique_ptr< Interface< dim > > create_gravity_model(ParameterHandler &prm)
void write_plugin_graph(std::ostream &output_stream)
static void declare_parameters(ParameterHandler &prm)