ASPECT
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
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 
29 #include <deal.II/grid/manifold.h>
30 #include <deal.II/base/function_lib.h>
31 
32 namespace aspect
33 {
34  namespace GeometryModel
35  {
60  template <int dim>
61  class TwoMergedChunks : public Interface<dim>, public SimulatorAccess<dim>
62  {
63  public:
64 
73  void initialize () override;
74 
78  void create_coarse_mesh (parallel::distributed::Triangulation<dim> &coarse_grid) const override;
79 
90  std::set<types::boundary_id>
91  get_used_boundary_indicators () const override;
92 
100  std::map<std::string,types::boundary_id>
101  get_symbolic_boundary_names_map () const override;
102 
113  double length_scale () const override;
114 
128  double depth(const Point<dim> &position) const override;
129 
134  double height_above_reference_surface(const Point<dim> &position) const override;
135 
139  Point<dim> representative_point(const double depth) const override;
140 
145  virtual
146  double west_longitude() const;
147 
152  virtual
153  double east_longitude() const;
154 
158  virtual
159  double longitude_range() const;
160 
165  virtual
166  double south_latitude() const;
167 
172  virtual
173  double north_latitude() const;
174 
178  virtual
179  double latitude_range() const;
180 
185  double maximal_depth() const override;
186 
190  virtual
191  double inner_radius() const;
192 
196  virtual
197  double outer_radius() const;
198 
199 
205  bool
206  has_curved_elements() const override;
207 
213  bool
214  point_is_in_domain(const Point<dim> &point) const override;
215 
221 
227  std::array<double,dim> cartesian_to_natural_coordinates(const Point<dim> &position) const override;
228 
234  Point<dim> natural_to_cartesian_coordinates(const std::array<double,dim> &position) const override;
235 
239  static
240  void
241  declare_parameters (ParameterHandler &prm);
242 
246  void
247  parse_parameters (ParameterHandler &prm) override;
248 
249  private:
258 
264  Point<dim> point1;
265 
271  Point<dim> point2;
272 
278  Point<dim> point3;
279 
285  Point<dim> point4;
286 
291  std::array<unsigned int, dim> lower_repetitions;
292  std::array<unsigned int, dim> upper_repetitions;
293 
305  std::unique_ptr<const internal::ChunkGeometry<dim>> manifold;
306 
310  static constexpr types::manifold_id my_manifold_id = 15;
311 
315  void set_boundary_indicators (parallel::distributed::Triangulation<dim> &triangulation) const;
316  };
317  }
318 }
319 
320 
321 #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
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