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  {
37  using namespace dealii;
38 
39 
64  template <int dim>
65  class TwoMergedChunks : public Interface<dim>, public SimulatorAccess<dim>
66  {
67  public:
68 
77  void initialize () override;
78 
82  void create_coarse_mesh (parallel::distributed::Triangulation<dim> &coarse_grid) const override;
83 
94  std::set<types::boundary_id>
95  get_used_boundary_indicators () const override;
96 
104  std::map<std::string,types::boundary_id>
105  get_symbolic_boundary_names_map () const override;
106 
117  double length_scale () const override;
118 
132  double depth(const Point<dim> &position) const override;
133 
138  double height_above_reference_surface(const Point<dim> &position) const override;
139 
143  Point<dim> representative_point(const double depth) const override;
144 
149  virtual
150  double west_longitude() const;
151 
156  virtual
157  double east_longitude() const;
158 
162  virtual
163  double longitude_range() const;
164 
169  virtual
170  double south_latitude() const;
171 
176  virtual
177  double north_latitude() const;
178 
182  virtual
183  double latitude_range() const;
184 
189  double maximal_depth() const override;
190 
194  virtual
195  double inner_radius() const;
196 
200  virtual
201  double outer_radius() const;
202 
203 
209  bool
210  has_curved_elements() const override;
211 
217  bool
218  point_is_in_domain(const Point<dim> &point) const override;
219 
224  aspect::Utilities::Coordinates::CoordinateSystem natural_coordinate_system() const override;
225 
231  std::array<double,dim> cartesian_to_natural_coordinates(const Point<dim> &position) const override;
232 
238  Point<dim> natural_to_cartesian_coordinates(const std::array<double,dim> &position) const override;
239 
243  static
244  void
245  declare_parameters (ParameterHandler &prm);
246 
250  void
251  parse_parameters (ParameterHandler &prm) override;
252 
253  private:
262 
268  Point<dim> point1;
269 
275  Point<dim> point2;
276 
282  Point<dim> point3;
283 
289  Point<dim> point4;
290 
295  std::array<unsigned int, dim> lower_repetitions;
296  std::array<unsigned int, dim> upper_repetitions;
297 
309  std::unique_ptr<const internal::ChunkGeometry<dim>> manifold;
310 
314  static constexpr types::manifold_id my_manifold_id = 15;
315 
319  void set_boundary_indicators (parallel::distributed::Triangulation<dim> &triangulation) const;
320  };
321  }
322 }
323 
324 
325 #endif
void declare_parameters(ParameterHandler &prm)
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
Definition: compat.h:42