ASPECT
two_merged_boxes.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_two_merged_boxes_h
23 #define _aspect_geometry_model_two_merged_boxes_h
24 
27 
28 namespace aspect
29 {
30  namespace GeometryModel
31  {
37  template <int dim>
38  class TwoMergedBoxes : public Interface<dim>, public SimulatorAccess<dim>
39  {
40  public:
41 
45  void create_coarse_mesh (parallel::distributed::Triangulation<dim> &coarse_grid) const override;
46 
51  Point<dim> get_extents () const;
52 
57  Point<dim> get_origin () const;
58 
65  double length_scale () const override;
66 
80  double depth(const Point<dim> &position) const override;
81 
86  double height_above_reference_surface(const Point<dim> &position) const override;
87 
91  Point<dim> representative_point(const double depth) const override;
92 
96  double maximal_depth() const override;
97 
109  std::set<types::boundary_id>
110  get_used_boundary_indicators () const override;
111 
126  std::map<std::string,types::boundary_id>
127  get_symbolic_boundary_names_map () const override;
128 
133  std::set<std::pair<std::pair<types::boundary_id, types::boundary_id>, unsigned int>>
134  get_periodic_boundary_pairs () const override;
135 
142  void
143  adjust_positions_for_periodicity (Point<dim> &position,
144  const ArrayView<Point<dim>> &connected_positions = {},
145  const ArrayView<Tensor<1, dim>> &connected_velocities = {}) const override;
146 
152  bool
153  has_curved_elements() const override;
154 
160  bool
161  point_is_in_domain(const Point<dim> &point) const override;
162 
168 
174  std::array<double,dim> cartesian_to_natural_coordinates(const Point<dim> &position) const override;
175 
181  Point<dim> natural_to_cartesian_coordinates(const std::array<double,dim> &position) const override;
182 
186  static
187  void
188  declare_parameters (ParameterHandler &prm);
189 
193  void
194  parse_parameters (ParameterHandler &prm) override;
195 
196  private:
205 
209  Point<dim> extents;
210 
214  Point<dim> lower_extents;
215 
219  Point<dim> lower_box_origin;
220 
224  Point<dim> upper_extents;
225 
229  Point<dim> upper_box_origin;
230 
235  bool periodic[dim+dim-1];
236 
240  std::array<unsigned int, dim> lower_repetitions;
241 
245  std::array<unsigned int, dim> upper_repetitions;
246 
251  double height_lith;
252 
256  void
257  set_boundary_indicators (parallel::distributed::Triangulation<dim> &triangulation) const;
258  };
259  }
260 }
261 
262 
263 #endif
void create_coarse_mesh(parallel::distributed::Triangulation< dim > &coarse_grid) const override
std::set< std::pair< std::pair< types::boundary_id, types::boundary_id >, unsigned int > > get_periodic_boundary_pairs() const override
double maximal_depth() const override
bool has_curved_elements() const override
double depth(const Point< dim > &position) const override
std::array< unsigned int, dim > upper_repetitions
void set_boundary_indicators(parallel::distributed::Triangulation< dim > &triangulation) const
std::set< types::boundary_id > get_used_boundary_indicators() const override
aspect::Utilities::Coordinates::CoordinateSystem natural_coordinate_system() const override
std::array< unsigned int, dim > lower_repetitions
double height_above_reference_surface(const Point< dim > &position) const override
std::map< std::string, types::boundary_id > get_symbolic_boundary_names_map() const override
void adjust_positions_for_periodicity(Point< dim > &position, const ArrayView< Point< dim >> &connected_positions={}, const ArrayView< Tensor< 1, dim >> &connected_velocities={}) const override
double length_scale() const override
Definition: compat.h:59
bool point_is_in_domain(const Point< dim > &point) const override
void parse_parameters(ParameterHandler &prm) override
Point< dim > natural_to_cartesian_coordinates(const std::array< double, dim > &position) const override
static void declare_parameters(ParameterHandler &prm)
Point< dim > representative_point(const double depth) const override
std::array< double, dim > cartesian_to_natural_coordinates(const Point< dim > &position) const override