ASPECT
Public Member Functions | Static Public Member Functions | Private Types | Private Member Functions | Private Attributes | Static Private Attributes | List of all members
aspect::GeometryModel::SphericalShell< dim > Class Template Reference
Inheritance diagram for aspect::GeometryModel::SphericalShell< dim >:
Inheritance graph
[legend]

Public Member Functions

 SphericalShell ()=default
 
void initialize () override
 
void create_coarse_mesh (parallel::distributed::Triangulation< dim > &coarse_grid) const override
 
double length_scale () const override
 
double depth (const Point< dim > &position) const override
 
double height_above_reference_surface (const Point< dim > &position) const override
 
Point< dim > representative_point (const double depth) const override
 
double maximal_depth () const override
 
std::set< types::boundary_id > get_used_boundary_indicators () const override
 
std::map< std::string, types::boundary_id > get_symbolic_boundary_names_map () const override
 
std::set< std::pair< std::pair< types::boundary_id, types::boundary_id >, unsigned int > > get_periodic_boundary_pairs () 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
 
bool has_curved_elements () const override
 
bool point_is_in_domain (const Point< dim > &point) 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
 
Point< dim > natural_to_cartesian_coordinates (const std::array< double, dim > &position) const override
 
void parse_parameters (ParameterHandler &prm) override
 
double inner_radius () const
 
double outer_radius () const
 
double opening_angle () const
 
void make_periodicity_constraints (const DoFHandler< dim > &dof_handler, AffineConstraints< double > &constraints) const override
 
- Public Member Functions inherited from aspect::GeometryModel::Interface< dim >
Utilities::NaturalCoordinate< dim > cartesian_to_other_coordinates (const Point< dim > &position, const Utilities::Coordinates::CoordinateSystem &coordinate_system) const
 
types::boundary_id translate_symbolic_boundary_name_to_id (const std::string &name) const
 
std::vector< types::boundary_id > translate_symbolic_boundary_names_to_ids (const std::vector< std::string > &names) const
 
std::string translate_id_to_symbol_name (const types::boundary_id boundary_id) const
 
- Public Member Functions inherited from aspect::Plugins::InterfaceBase
virtual ~InterfaceBase ()=default
 
virtual void update ()
 
- Public Member Functions inherited from aspect::SimulatorAccess< dim >
 SimulatorAccess ()
 
 SimulatorAccess (const Simulator< dim > &simulator_object)
 
virtual ~SimulatorAccess ()=default
 
virtual void initialize_simulator (const Simulator< dim > &simulator_object)
 
const Introspection< dim > & introspection () const
 
const Simulator< dim > & get_simulator () const
 
const Parameters< dim > & get_parameters () const
 
SimulatorSignals< dim > & get_signals () const
 
MPI_Comm get_mpi_communicator () const
 
TimerOutput & get_computing_timer () const
 
const ConditionalOStream & get_pcout () const
 
double get_time () const
 
double get_timestep () const
 
double get_old_timestep () const
 
unsigned int get_timestep_number () const
 
const TimeStepping::Manager< dim > & get_timestepping_manager () const
 
unsigned int get_nonlinear_iteration () const
 
const parallel::distributed::Triangulation< dim > & get_triangulation () const
 
double get_volume () const
 
const Mapping< dim > & get_mapping () const
 
std::string get_output_directory () const
 
bool include_adiabatic_heating () const
 
bool include_latent_heat () const
 
bool include_melt_transport () const
 
int get_stokes_velocity_degree () const
 
double get_adiabatic_surface_temperature () const
 
double get_surface_pressure () const
 
bool convert_output_to_years () const
 
unsigned int get_pre_refinement_step () const
 
unsigned int n_compositional_fields () const
 
double get_end_time () const
 
void get_refinement_criteria (Vector< float > &estimated_error_per_cell) const
 
void get_artificial_viscosity (Vector< float > &viscosity_per_cell, const bool skip_interior_cells=false) const
 
void get_artificial_viscosity_composition (Vector< float > &viscosity_per_cell, const unsigned int compositional_variable) const
 
const LinearAlgebra::BlockVectorget_current_linearization_point () const
 
const LinearAlgebra::BlockVectorget_solution () const
 
const LinearAlgebra::BlockVectorget_old_solution () const
 
const LinearAlgebra::BlockVectorget_old_old_solution () const
 
const LinearAlgebra::BlockVectorget_reaction_vector () const
 
const LinearAlgebra::BlockVectorget_mesh_velocity () const
 
const DoFHandler< dim > & get_dof_handler () const
 
const FiniteElement< dim > & get_fe () const
 
