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 #ifndef _aspect_geometry_model_interface_h
23 #define _aspect_geometry_model_interface_h
24 
25 #include <aspect/plugins.h>
26 #include <deal.II/base/parameter_handler.h>
27 #include <deal.II/distributed/tria.h>
28 #include <array>
29 #include <aspect/utilities.h>
31 
32 #include <set>
33 
34 
35 namespace aspect
36 {
47  namespace GeometryModel
48  {
49  using namespace dealii;
50 
58  template <int dim>
59  class Interface : public Plugins::InterfaceBase
60  {
61  public:
65  virtual
66  void create_coarse_mesh (parallel::distributed::Triangulation<dim> &coarse_grid) const = 0;
67 
81  virtual
82  double length_scale () const = 0;
83 
121  virtual
122  double depth(const Point<dim> &position) const = 0;
123 
133  virtual
134  double height_above_reference_surface(const Point<dim> &position) const = 0;
135 
141  Utilities::NaturalCoordinate<dim>
142  cartesian_to_other_coordinates(const Point<dim> &position,
143  const Utilities::Coordinates::CoordinateSystem &coordinate_system) const;
144 
148  virtual
149  aspect::Utilities::Coordinates::CoordinateSystem natural_coordinate_system() const = 0;
150 
158  virtual
159  std::array<double,dim> cartesian_to_natural_coordinates(const Point<dim> &position) const;
160 
166  virtual
167  Point<dim> natural_to_cartesian_coordinates(const std::array<double,dim> &position) const;
168 
186  virtual
187  Point<dim> representative_point(const double depth) const = 0;
188 
197  virtual
198  double maximal_depth() const = 0;
199 
200 
206  virtual
207  std::set<types::boundary_id>
208  get_used_boundary_indicators () const = 0;
209 
245  virtual
246  std::map<std::string,types::boundary_id>
247  get_symbolic_boundary_names_map () const;
248 
263  types::boundary_id
264  translate_symbolic_boundary_name_to_id (const std::string &name) const;
265 
280  std::vector<types::boundary_id>
281  translate_symbolic_boundary_names_to_ids (const std::vector<std::string> &names) const;
282 
292  std::string
293  translate_id_to_symbol_name (const types::boundary_id boundary_id) const;
294 
301  virtual
302  std::set<std::pair<std::pair<types::boundary_id, types::boundary_id>, unsigned int>>
303  get_periodic_boundary_pairs () const;
304 
326  virtual
327  void
328  adjust_positions_for_periodicity (Point<dim> &position,
329  const ArrayView<Point<dim>> &connected_positions = {},
330  const ArrayView<Tensor<1,dim>> &connected_velocities = {}) const;
331 
338  virtual
339  bool
340  has_curved_elements() const;
341 
346  virtual
347  bool
348  point_is_in_domain(const Point<dim> &p) const = 0;
349 
356  virtual
357  void
358  make_periodicity_constraints(const DoFHandler<dim> &dof_handler,
359  AffineConstraints<double> &constraints) const;
360  };
361 
362 
363 
379  template <int dim>
380  void
381  register_geometry_model (const std::string &name,
382  const std::string &description,
383  void (*declare_parameters_function) (ParameterHandler &),
384  std::unique_ptr<Interface<dim>> (*factory_function) ());
385 
396  template <int dim>
397  std::unique_ptr<Interface<dim>>
398  create_geometry_model (ParameterHandler &prm);
399 
405  template <int dim>
406  void
407  declare_parameters (ParameterHandler &prm);
408 
418  template <int dim>
419  void
420  write_plugin_graph (std::ostream &output_stream);
421 
422 
423 
431 #define ASPECT_REGISTER_GEOMETRY_MODEL(classname,name,description) \
432  template class classname<2>; \
433  template class classname<3>; \
434  namespace ASPECT_REGISTER_GEOMETRY_MODEL_ ## classname \
435  { \
436  aspect::internal::Plugins::RegisterHelper<aspect::GeometryModel::Interface<2>,classname<2>> \
437  dummy_ ## classname ## _2d (&aspect::GeometryModel::register_geometry_model<2>, \
438  name, description); \
439  aspect::internal::Plugins::RegisterHelper<aspect::GeometryModel::Interface<3>,classname<3>> \
440  dummy_ ## classname ## _3d (&aspect::GeometryModel::register_geometry_model<3>, \
441  name, description); \
442  }
443  }
444 }
445 
446 
447 #endif
std::unique_ptr< Interface< dim > > create_geometry_model(ParameterHandler &prm)
void declare_parameters(ParameterHandler &prm)
Definition: compat.h:59
void write_plugin_graph(std::ostream &output_stream)
void register_geometry_model(const std::string &name, const std::string &description, void(*declare_parameters_function)(ParameterHandler &), std::unique_ptr< Interface< dim >>(*factory_function)())
Definition: compat.h:42