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 
111  virtual Tensor<1, dim>
112  normal_vector(
113  const typename Triangulation<dim, dim>::face_iterator &face,
114  const Point<dim> &p) const override;
115 
119  virtual void
120  get_normals_at_vertices(
121  const typename Triangulation<dim, dim>::face_iterator &face,
122  typename Manifold<dim, dim>::FaceVertexNormals &face_vertex_normals)
123  const override;
124 
125 
140  virtual void
141  get_new_points(const ArrayView<const Point<dim>> &surrounding_points,
142  const Table<2, double> &weights,
143  ArrayView<Point<dim>> new_points) const override;
144 
149  virtual Point<dim>
150  get_new_point(const ArrayView<const Point<dim>> &vertices,
151  const ArrayView<const double> &weights) const override;
152 
153  private:
158 
162  const double R0, R1;
163 
168  double topography_for_point (const Point<dim> &x_y_z) const;
169  };
170 
171  }
172 
183  template <int dim>
184  class SphericalShell : public Interface<dim>, public SimulatorAccess<dim>
185  {
186  public:
190  SphericalShell() = default;
191 
200  void initialize () override;
201 
205  void create_coarse_mesh (parallel::distributed::Triangulation<dim> &coarse_grid) const override;
206 
214  double length_scale () const override;
215 
229  double depth(const Point<dim> &position) const override;
230 
235  double height_above_reference_surface(const Point<dim> &position) const override;
236 
240  Point<dim> representative_point(const double depth) const override;
241 
245  double maximal_depth() const override;
246 
268  std::set<types::boundary_id>
269  get_used_boundary_indicators () const override;
270 
276  std::map<std::string,types::boundary_id>
277  get_symbolic_boundary_names_map () const override;
278 
283  std::set<std::pair<std::pair<types::boundary_id, types::boundary_id>, unsigned int>>
284  get_periodic_boundary_pairs () const override;
285 
292  void
293  adjust_positions_for_periodicity (Point<dim> &position,
294  const ArrayView<Point<dim>> &connected_positions = {},
295  const ArrayView<Tensor<1, dim>> &connected_velocities = {}) const override;
296 
302  bool
303  has_curved_elements() const override;
304 
310  bool
311  point_is_in_domain(const Point<dim> &point) const override;
312 
317  aspect::Utilities::Coordinates::CoordinateSystem natural_coordinate_system() const override;
318 
324  std::array<double,dim> cartesian_to_natural_coordinates(const Point<dim> &position) const override;
325 
331  Point<dim> natural_to_cartesian_coordinates(const std::array<double,dim> &position) const override;
332 
333 
340  static
341  void
342  declare_parameters (ParameterHandler &prm);
343 
350  void
351  parse_parameters (ParameterHandler &prm) override;
352 
356  double
357  inner_radius () const;
358 
362  double
363  outer_radius () const;
364 
368  double
369  opening_angle () const;
370 
375  void
376  make_periodicity_constraints(const DoFHandler<dim> &dof_handler,
377  AffineConstraints<double> &constraints) const override;
378 
379  private:
385  {
388  slices
389  } custom_mesh;
390 
395 
399  unsigned int n_slices;
400 
404  std::vector<double> R_values_list;
405 
409  double R0, R1;
410 
414  double phi;
415 
420 
425  void set_manifold_ids (parallel::distributed::Triangulation<dim> &triangulation) const;
426 
430  bool periodic;
431 
443  std::unique_ptr<const internal::SphericalManifoldWithTopography<dim>> manifold;
444 
448  static const types::manifold_id my_manifold_id = 99;
449  };
450  }
451 }
452 
453 
454 #endif
std::unique_ptr< const internal::SphericalManifoldWithTopography< dim > > manifold
const InitialTopographyModel::Interface< dim > * topo
void declare_parameters(ParameterHandler &prm)
Definition: compat.h:59
Definition: compat.h:42