ASPECT
ellipsoidal_chunk.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_ellipsoidal_chunk_h
23 #define _aspect_geometry_model_ellipsoidal_chunk_h
24 
28 #include <deal.II/grid/manifold.h>
29 
30 
31 namespace aspect
32 {
33  namespace GeometryModel
34  {
35  namespace internal
36  {
40  template <int dim>
41  class EllipsoidalChunkGeometry : public ChartManifold<dim,3,3>
42  {
43  public:
48  const double para_semi_major_axis_a,
49  const double para_eccentricity,
50  const double para_semi_minor_axis_b,
51  const double para_bottom_depth,
52  const std::vector<Point<2>> &para_corners);
53 
57  EllipsoidalChunkGeometry(const EllipsoidalChunkGeometry &other) = default;
58 
65  Point<3>
66  pull_back(const Point<3> &space_point) const override;
67 
72  virtual
73  Point<2>
74  pull_back(const Point<2> &space_point) const;
75 
82  Point<3>
83  push_forward(const Point<3> &chart_point) const override;
84 
89  Point<3> pull_back_ellipsoid (const Point<3> &x, const double semi_major_axis_a, const double eccentricity) const;
90 
95  Point<3> push_forward_ellipsoid (const Point<3> &phi_theta_d, const double semi_major_axis_a, const double eccentricity) const;
96 
100  std::unique_ptr<Manifold<dim,3>>
101  clone() const override;
102 
103  private:
108  Point<3> push_forward_topography (const Point<3> &phi_theta_d_hat) const;
109 
114  Point<3> pull_back_topography (const Point<3> &phi_theta_d) const;
115 
120 
121  const double semi_major_axis_a;
122  const double eccentricity;
123  const double semi_minor_axis_b;
124  const double bottom_depth;
125  const std::vector<Point<2>> corners;
126  };
127  }
128 
137  template <int dim>
138  class EllipsoidalChunk : public Interface<dim>, public SimulatorAccess<dim>
139  {
140  public:
144  void
145  initialize () override;
146 
147 
151  void
152  create_coarse_mesh(parallel::distributed::Triangulation<dim> &coarse_grid) const override;
153 
158  double
159  length_scale() const override;
160 
174  double
175  depth(const Point<dim> &position) const override;
176 
181  double
182  height_above_reference_surface(const Point<dim> &position) const override;
183 
187  Point<dim>
188  representative_point(const double depth) const override;
189 
195  bool
196  point_is_in_domain(const Point<dim> &point) const override;
197 
202  double
203  maximal_depth() const override;
204 
214  std::set<types::boundary_id>
215  get_used_boundary_indicators() const override;
216 
220  std::map<std::string,types::boundary_id>
221  get_symbolic_boundary_names_map() const override;
222 
227  aspect::Utilities::Coordinates::CoordinateSystem natural_coordinate_system() const override;
228 
236  std::array<double,dim> cartesian_to_natural_coordinates(const Point<dim> &position) const override;
237 
243  Point<dim> natural_to_cartesian_coordinates(const std::array<double,dim> &position) const override;
244 
248  static
249  void
250  declare_parameters(ParameterHandler &prm);
251 
256  void
257  parse_parameters(ParameterHandler &prm) override;
258 
262  double
263  get_radius(const Point<dim> &position) const;
264 
268  double
269  get_semi_minor_axis_b() const;
270 
274  double
275  get_semi_major_axis_a() const;
276 
277 
281  double
282  get_eccentricity() const;
283 
284 
292  const std::vector<Point<2>> &
293  get_corners() const;
294 
295 
300  get_manifold() const;
301 
302  private:
303  // Declare variables for reading in coordinates of the region of interest.
304  std::vector<Point<2>> corners;
306  double eccentricity;
312  double bottom_depth;
317  // Declare variables for subdividing
318  unsigned int EW_subdiv;
319  unsigned int NS_subdiv;
320  unsigned int depth_subdiv;
321 
322 
323 
335  std::unique_ptr<const internal::EllipsoidalChunkGeometry<dim>> manifold;
336 
337  void
338  set_boundary_ids(parallel::distributed::Triangulation<dim> &coarse_grid) const;
339  };
340  }
341 }
342 
343 
344 #endif
Point< 3 > push_forward(const Point< 3 > &chart_point) const override
Point< 3 > push_forward_topography(const Point< 3 > &phi_theta_d_hat) const
Point< 3 > pull_back(const Point< 3 > &space_point) const override
void declare_parameters(ParameterHandler &prm)
Point< 3 > pull_back_ellipsoid(const Point< 3 > &x, const double semi_major_axis_a, const double eccentricity) const
const InitialTopographyModel::Interface< dim > * topography
EllipsoidalChunkGeometry(const InitialTopographyModel::Interface< dim > &topography, const double para_semi_major_axis_a, const double para_eccentricity, const double para_semi_minor_axis_b, const double para_bottom_depth, const std::vector< Point< 2 >> &para_corners)
std::unique_ptr< const internal::EllipsoidalChunkGeometry< dim > > manifold
Definition: compat.h:59
std::unique_ptr< Manifold< dim, 3 > > clone() const override
Point< 3 > push_forward_ellipsoid(const Point< 3 > &phi_theta_d, const double semi_major_axis_a, const double eccentricity) const
Point< 3 > pull_back_topography(const Point< 3 > &phi_theta_d) const