ASPECT
spherical_shell.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_spherical_shell_h
23 #define _aspect_geometry_model_spherical_shell_h
24 
27 
28 namespace aspect
29 {
30  namespace GeometryModel
31  {
32  using namespace dealii;
33 
34  namespace internal
35  {
40  template <int dim>
42  {
43  public:
48  const double inner_radius,
49  const double outer_radius);
50 
55 
59  virtual std::unique_ptr<Manifold<dim, dim>>
60  clone() const override;
61 
66  Point<dim>
67  push_forward_from_sphere (const Point<dim> &p) const;
68 
73  Point<dim>
74  pull_back_to_sphere (const Point<dim> &p) const;
75 
83  virtual Point<dim>
84  get_intermediate_point(const Point<dim> &p1,
85  const Point<dim> &p2,
86  const double w) const override;
87 
92  virtual Tensor<1, dim>
93  get_tangent_vector(const Point<dim> &x1,
94  const Point<dim> &x2) const override;
95 
99  virtual Tensor<1, dim>
100  normal_vector(
101  const typename Triangulation<dim, dim>::face_iterator &face,
102  const Point<dim> &p) const override;
103 
107  virtual void
108  get_normals_at_vertices(
109  const typename Triangulation<dim, dim>::face_iterator &face,
110  typename Manifold<dim, dim>::FaceVertexNormals &face_vertex_normals)
111  const override;
112 
113 
128  virtual void
129  get_new_points(const ArrayView<const Point<dim>> &surrounding_points,
130  const Table<2, double> &weights,
131  ArrayView<Point<dim>> new_points) const override;
132 
137  virtual Point<dim>
138  get_new_point(const ArrayView<const Point<dim>> &vertices,
139  const ArrayView<const double> &weights) const override;
140 
141  private:
146 
150  const double R0, R1;
151 
156  double topography_for_point (const Point<dim> &x_y_z) const;
157  };
158 
159  }
160 
171  template <int dim>
172  class SphericalShell : public Interface<dim>, public SimulatorAccess<dim>
173  {
174  public:
178  SphericalShell() = default;
179 
188  void initialize () override;
189 
193  void create_coarse_mesh (parallel::distributed::Triangulation<dim> &coarse_grid) const override;
194 
202  double length_scale () const override;
203 
217  double depth(const Point<dim> &position) const override;
218 
223  double height_above_reference_surface(const Point<dim> &position) const override;
224 
228  Point<dim> representative_point(const double depth) const override;
229 
233  double maximal_depth() const override;
234 
256  std::set<types::boundary_id>
257  get_used_boundary_indicators () const override;
258 
264  std::map<std::string,types::boundary_id>
265  get_symbolic_boundary_names_map () const override;
266 
271  std::set<std::pair<std::pair<types::boundary_id, types::boundary_id>, unsigned int>>
272  get_periodic_boundary_pairs () const override;
273 
280  void
281  adjust_positions_for_periodicity (Point<dim> &position,
282  const ArrayView<Point<dim>> &connected_positions = {},
283  const ArrayView<Tensor<1, dim>> &connected_velocities = {}) const override;
284 
290  bool
291  has_curved_elements() const override;
292 
298  bool
299  point_is_in_domain(const Point<dim> &point) const override;
300 
305  aspect::Utilities::Coordinates::CoordinateSystem natural_coordinate_system() const override;
306 
312  std::array<double,dim> cartesian_to_natural_coordinates(const Point<dim> &position) const override;
313 
319  Point<dim> natural_to_cartesian_coordinates(const std::array<double,dim> &position) const override;
320 
321 
328  static
329  void
330  declare_parameters (ParameterHandler &prm);
331 
338  void
339  parse_parameters (ParameterHandler &prm) override;
340 
344  double
345  inner_radius () const;
346 
350  double
351  outer_radius () const;
352 
356  double
357  opening_angle () const;
358 
363  void
364  make_periodicity_constraints(const DoFHandler<dim> &dof_handler,
365  AffineConstraints<double> &constraints) const override;
366 
367  private:
373  {
376  slices
377  } custom_mesh;
378 
383 
387  unsigned int n_slices;
388 
392  std::vector<double> R_values_list;
393 
397  double R0, R1;
398 
402  double phi;
403 
408 
413  void set_manifold_ids (parallel::distributed::Triangulation<dim> &triangulation) const;
414 
418  bool periodic;
419 
431  std::unique_ptr<const internal::SphericalManifoldWithTopography<dim>> manifold;
432 
436  static const types::manifold_id my_manifold_id = 99;
437  };
438  }
439 }
440 
441 
442 #endif
std::unique_ptr< const internal::SphericalManifoldWithTopography< dim > > manifold
const InitialTopographyModel::Interface< dim > * topo
void declare_parameters(ParameterHandler &prm)
Definition: compat.h:42