22 #ifndef _aspect_geometry_model_spherical_shell_h 23 #define _aspect_geometry_model_spherical_shell_h 30 namespace GeometryModel
46 const double inner_radius,
47 const double outer_radius);
57 virtual std::unique_ptr<Manifold<dim, dim>>
58 clone()
const override;
84 const double w)
const override;
90 virtual Tensor<1, dim>
92 const Point<dim> &x2)
const override;
109 virtual Tensor<1, dim>
111 const typename Triangulation<dim, dim>::face_iterator &face,
112 const Point<dim> &p)
const override;
119 const typename Triangulation<dim, dim>::face_iterator &face,
120 typename Manifold<dim, dim>::FaceVertexNormals &face_vertex_normals)
139 get_new_points(
const ArrayView<
const Point<dim>> &surrounding_points,
140 const Table<2, double> &weights,
141 ArrayView<Point<dim>> new_points)
const override;
149 const ArrayView<const double> &weights)
const override;
160 const std::shared_ptr<const InitialTopographyModel::Interface<dim>>
topo;
202 void create_coarse_mesh (parallel::distributed::Triangulation<dim> &coarse_grid)
const override;
211 double length_scale ()
const override;
226 double depth(
const Point<dim> &position)
const override;
232 double height_above_reference_surface(
const Point<dim> &position)
const override;
237 Point<dim> representative_point(
const double depth)
const override;
242 double maximal_depth()
const override;
265 std::set<types::boundary_id>
266 get_used_boundary_indicators ()
const override;
273 std::map<std::string,types::boundary_id>
274 get_symbolic_boundary_names_map ()
const override;
280 std::set<std::pair<std::pair<types::boundary_id, types::boundary_id>,
unsigned int>>
281 get_periodic_boundary_pairs ()
const override;
290 adjust_positions_for_periodicity (Point<dim> &position,
291 const ArrayView<Point<dim>> &connected_positions = {},
292 const ArrayView<Tensor<1, dim>> &connected_velocities = {})
const override;
300 has_curved_elements()
const override;
308 point_is_in_domain(
const Point<dim> &point)
const override;
321 std::array<double,dim> cartesian_to_natural_coordinates(
const Point<dim> &position)
const override;
328 Point<dim> natural_to_cartesian_coordinates(
const std::array<double,dim> &position)
const override;
342 parse_parameters (ParameterHandler &prm)
override;
348 inner_radius ()
const;
354 outer_radius ()
const;
360 opening_angle ()
const;
367 make_periodicity_constraints(
const DoFHandler<dim> &dof_handler,
368 AffineConstraints<double> &constraints)
const override;
416 void set_manifold_ids (parallel::distributed::Triangulation<dim> &triangulation)
const;
434 std::unique_ptr<const internal::SphericalManifoldWithTopography<dim>>
manifold;
439 static constexpr types::manifold_id my_manifold_id = 99;
virtual Point< dim > get_new_point(const ArrayView< const Point< dim >> &vertices, const ArrayView< const double > &weights) const override
std::unique_ptr< const internal::SphericalManifoldWithTopography< dim > > manifold
virtual Point< dim > get_intermediate_point(const Point< dim > &p1, const Point< dim > &p2, const double w) const override
const std::shared_ptr< const InitialTopographyModel::Interface< dim > > topo
virtual std::unique_ptr< Manifold< dim, dim > > clone() const override
SphericalManifoldWithTopography(const std::shared_ptr< const InitialTopographyModel::Interface< dim >> topography, const double inner_radius, const double outer_radius)
std::vector< double > R_values_list
void declare_parameters(ParameterHandler &prm)
Point< dim > push_forward_from_sphere(const Point< dim > &p) const
int n_cells_along_circumference
virtual void get_new_points(const ArrayView< const Point< dim >> &surrounding_points, const Table< 2, double > &weights, ArrayView< Point< dim >> new_points) const override
virtual void get_normals_at_vertices(const typename Triangulation< dim, dim >::face_iterator &face, typename Manifold< dim, dim >::FaceVertexNormals &face_vertex_normals) const override
CustomMeshRadialSubdivision
Point< dim > pull_back_to_sphere(const Point< dim > &p) const
virtual Tensor< 1, dim > normal_vector(const typename Triangulation< dim, dim >::face_iterator &face, const Point< dim > &p) const override
double topography_for_point(const Point< dim > &x_y_z) const
unsigned int initial_lateral_refinement
virtual Tensor< 1, dim > get_tangent_vector(const Point< dim > &x1, const Point< dim > &x2) const override