ASPECT
chunk.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2011 - 2023 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:
60  ChunkGeometry();
61 
65  ChunkGeometry(const ChunkGeometry &other);
66 
67  /*
68  * An initialization function to make sure that the
69  * manifold has access to the topography plugins.
70  */
71  void
73 
79  Point<dim>
80  pull_back(const Point<dim> &space_point) const override;
81 
88  Point<dim>
89  push_forward(const Point<dim> &chart_point) const override;
90 
96  DerivativeForm<1, dim, dim>
97  push_forward_gradient(const Point<dim> &chart_point) const override;
98 
103  Point<dim>
104  pull_back_sphere(const Point<dim> &space_point) const;
105 
111  Point<dim>
112  push_forward_sphere(const Point<dim> &chart_point) const;
113 
117  virtual Tensor<1, dim>
118  normal_vector(
119  const typename Triangulation<dim>::face_iterator &face,
120  const Point<dim> &p) const override;
121 
128  double
129  get_radius(const Point<dim> &space_point) const;
130 
137  virtual
138  void
139  set_min_longitude(const double p1_lon);
140 
144  std::unique_ptr<Manifold<dim,dim>>
145  clone() const override;
146 
150  void
151  set_min_radius(const double p1_rad);
152 
156  void
157  set_max_depth(const double p2_rad_minus_p1_rad);
158 
159  private:
163  double point1_lon;
164 
168  double inner_radius;
169 
174  double max_depth;
175 
181  virtual
182  Point<dim>
183  pull_back_topo(const Point<dim> &space_point) const;
184 
190  virtual
191  Point<dim>
192  push_forward_topo(const Point<dim> &chart_point) const;
193 
198  };
199  }
200 
218  template <int dim>
219  class Chunk : public Interface<dim>, public SimulatorAccess<dim>
220  {
221  public:
222 
231  void initialize () override;
232 
237  void set_topography_model (const InitialTopographyModel::Interface<dim> *topo_pointer);
238 
242  void create_coarse_mesh (parallel::distributed::Triangulation<dim> &coarse_grid) const override;
243 
244 
255  std::set<types::boundary_id>
256  get_used_boundary_indicators () const override;
257 
270  std::map<std::string,types::boundary_id>
271  get_symbolic_boundary_names_map () const override;
272 
273 
284  double length_scale () const override;
285 
299  double depth(const Point<dim> &position) const override;
300 
305  double height_above_reference_surface(const Point<dim> &position) const override;
306 
310  Point<dim> representative_point(const double depth) const override;
311 
325  double depth_wrt_topo(const Point<dim> &position) const;
326 
331  virtual
332  double west_longitude() const;
333 
338  virtual
339  double east_longitude() const;
340 
344  virtual
345  double longitude_range() const;
346 
351  virtual
352  double south_latitude() const;
353 
358  virtual
359  double north_latitude() const;
360 
364  virtual
365  double latitude_range() const;
366 
371  double maximal_depth() const override;
372 
376  virtual
377  double inner_radius() const;
378 
382  virtual
383  double outer_radius() const;
384 
385 
391  bool
392  has_curved_elements() const override;
393 
399  bool
400  point_is_in_domain(const Point<dim> &point) const override;
401 
406  aspect::Utilities::Coordinates::CoordinateSystem natural_coordinate_system() const override;
407 
413  std::array<double,dim> cartesian_to_natural_coordinates(const Point<dim> &position) const override;
414 
420  Point<dim> natural_to_cartesian_coordinates(const std::array<double,dim> &position) const override;
421 
425  static
426  void
427  declare_parameters (ParameterHandler &prm);
428 
432  void
433  parse_parameters (ParameterHandler &prm) override;
434 
435  private:
440  Point<dim> point1;
441 
446  Point<dim> point2;
447 
451  std::array<unsigned int, dim> repetitions;
452 
457  };
458  }
459 }
460 
461 
462 #endif
std::array< unsigned int, dim > repetitions
Definition: chunk.h:451
const InitialTopographyModel::Interface< dim > * topo
Definition: chunk.h:197
internal::ChunkGeometry< dim > manifold
Definition: chunk.h:456
void declare_parameters(ParameterHandler &prm)
Definition: compat.h:42