const LinearAlgebra::BlockSparseMatrixget_system_matrix () const
 
const LinearAlgebra::BlockSparseMatrixget_system_preconditioner_matrix () const
 
const MaterialModel::Interface< dim > & get_material_model () const
 
const GravityModel::Interface< dim > & get_gravity_model () const
 
const InitialTopographyModel::Interface< dim > & get_initial_topography_model () const
 
const GeometryModel::Interface< dim > & get_geometry_model () const
 
const AdiabaticConditions::Interface< dim > & get_adiabatic_conditions () const
 
bool has_boundary_temperature () const
 
const BoundaryTemperature::Manager< dim > & get_boundary_temperature_manager () const
 
const BoundaryHeatFlux::Interface< dim > & get_boundary_heat_flux () const
 
bool has_boundary_composition () const
 
const BoundaryComposition::Manager< dim > & get_boundary_composition_manager () const
 
const BoundaryTraction::Manager< dim > & get_boundary_traction_manager () const
 
std::shared_ptr< const InitialTemperature::Manager< dim > > get_initial_temperature_manager_pointer () const
 
const InitialTemperature::Manager< dim > & get_initial_temperature_manager () const
 
std::shared_ptr< const InitialComposition::Manager< dim > > get_initial_composition_manager_pointer () const
 
const InitialComposition::Manager< dim > & get_initial_composition_manager () const
 
const std::set< types::boundary_id > & get_fixed_temperature_boundary_indicators () const
 
const std::set< types::boundary_id > & get_fixed_heat_flux_boundary_indicators () const
 
const std::set< types::boundary_id > & get_fixed_composition_boundary_indicators () const
 
const std::set< types::boundary_id > & get_mesh_deformation_boundary_indicators () const
 
const BoundaryVelocity::Manager< dim > & get_boundary_velocity_manager () const
 
const HeatingModel::Manager< dim > & get_heating_model_manager () const
 
const MeshRefinement::Manager< dim > & get_mesh_refinement_manager () const
 
const MeltHandler< dim > & get_melt_handler () const
 
const VolumeOfFluidHandler< dim > & get_volume_of_fluid_handler () const
 
const NewtonHandler< dim > & get_newton_handler () const
 
const MeshDeformation::MeshDeformationHandler< dim > & get_mesh_deformation_handler () const
 
const LateralAveraging< dim > & get_lateral_averaging () const
 
const AffineConstraints< double > & get_current_constraints () const
 
bool simulator_is_past_initialization () const
 
double get_pressure_scaling () const
 
bool pressure_rhs_needs_compatibility_modification () const
 
bool model_has_prescribed_stokes_solution () const
 
TableHandler & get_statistics_object () const
 
const Postprocess::Manager< dim > & get_postprocess_manager () const
 
unsigned int n_particle_managers () const
 
const Particle::Manager< dim > & get_particle_manager (const unsigned int particle_manager_index) const
 
Particle::Manager< dim > & get_particle_manager (const unsigned int particle_manager_index)
 
bool is_stokes_matrix_free ()
 
const StokesMatrixFreeHandler< dim > & get_stokes_matrix_free () const
 
RotationProperties< dim > compute_net_angular_momentum (const bool use_constant_density, const LinearAlgebra::BlockVector &solution, const bool limit_to_top_faces=false) const
 

Static Public Member Functions

static void declare_parameters (ParameterHandler &prm)
 
- Static Public Member Functions inherited from aspect::Plugins::InterfaceBase
static void declare_parameters (ParameterHandler &prm)
 
- Static Public Member Functions inherited from aspect::SimulatorAccess< dim >
static void get_composition_values_at_q_point (const std::vector< std::vector< double >> &composition_values, const unsigned int q, std::vector< double > &composition_values_at_q_point)
 

Private Types

enum  CustomMeshRadialSubdivision { none, list, slices }
 

Private Member Functions

void set_manifold_ids (parallel::distributed::Triangulation< dim > &triangulation) const
 

Private Attributes

enum aspect::GeometryModel::SphericalShell::CustomMeshRadialSubdivision custom_mesh
 
unsigned int initial_lateral_refinement
 
unsigned int n_slices
 
std::vector< double > R_values_list
 
double R0
 
double R1
 
double phi
 
int n_cells_along_circumference
 
bool periodic
 
std::unique_ptr< const internal::SphericalManifoldWithTopography< dim > > manifold
 

Static Private Attributes

static constexpr types::manifold_id my_manifold_id = 99
 

Detailed Description

template<int dim>
class aspect::GeometryModel::SphericalShell< dim >

