ASPECT
two_merged_chunks.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_two_merged_chunks_h
23 #define _aspect_geometry_model_two_merged_chunks_h
24 
28 #include <aspect/compat.h>
29 
30 #include <deal.II/grid/manifold.h>
31 #include <deal.II/base/function_lib.h>
32 
33 namespace aspect
34 {
35  namespace GeometryModel
36  {
61  template <int dim>
62  class TwoMergedChunks : public Interface<dim>, public SimulatorAccess<dim>
63  {
64  public:
65 
74  void initialize () override;
75 
79  void create_coarse_mesh (parallel::distributed::Triangulation<dim> &coarse_grid) const override;
80 
91  std::set<types::boundary_id>
92  get_used_boundary_indicators () const override;
93 
101  std::map<std::string,types::boundary_id>
102  get_symbolic_boundary_names_map () const override;
103 
114  double length_scale () const override;
115 
129  double depth(const Point<dim> &position) const override;
130 
135  double height_above_reference_surface(const Point<dim> &position) const override;
136 
140  Point<dim> representative_point(const double depth) const override;
141 
146  virtual
147  double west_longitude() const;
148 
153  virtual
154  double east_longitude() const;
155 
159  virtual
160  double longitude_range() const;
161 
166  virtual
167  double south_latitude() const;
168 
173  virtual
174  double north_latitude() const;
175 
179  virtual
180  double latitude_range() const;
181 
186  double maximal_depth() const override;
187 
191  virtual
192  double inner_radius() const;
193 
197  virtual
198  double outer_radius() const;
199 
200 
206  bool
207  has_curved_elements() const override;
208 
214  bool
215  point_is_in_domain(const Point<dim> &point) const override;
216 
222 
228  std::array<double,dim> cartesian_to_natural_coordinates(const Point<dim> &position) const override;
229 
235  Point<dim> natural_to_cartesian_coordinates(const std::array<double,dim> &position) const override;
236 
240  static
241  void
242  declare_parameters (ParameterHandler &prm);
243 
247  void
248  parse_parameters (ParameterHandler &prm) override;
249 
250  private:
259 
265  Point<dim> point1;
266 
272  Point<dim> point2;
273 
279  Point<dim> point3;
280 
286  Point<dim> point4;
287 
292  std::array<unsigned int, dim> lower_repetitions;
293  std::array<unsigned int, dim> upper_repetitions;
294 
306  std::unique_ptr<const internal::ChunkGeometry<dim>> manifold;
307 
311  static constexpr types::manifold_id my_manifold_id = 15;
312 
316  void set_boundary_indicators (parallel::distributed::Triangulation<dim> &triangulation) const;
317  };
318  }
319 }
320 
321 
322 #endif
virtual double latitude_range() const
virtual double outer_radius() const
Point< dim > representative_point(const double depth) const override
aspect::Utilities::Coordinates::CoordinateSystem natural_coordinate_system() const override
std::array< double, dim > cartesian_to_natural_coordinates(const Point< dim > &position) const override
bool point_is_in_domain(const Point< dim > &point) const override
double length_scale() const override
bool has_curved_elements() const override
virtual double south_latitude() const
static void declare_parameters(ParameterHandler &prm)
virtual double longitude_range() const
virtual double inner_radius() const
virtual double east_longitude() const
virtual double north_latitude() const
static constexpr types::manifold_id my_manifold_id
double maximal_depth() const override
void parse_parameters(ParameterHandler &prm) override
Definition: compat.h:59
std::array< unsigned int, dim > lower_repetitions
std::unique_ptr< const internal::ChunkGeometry< dim > > manifold
std::array< unsigned int, dim > upper_repetitions
void create_coarse_mesh(parallel::distributed::Triangulation< dim > &coarse_grid) const override
std::map< std::string, types::boundary_id > get_symbolic_boundary_names_map() const override
double depth(const Point< dim > &position) const override
double height_above_reference_surface(const Point< dim > &position) const override
void set_boundary_indicators(parallel::distributed::Triangulation< dim > &triangulation) const
Point< dim > natural_to_cartesian_coordinates(const std::array< double, dim > &position) const override
std::set< types::boundary_id > get_used_boundary_indicators() const override
virtual double west_longitude() const