ASPECT
interface.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2011 - 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 
23 #ifndef _aspect_prescribed_stokes_solution_interface_h
24 #define _aspect_prescribed_stokes_solution_interface_h
25 
26 #include <aspect/plugins.h>
27 
28 #include <deal.II/base/point.h>
29 #include <deal.II/base/parameter_handler.h>
30 
31 #include <deal.II/base/function.h>
32 
33 namespace aspect
34 {
41  namespace PrescribedStokesSolution
42  {
53  template <int dim>
55  {
56  public:
63  virtual
64  void stokes_solution (const Point<dim> &p, Vector<double> &value) const = 0;
65  };
66 
67 
68 
69 
86  template <int dim>
87  void
88  register_prescribed_stokes_solution_model (const std::string &name,
89  const std::string &description,
90  void (*declare_parameters_function) (ParameterHandler &),
91  std::unique_ptr<Interface<dim>> (*factory_function) ());
92 
103  template <int dim>
104  std::unique_ptr<Interface<dim>>
105  create_prescribed_stokes_solution (ParameterHandler &prm);
106 
107 
114  template <int dim>
115  void
116  declare_parameters (ParameterHandler &prm);
117 
118 
128  template <int dim>
129  void
130  write_plugin_graph (std::ostream &output_stream);
131 
132 
133 
141 #define ASPECT_REGISTER_PRESCRIBED_STOKES_SOLUTION(classname,name,description) \
142  template class classname<2>; \
143  template class classname<3>; \
144  namespace ASPECT_REGISTER_PRESCRIBED_STOKES_SOLUTION_ ## classname \
145  { \
146  aspect::internal::Plugins::RegisterHelper<aspect::PrescribedStokesSolution::Interface<2>,classname<2>> \
147  dummy_ ## classname ## _2d (&aspect::PrescribedStokesSolution::register_prescribed_stokes_solution_model<2>, \
148  name, description); \
149  aspect::internal::Plugins::RegisterHelper<aspect::PrescribedStokesSolution::Interface<3>,classname<3>> \
150  dummy_ ## classname ## _3d (&aspect::PrescribedStokesSolution::register_prescribed_stokes_solution_model<3>, \
151  name, description); \
152  }
153  }
154 }
155 
156 
157 #endif
void register_prescribed_stokes_solution_model(const std::string &name, const std::string &description, void(*declare_parameters_function)(ParameterHandler &), std::unique_ptr< Interface< dim >>(*factory_function)())
Definition: compat.h:59
virtual void stokes_solution(const Point< dim > &p, Vector< double > &value) const =0
std::unique_ptr< Interface< dim > > create_prescribed_stokes_solution(ParameterHandler &prm)
void write_plugin_graph(std::ostream &output_stream)
static void declare_parameters(ParameterHandler &prm)