A class that describes the geometry as a spherical shell. To be more precise, at least in 2d this class can also simulate just a sector of the spherical shell geometry, in particular a half ring and a quarter ring.

The parameters that describe this geometry and that are read from the input file are the inner and outer radii of the shell and the opening angle of the section of the shell we want to build.

Definition at line 184 of file spherical_shell.h.

Member Enumeration Documentation

§ CustomMeshRadialSubdivision

Specify the radial subdivision of the spherical shell mesh.

Enumerator
none 
list 
slices 

Definition at line 384 of file spherical_shell.h.

Constructor & Destructor Documentation

§ SphericalShell()

template<int dim>
aspect::GeometryModel::SphericalShell< dim >::SphericalShell ( )
default

Constructor.

Member Function Documentation

§ initialize()

template<int dim>
void aspect::GeometryModel::SphericalShell< dim >::initialize ( )
overridevirtual

Initialization function. This function is called once at the beginning of the program after parse_parameters is run and after the SimulatorAccess (if applicable) is initialized. This function calls the initialize function of the manifold with a pointer to the initial topography model obtained from SimulatorAccess.

Reimplemented from aspect::Plugins::InterfaceBase.

§ create_coarse_mesh()

template<int dim>
void aspect::GeometryModel::SphericalShell< dim >::create_coarse_mesh ( parallel::distributed::Triangulation< dim > &  coarse_grid) const
overridevirtual

Generate a coarse mesh for the geometry described by this class.

Implements aspect::GeometryModel::Interface< dim >.

§ length_scale()

template<int dim>
double aspect::GeometryModel::SphericalShell< dim >::length_scale ( ) const
overridevirtual

Return the typical length scale one would expect of features in this geometry, assuming realistic parameters.

As discussed in the step-32 tutorial program, an appropriate length scale for this geometry is 10km, so we return \(10^4\) here.

Implements aspect::GeometryModel::Interface< dim >.

§ depth()

template<int dim>
double aspect::GeometryModel::SphericalShell< dim >::depth ( const Point< dim > &  position) const
overridevirtual

Return the depth that corresponds to the given position. The documentation of the base class (see GeometryModel::Interface::depth()) describes in detail how "depth" is interpreted in general.

Computing a depth requires a geometry model to define a "vertical" direction. The current class considers the radial vector away from the origin as vertical and considers the "outer" boundary as the "surface". In almost all cases one will use a gravity model that also matches these definitions.

Implements aspect::GeometryModel::Interface< dim >.

§ height_above_reference_surface()

template<int dim>
double aspect::GeometryModel::SphericalShell< dim >::height_above_reference_surface ( const Point< dim > &  position) const
overridevirtual

Return the height of the given position relative to the outer radius of the shell.

Implements aspect::GeometryModel::Interface< dim >.

§ representative_point()

template<int dim>
Point<dim> aspect::GeometryModel::SphericalShell< dim >::representative_point ( const double  depth) const
overridevirtual

Returns a representative point for a given depth. Such a point must lie inside the domain for sure (assuming the given depth is valid), but it is not important where exactly in the domain it is. A good choice would be on a middle axis, for example.

The function is used, for example, in computing an initial adiabatic profile. For this, we need to be able to query the density as a function of (adiabatic) pressure and temperature. However, the functions returning the density from the material model also depend on location. Since we are interested only in computing a depth-dependent adiabatic profile without lateral variation, we need to be able to query the density at "representative points" at given depths, without caring too much where exactly that is – but at points that we know for sure are inside the domain.

Implements aspect::GeometryModel::Interface< dim >.

§ maximal_depth()

template<int dim>
double aspect::GeometryModel::SphericalShell< dim >::maximal_depth ( ) const
overridevirtual

Returns the maximal depth of this geometry. For a definition of how "depth" is to be interpreted, see the documentation of the depth() member function. In particular, the maximal depth of a geometry model only represents the depth computed with regard to some reference configuration, ignoring any dynamic or initial topography.

Implements aspect::GeometryModel::Interface< dim >.

§ get_used_boundary_indicators()

template<int dim>
std::set<types::boundary_id> aspect::GeometryModel::SphericalShell< dim >::get_used_boundary_indicators ( ) const
overridevirtual

Return the set of boundary indicators that are used by this model. This information is used to determine what boundary indicators can be used in the input file.

The spherical shell may be generated as per in the original code (with respect to the inner and outer radius, and an initial number of cells along circumference) or following a custom mesh scheme: list of radial values or number of slices. A surface mesh is first generated and refined as desired, before it is extruded radially. A list of radial values subdivides the spherical shell at specified radii. The number of slices subdivides the spherical shell into N slices of equal thickness. The custom spherical shell only works with an opening angle of 360 degrees.

