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 #ifndef _aspect_geometry_model_initial_topography_model_interface_h
23 #define _aspect_geometry_model_initial_topography_model_interface_h
24 
25 #include <aspect/plugins.h>
26 #include <deal.II/base/parameter_handler.h>
27 #include <deal.II/base/point.h>
28 
29 #include <set>
30 
31 
32 namespace aspect
33 {
42  namespace InitialTopographyModel
43  {
44  using namespace dealii;
45 
52  template <int dim>
53  class Interface
54  {
55  public:
60  virtual ~Interface() = default;
61 
67  virtual void initialize ();
68 
79  virtual
80  double value (const Point<dim-1> &surface_point) const = 0;
81 
85  virtual
86  double max_topography () const = 0;
87 
94  static
95  void
96  declare_parameters (ParameterHandler &prm);
97 
104  virtual
105  void
106  parse_parameters (ParameterHandler &prm);
107  };
108 
109 
110 
126  template <int dim>
127  void
128  register_initial_topography_model (const std::string &name,
129  const std::string &description,
130  void (*declare_parameters_function) (ParameterHandler &),
131  std::unique_ptr<Interface<dim>> (*factory_function) ());
132 
143  template <int dim>
144  std::unique_ptr<Interface<dim>>
145  create_initial_topography_model (ParameterHandler &prm);
146 
147 
154  template <int dim>
155  void
156  declare_parameters (ParameterHandler &prm);
157 
158 
168  template <int dim>
169  void
170  write_plugin_graph (std::ostream &output_stream);
171 
172 
180 #define ASPECT_REGISTER_INITIAL_TOPOGRAPHY_MODEL(classname,name,description) \
181  template class classname<2>; \
182  template class classname<3>; \
183  namespace ASPECT_REGISTER_INITIAL_TOPOGRAPHY_MODEL_ ## classname \
184  { \
185  aspect::internal::Plugins::RegisterHelper<aspect::InitialTopographyModel::Interface<2>,classname<2>> \
186  dummy_ ## classname ## _2d (&aspect::InitialTopographyModel::register_initial_topography_model<2>, \
187  name, description); \
188  aspect::internal::Plugins::RegisterHelper<aspect::InitialTopographyModel::Interface<3>,classname<3>> \
189  dummy_ ## classname ## _3d (&aspect::InitialTopographyModel::register_initial_topography_model<3>, \
190  name, description); \
191  }
192  }
193 }
194 
195 
196 #endif
void register_initial_topography_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_initial_topography_model(ParameterHandler &prm)
void declare_parameters(ParameterHandler &prm)
void write_plugin_graph(std::ostream &output_stream)
Definition: compat.h:42