ASPECT
interface.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2011 - 2022 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  {
43  using namespace dealii;
44 
55  template <int dim>
56  class Interface
57  {
58  public:
63  virtual ~Interface() = default;
64 
70  virtual
71  void
72  initialize ();
73 
84  virtual
85  void
86  update ();
87 
94  virtual
95  void stokes_solution (const Point<dim> &p, Vector<double> &value) const = 0;
96 
103  static
104  void
105  declare_parameters (ParameterHandler &prm);
106 
113  virtual
114  void
115  parse_parameters (ParameterHandler &prm);
116 
117  };
118 
119 
120 
121 
138  template <int dim>
139  void
140  register_prescribed_stokes_solution_model (const std::string &name,
141  const std::string &description,
142  void (*declare_parameters_function) (ParameterHandler &),
143  std::unique_ptr<Interface<dim>> (*factory_function) ());
144 
155  template <int dim>
156  std::unique_ptr<Interface<dim>>
157  create_prescribed_stokes_solution (ParameterHandler &prm);
158 
159 
166  template <int dim>
167  void
168  declare_parameters (ParameterHandler &prm);
169 
170 
180  template <int dim>
181  void
182  write_plugin_graph (std::ostream &output_stream);
183 
184 
185 
193 #define ASPECT_REGISTER_PRESCRIBED_STOKES_SOLUTION(classname,name,description) \
194  template class classname<2>; \
195  template class classname<3>; \
196  namespace ASPECT_REGISTER_PRESCRIBED_STOKES_SOLUTION_ ## classname \
197  { \
198  aspect::internal::Plugins::RegisterHelper<aspect::PrescribedStokesSolution::Interface<2>,classname<2>> \
199  dummy_ ## classname ## _2d (&aspect::PrescribedStokesSolution::register_prescribed_stokes_solution_model<2>, \
200  name, description); \
201  aspect::internal::Plugins::RegisterHelper<aspect::PrescribedStokesSolution::Interface<3>,classname<3>> \
202  dummy_ ## classname ## _3d (&aspect::PrescribedStokesSolution::register_prescribed_stokes_solution_model<3>, \
203  name, description); \
204  }
205  }
206 }
207 
208 
209 #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)())
void declare_parameters(ParameterHandler &prm)
std::unique_ptr< Interface< dim > > create_prescribed_stokes_solution(ParameterHandler &prm)
void write_plugin_graph(std::ostream &output_stream)
Definition: compat.h:42