The spherical shell model uses boundary indicators zero and one, with zero corresponding to the inner surface and one corresponding to the outer surface. In 2d, if the geometry is only a slice of the shell, boundary indicators 2 and 3 indicate the left and right radial bounding lines.

Implements aspect::GeometryModel::Interface< dim >.

§ get_symbolic_boundary_names_map()

template<int dim>
std::map<std::string,types::boundary_id> aspect::GeometryModel::SphericalShell< dim >::get_symbolic_boundary_names_map ( ) const
overridevirtual

Return symbolic names for all boundary components. Their names are described in the documentation of this plugin, at the bottom of the .cc file.

Reimplemented from aspect::GeometryModel::Interface< dim >.

§ get_periodic_boundary_pairs()

template<int dim>
std::set<std::pair<std::pair<types::boundary_id, types::boundary_id>, unsigned int> > aspect::GeometryModel::SphericalShell< dim >::get_periodic_boundary_pairs ( ) const
overridevirtual

Return the set of periodic boundaries as described in the input file.

Reimplemented from aspect::GeometryModel::Interface< dim >.

§ adjust_positions_for_periodicity()

template<int dim>
void aspect::GeometryModel::SphericalShell< dim >::adjust_positions_for_periodicity ( Point< dim > &  position,
const ArrayView< Point< dim >> &  connected_positions = {},
const ArrayView< Tensor< 1, dim >> &  connected_velocities = {} 
) const
overridevirtual

Adjust positions to be inside the domain considering periodic boundary conditions.

This function checks if position is outside the domain and if it could have reasonably crossed a periodic boundary. If so, it will adjust the position as described by the periodic boundary (e.g. a translation in a box, or a rotation in a spherical shell). Afterwards, if it adjusted position, it will also adjust all locations given in connected_positions (if any where given) and all velocities given in connected_velocities in the same way. Adjusting both of these allows to adjust related temporary variables, e.g. the intermediate results of an ordinary differential equation solver that are used to compute differences/directions between points. The reason positions and velocities have to be treated separately is that some geometries have to adjust them differently across a periodic boundary, e.g. a box has to translate positions but not velocities, while a spherical shell has to translate positions (by rotation around the center) and rotate velocities.

This function does not check that the position after the adjustment is inside the domain; to check this is the responsibility of the calling function. A common application of this function are particles that crossed a periodic boundary.

Apply a rotation to all points outside of the domain to account for periodicity.

Reimplemented from aspect::GeometryModel::Interface< dim >.

§ has_curved_elements()

template<int dim>
bool aspect::GeometryModel::SphericalShell< dim >::has_curved_elements ( ) const
overridevirtual

If true, the geometry contains cells with boundaries that are not straight and have a deal.II boundary object attached to it. If the return value is false, certain operation can be optimized.The default implementation of this function will return true.

Return true because we have a curved boundary.

Reimplemented from aspect::GeometryModel::Interface< dim >.

§ point_is_in_domain()

template<int dim>
bool aspect::GeometryModel::SphericalShell< dim >::point_is_in_domain ( const Point< dim > &  point) const
overridevirtual

Return whether the given point lies within the domain specified by the geometry. This function does not take into account initial or dynamic surface topography.

Implements aspect::GeometryModel::Interface< dim >.

§ natural_coordinate_system()

template<int dim>
aspect::Utilities::Coordinates::CoordinateSystem aspect::GeometryModel::SphericalShell< dim >::natural_coordinate_system ( ) const
overridevirtual

Returns what the natural coordinate system for this geometry model is, which for a spherical shell is Spherical.

Implements aspect::GeometryModel::Interface< dim >.

§ cartesian_to_natural_coordinates()

template<int dim>
std::array<double,dim> aspect::GeometryModel::SphericalShell< dim >::cartesian_to_natural_coordinates ( const Point< dim > &  position) const
overridevirtual

Takes the Cartesian points (x,z or x,y,z) and returns standardized coordinates which are most 'natural' to the geometry model. For a spherical shell this is (radius, longitude) in 2d and (radius, longitude, latitude) in 3d.

Reimplemented from aspect::GeometryModel::Interface< dim >.

§ natural_to_cartesian_coordinates()

template<int dim>
Point<dim> aspect::GeometryModel::SphericalShell< dim >::natural_to_cartesian_coordinates ( const std::array< double, dim > &  position) const
overridevirtual

Undoes the action of cartesian_to_natural_coordinates, and turns the coordinate system which is most 'natural' to the geometry model into Cartesian coordinates.

