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_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
60  {
61  public:
66  virtual ~Interface() = default;
67 
73  virtual void initialize ();
74 
78  virtual
79  void create_coarse_mesh (parallel::distributed::Triangulation<dim> &coarse_grid) const = 0;
80 
94  virtual
95  double length_scale () const = 0;
96 
134  virtual
135  double depth(const Point<dim> &position) const = 0;
136 
146  virtual
147  double height_above_reference_surface(const Point<dim> &position) const = 0;
148 
154  Utilities::NaturalCoordinate<dim>
155  cartesian_to_other_coordinates(const Point<dim> &position,
156  const Utilities::Coordinates::CoordinateSystem &coordinate_system) const;
157 
161  virtual
162  aspect::Utilities::Coordinates::CoordinateSystem natural_coordinate_system() const = 0;
163 
171  virtual
172  std::array<double,dim> cartesian_to_natural_coordinates(const Point<dim> &position) const;
173 
179  virtual
180  Point<dim> natural_to_cartesian_coordinates(const std::array<double,dim> &position) const;
181 
199  virtual
200  Point<dim> representative_point(const double depth) const = 0;
201 
210  virtual
211  double maximal_depth() const = 0;
212 
213 
219  virtual
220  std::set<types::boundary_id>
221  get_used_boundary_indicators () const = 0;
222 
258  virtual
259  std::map<std::string,types::boundary_id>
260  get_symbolic_boundary_names_map () const;
261 
276  types::boundary_id
277  translate_symbolic_boundary_name_to_id (const std::string &name) const;
278 
293  std::vector<types::boundary_id>
294  translate_symbolic_boundary_names_to_ids (const std::vector<std::string> &names) const;
295 
305  std::string
306  translate_id_to_symbol_name (const types::boundary_id boundary_id) const;
307 
314  virtual
315  std::set<std::pair<std::pair<types::boundary_id, types::boundary_id>, unsigned int>>
316  get_periodic_boundary_pairs () const;
317 
339  virtual
340  void
341  adjust_positions_for_periodicity (Point<dim> &position,
342  const ArrayView<Point<dim>> &connected_positions = {},
343  const ArrayView<Tensor<1,dim>> &connected_velocities = {}) const;
344 
351  virtual
352  bool
353  has_curved_elements() const;
354 
359  virtual
360  bool
361  point_is_in_domain(const Point<dim> &p) const = 0;
362 
369  static
370  void
371  declare_parameters (ParameterHandler &prm);
372 
379  virtual
380  void
381  parse_parameters (ParameterHandler &prm);
382 
389  virtual
390  void
391  make_periodicity_constraints(const DoFHandler<dim> &dof_handler,
392  AffineConstraints<double> &constraints) const;
393  };
394 
395 
396 
412  template <int dim>
413  void
414  register_geometry_model (const std::string &name,
415  const std::string &description,
416  void (*declare_parameters_function) (ParameterHandler &),
417  std::unique_ptr<Interface<dim>> (*factory_function) ());
418 
429  template <int dim>
430  std::unique_ptr<Interface<dim>>
431  create_geometry_model (ParameterHandler &prm);
432 
438  template <int dim>
439  void
440  declare_parameters (ParameterHandler &prm);
441 
451  template <int dim>
452  void
453  write_plugin_graph (std::ostream &output_stream);
454 
455 
456 
464 #define ASPECT_REGISTER_GEOMETRY_MODEL(classname,name,description) \
465  template class classname<2>; \
466  template class classname<3>; \
467  namespace ASPECT_REGISTER_GEOMETRY_MODEL_ ## classname \
468  { \
469  aspect::internal::Plugins::RegisterHelper<aspect::GeometryModel::Interface<2>,classname<2>> \
470  dummy_ ## classname ## _2d (&aspect::GeometryModel::register_geometry_model<2>, \
471  name, description); \
472  aspect::internal::Plugins::RegisterHelper<aspect::GeometryModel::Interface<3>,classname<3>> \
473  dummy_ ## classname ## _3d (&aspect::GeometryModel::register_geometry_model<3>, \
474  name, description); \
475  }
476  }
477 }
478 
479 
480 #endif
std::unique_ptr< Interface< dim > > create_geometry_model(ParameterHandler &prm)
void declare_parameters(ParameterHandler &prm)
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