ASPECT
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
spherical_shell.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_spherical_shell_h
23 #define _aspect_geometry_model_spherical_shell_h
24 
27 
28 namespace aspect
29 {
30  namespace GeometryModel
31  {
32  namespace internal
33  {
38  template <int dim>
40  {
41  public:
45  SphericalManifoldWithTopography(const std::shared_ptr<const InitialTopographyModel::Interface<dim>> topography,
46  const double inner_radius,
47  const double outer_radius);
48 
53 
57  virtual std::unique_ptr<Manifold<dim, dim>>
58  clone() const override;
59 
64  Point<dim>
65  push_forward_from_sphere (const Point<dim> &p) const;
66 
71  Point<dim>
72  pull_back_to_sphere (const Point<dim> &p) const;
73 
81  virtual Point<dim>
82  get_intermediate_point(const Point<dim> &p1,
83  const Point<dim> &p2,
84  const double w) const override;
85 
90  virtual Tensor<1, dim>
91  get_tangent_vector(const Point<dim> &x1,
92  const Point<dim> &x2) const override;
93 
109  virtual Tensor<1, dim>
111  const typename Triangulation<dim, dim>::face_iterator &face,
112  const Point<dim> &p) const override;
113 
117  virtual void
119  const typename Triangulation<dim, dim>::face_iterator &face,
120  typename Manifold<dim, dim>::FaceVertexNormals &face_vertex_normals)
121  const override;
122 
123 
138  virtual void
139  get_new_points(const ArrayView<const Point<dim>> &surrounding_points,
140  const Table<2, double> &weights,
141  ArrayView<Point<dim>> new_points) const override;
142 
147  virtual Point<dim>
148  get_new_point(const ArrayView<const Point<dim>> &vertices,
149  const ArrayView<const double> &weights) const override;
154  double topography_for_point (const Point<dim> &x_y_z) const;
155 
156  private:
160  const std::shared_ptr<const InitialTopographyModel::Interface<dim>> topo;
161 
165  const double R0, R1;
166  };
167 
168  }
169 
180  template <int dim>
181  class SphericalShell : public Interface<dim>, public SimulatorAccess<dim>
182  {
183  public:
187  SphericalShell() = default;
188 
197  void initialize () override;
198 
202  void create_coarse_mesh (parallel::distributed::Triangulation<dim> &coarse_grid) const override;
203 
211  double length_scale () const override;
212 
226  double depth(const Point<dim> &position) const override;
227 
232  double height_above_reference_surface(const Point<dim> &position) const override;
233 
237  Point<dim> representative_point(const double depth) const override;
238 
242  double maximal_depth() const override;
243 
265  std::set<types::boundary_id>
266  get_used_boundary_indicators () const override;
267 
273  std::map<std::string,types::boundary_id>
274  get_symbolic_boundary_names_map () const override;
275 
280  std::set<std::pair<std::pair<types::boundary_id, types::boundary_id>, unsigned int>>
281  get_periodic_boundary_pairs () const override;
282 
289  void
290  adjust_positions_for_periodicity (Point<dim> &position,
291  const ArrayView<Point<dim>> &connected_positions = {},
292  const ArrayView<Tensor<1, dim>> &connected_velocities = {}) const override;
293 
299  bool
300  has_curved_elements() const override;
301 
307  bool
308  point_is_in_domain(const Point<dim> &point) const override;
309 
314  aspect::Utilities::Coordinates::CoordinateSystem natural_coordinate_system() const override;
315 
321  std::array<double,dim> cartesian_to_natural_coordinates(const Point<dim> &position) const override;
322 
328  Point<dim> natural_to_cartesian_coordinates(const std::array<double,dim> &position) const override;
329 
330 
334  static
335  void
336  declare_parameters (ParameterHandler &prm);
337 
341  void
342  parse_parameters (ParameterHandler &prm) override;
343 
347  double
348  inner_radius () const;
349 
353  double
354  outer_radius () const;
355 
359  double
360  opening_angle () const;
361 
366  void
367  make_periodicity_constraints(const DoFHandler<dim> &dof_handler,
368  AffineConstraints<double> &constraints) const override;
369 
370  private:
376  {
379  slices
380  } custom_mesh;
381 
386 
390  unsigned int n_slices;
391 
395  std::vector<double> R_values_list;
396 
400  double R0, R1;
401 
405  double phi;
406 
411 
416  void set_manifold_ids (parallel::distributed::Triangulation<dim> &triangulation) const;
417 
421  bool periodic;
422 
434  std::unique_ptr<const internal::SphericalManifoldWithTopography<dim>> manifold;
435 
439  static constexpr types::manifold_id my_manifold_id = 99;
440  };
441  }
442 }
443 
444 
445 #endif
virtual Point< dim > get_new_point(const ArrayView< const Point< dim >> &vertices, const ArrayView< const double > &weights) const override
std::unique_ptr< const internal::SphericalManifoldWithTopography< dim > > manifold
virtual Point< dim > get_intermediate_point(const Point< dim > &p1, const Point< dim > &p2, const double w) const override
const std::shared_ptr< const InitialTopographyModel::Interface< dim > > topo
virtual std::unique_ptr< Manifold< dim, dim > > clone() const override
SphericalManifoldWithTopography(const std::shared_ptr< const InitialTopographyModel::Interface< dim >> topography, const double inner_radius, const double outer_radius)
void declare_parameters(ParameterHandler &prm)
Point< dim > push_forward_from_sphere(const Point< dim > &p) const
virtual void get_new_points(const ArrayView< const Point< dim >> &surrounding_points, const Table< 2, double > &weights, ArrayView< Point< dim >> new_points) const override
virtual void get_normals_at_vertices(const typename Triangulation< dim, dim >::face_iterator &face, typename Manifold< dim, dim >::FaceVertexNormals &face_vertex_normals) const override
Point< dim > pull_back_to_sphere(const Point< dim > &p) const
virtual Tensor< 1, dim > normal_vector(const typename Triangulation< dim, dim >::face_iterator &face, const Point< dim > &p) const override
double topography_for_point(const Point< dim > &x_y_z) const
virtual Tensor< 1, dim > get_tangent_vector(const Point< dim > &x1, const Point< dim > &x2) const override