Reimplemented from aspect::GeometryModel::Interface< dim >.

§ declare_parameters()

template<int dim>
static void aspect::GeometryModel::SphericalShell< dim >::declare_parameters ( ParameterHandler &  prm)
static

Declare the parameters this class takes through input files. The default implementation of this function does not describe any parameters. Consequently, derived classes do not have to overload this function if they do not take any runtime parameters.

§ parse_parameters()

template<int dim>
void aspect::GeometryModel::SphericalShell< dim >::parse_parameters ( ParameterHandler &  prm)
overridevirtual

Read the parameters this class declares from the parameter file. The default implementation of this function does not read any parameters. Consequently, derived classes do not have to overload this function if they do not take any runtime parameters.

Reimplemented from aspect::Plugins::InterfaceBase.

§ inner_radius()

template<int dim>
double aspect::GeometryModel::SphericalShell< dim >::inner_radius ( ) const

Return the inner radius of the shell.

§ outer_radius()

template<int dim>
double aspect::GeometryModel::SphericalShell< dim >::outer_radius ( ) const

Return the outer radius of the shell.

§ opening_angle()

template<int dim>
double aspect::GeometryModel::SphericalShell< dim >::opening_angle ( ) const

Return the opening angle of the shell sector.

§ make_periodicity_constraints()

template<int dim>
void aspect::GeometryModel::SphericalShell< dim >::make_periodicity_constraints ( const DoFHandler< dim > &  dof_handler,
AffineConstraints< double > &  constraints 
) const
overridevirtual

Collects periodic boundaries constraints for the given geometry, which will be added to the existing constraints.

Reimplemented from aspect::GeometryModel::Interface< dim >.

§ set_manifold_ids()

template<int dim>
void aspect::GeometryModel::SphericalShell< dim >::set_manifold_ids ( parallel::distributed::Triangulation< dim > &  triangulation) const
private

Set the manifold ids on all cells (also boundaries) before refinement to generate well shaped cells.

Member Data Documentation

§ custom_mesh

§ initial_lateral_refinement

template<int dim>
unsigned int aspect::GeometryModel::SphericalShell< dim >::initial_lateral_refinement
private

Initial surface refinement for the custom mesh cases.

Definition at line 394 of file spherical_shell.h.

§ n_slices

template<int dim>
unsigned int aspect::GeometryModel::SphericalShell< dim >::n_slices
private

Initial surface refinement for the custom mesh cases.

Definition at line 399 of file spherical_shell.h.

§ R_values_list

template<int dim>
std::vector<double> aspect::GeometryModel::SphericalShell< dim >::R_values_list
private

List of radial values for the list custom mesh.

Definition at line 404 of file spherical_shell.h.

§ R0

template<int dim>
double aspect::GeometryModel::SphericalShell< dim >::R0
private

Inner and outer radii of the spherical shell.

Definition at line 409 of file spherical_shell.h.

§ R1

template<int dim>
double aspect::GeometryModel::SphericalShell< dim >::R1
private

Definition at line 409 of file spherical_shell.h.

§ phi

template<int dim>
double aspect::GeometryModel::SphericalShell< dim >::phi
private

Opening angle of the section of the shell that we simulate.

Definition at line 414 of file spherical_shell.h.

§ n_cells_along_circumference

template<int dim>
int aspect::GeometryModel::SphericalShell< dim >::n_cells_along_circumference
private

Number of tangential mesh cells in the initial, coarse mesh.

Definition at line 419 of file spherical_shell.h.

§ periodic

template<int dim>
bool aspect::GeometryModel::SphericalShell< dim >::periodic
private

Flag whether the 2D quarter shell is periodic in phi.

Definition at line 430 of file spherical_shell.h.

§ manifold

template<int dim>
std::unique_ptr<const internal::SphericalManifoldWithTopography<dim> > aspect::GeometryModel::SphericalShell< dim >::manifold
private

An object that describes the geometry. This pointer is initialized in the initialize() function, and serves as the manifold object that the triangulation is later given in create_coarse_mesh() where the triangulation clones it.

The object is marked as 'const' to make it clear that it should not be modified once created. That is because the triangulation copies it, and modifying the current object will not have any impact on the manifold used by the triangulation.

Definition at line 443 of file spherical_shell.h.

§ my_manifold_id

template<int dim>
constexpr types::manifold_id aspect::GeometryModel::SphericalShell< dim >::my_manifold_id = 99
staticprivate

Give a symbolic name to the manifold id to be used by this class.

Definition at line 448 of file spherical_shell.h.


The documentation for this class was generated from the following file: