ASPECT
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_chunk_h
23 #define _aspect_geometry_model_chunk_h
24 
27 #include <aspect/compat.h>
28 
29 #include <deal.II/grid/manifold.h>
30 #include <deal.II/base/function_lib.h>
31 #include <deal.II/grid/grid_out.h>
32 
33 namespace aspect
34 {
35  namespace GeometryModel
36  {
37  using namespace dealii;
38 
39  namespace internal
40  {
53  template <int dim>
54  class ChunkGeometry : public ChartManifold<dim,dim>
55  {
56  public:
61  const double min_longitude,
62  const double min_radius,
63  const double max_depth);
64 
68  ChunkGeometry(const ChunkGeometry &other) = default;
69 
75  Point<dim>
76  pull_back(const Point<dim> &space_point) const override;
77 
84  Point<dim>
85  push_forward(const Point<dim> &chart_point) const override;
86 
92  DerivativeForm<1, dim, dim>
93  push_forward_gradient(const Point<dim> &chart_point) const override;
94 
99  Point<dim>
100  pull_back_sphere(const Point<dim> &space_point) const;
101 
107  Point<dim>
108  push_forward_sphere(const Point<dim> &chart_point) const;
109 
113  virtual Tensor<1, dim>
114  normal_vector(
115  const typename Triangulation<dim>::face_iterator &face,
116  const Point<dim> &p) const override;
117 
124  double
125  get_radius(const Point<dim> &space_point) const;
126 
130  std::unique_ptr<Manifold<dim,dim>>
131  clone() const override;
132 
133  private:
138 
142  double point1_lon;
143 
147  double inner_radius;
148 
153  double max_depth;
154 
160  virtual
161  Point<dim>
162  pull_back_topo(const Point<dim> &space_point) const;
163 
169  virtual
170  Point<dim>
171  push_forward_topo(const Point<dim> &chart_point) const;
172  };
173  }
174 
192  template <int dim>
193  class Chunk : public Interface<dim>, public SimulatorAccess<dim>
194  {
195  public:
196 
205  void initialize () override;
206 
210  void create_coarse_mesh (parallel::distributed::Triangulation<dim> &coarse_grid) const override;
211 
212 
223  std::set<types::boundary_id>
224  get_used_boundary_indicators () const override;
225 
238  std::map<std::string,types::boundary_id>
239  get_symbolic_boundary_names_map () const override;
240 
241 
252  double length_scale () const override;
253 
267  double depth(const Point<dim> &position) const override;
268 
273  double height_above_reference_surface(const Point<dim> &position) const override;
274 
278  Point<dim> representative_point(const double depth) const override;
279 
284  virtual
285  double west_longitude() const;
286 
291  virtual
292  double east_longitude() const;
293 
297  virtual
298  double longitude_range() const;
299 
304  virtual
305  double south_latitude() const;
306 
311  virtual
312  double north_latitude() const;
313 
317  virtual
318  double latitude_range() const;
319 
324  double maximal_depth() const override;
325 
329  virtual
330  double inner_radius() const;
331 
335  virtual
336  double outer_radius() const;
337 
338 
344  bool
345  has_curved_elements() const override;
346 
352  bool
353  point_is_in_domain(const Point<dim> &point) const override;
354 
359  aspect::Utilities::Coordinates::CoordinateSystem natural_coordinate_system() const override;
360 
366  std::array<double,dim> cartesian_to_natural_coordinates(const Point<dim> &position) const override;
367 
373  Point<dim> natural_to_cartesian_coordinates(const std::array<double,dim> &position) const override;
374 
378  static
379  void
380  declare_parameters (ParameterHandler &prm);
381 
385  void
386  parse_parameters (ParameterHandler &prm) override;
387 
388  private:
393  Point<dim> point1;
394 
399  Point<dim> point2;
400 
404  std::array<unsigned int, dim> repetitions;
405 
417  std::unique_ptr<const internal::ChunkGeometry<dim>> manifold;
418 
422  static constexpr types::manifold_id my_manifold_id = 15;
423  };
424  }
425 }
426 
427 
428 #endif
std::array< unsigned int, dim > repetitions
Definition: chunk.h:404
const InitialTopographyModel::Interface< dim > * topo
Definition: chunk.h:137
std::unique_ptr< const internal::ChunkGeometry< dim > > manifold
Definition: chunk.h:417
void declare_parameters(ParameterHandler &prm)
Definition: compat.h:59
Definition: compat.h:42