ASPECT
ellipsoidal_chunk.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_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  using namespace dealii;
36 
45  template <int dim>
46  class EllipsoidalChunk : public Interface<dim>, public SimulatorAccess<dim>
47  {
48  public:
52  class EllipsoidalChunkGeometry : public ChartManifold<dim,3,3>
53  {
54  public:
59 
64 
69  void
71 
75  void
76  set_manifold_parameters(const double para_semi_major_axis_a,
77  const double para_eccentricity,
78  const double para_semi_minor_axis_b,
79  const double para_bottom_depth,
80  const std::vector<Point<2>> &para_corners);
81 
88  Point<3>
89  pull_back(const Point<3> &space_point) const override;
90 
95  virtual
96  Point<2>
97  pull_back(const Point<2> &space_point) const;
98 
105  Point<3>
106  push_forward(const Point<3> &chart_point) const override;
107 
112  Point<3> pull_back_ellipsoid (const Point<3> &x, const double semi_major_axis_a, const double eccentricity) const;
113 
118  Point<3> push_forward_ellipsoid (const Point<3> &phi_theta_d, const double semi_major_axis_a, const double eccentricity) const;
119 
123  std::unique_ptr<Manifold<dim,3>>
124  clone() const override;
125 
126  private:
131  Point<3> push_forward_topography (const Point<3> &phi_theta_d_hat) const;
132 
137  Point<3> pull_back_topography (const Point<3> &phi_theta_d) const;
138 
139 
141  double eccentricity;
143  double bottom_depth;
144  std::vector<Point<2>> corners;
146  };
147 
151  void
152  initialize () override;
153 
154 
158  void
159  create_coarse_mesh(parallel::distributed::Triangulation<dim> &coarse_grid) const override;
160 
165  double
166  length_scale() const override;
167 
181  double
182  depth(const Point<dim> &position) const override;
183 
188  double
189  height_above_reference_surface(const Point<dim> &position) const override;
190 
194  Point<dim>
195  representative_point(const double depth) const override;
196 
202  bool
203  point_is_in_domain(const Point<dim> &point) const override;
204 
209  double
210  maximal_depth() const override;
211 
221  std::set<types::boundary_id>
222  get_used_boundary_indicators() const override;
223 
227  std::map<std::string,types::boundary_id>
228  get_symbolic_boundary_names_map() const override;
229 
234  aspect::Utilities::Coordinates::CoordinateSystem natural_coordinate_system() const override;
235 
243  std::array<double,dim> cartesian_to_natural_coordinates(const Point<dim> &position) const override;
244 
250  Point<dim> natural_to_cartesian_coordinates(const std::array<double,dim> &position) const override;
251 
255  static
256  void
257  declare_parameters(ParameterHandler &prm);
258 
263  void
264  parse_parameters(ParameterHandler &prm) override;
265 
269  double
270  get_radius(const Point<dim> &position) const;
271 
275  double
276  get_semi_minor_axis_b() const;
277 
281  double
282  get_semi_major_axis_a() const;
283 
284 
288  double
289  get_eccentricity() const;
290 
291 
299  const std::vector<Point<2>> &
300  get_corners() const;
301 
302 
307  get_manifold() const;
308 
309  private:
310  // Declare variables for reading in coordinates of the region of interest.
311  std::vector<Point<2>> corners;
313  double eccentricity;
319  double bottom_depth;
324  // Declare variables for subdividing
325  unsigned int EW_subdiv;
326  unsigned int NS_subdiv;
327  unsigned int depth_subdiv;
328 
329 
330 
335 
336  void
337  set_boundary_ids(parallel::distributed::Triangulation<dim> &coarse_grid) const;
338  };
339  }
340 }
341 
342 
343 #endif
const InitialTopographyModel::Interface< dim > * topography
void declare_parameters(ParameterHandler &prm)
Definition: compat.h:42