![]() |
ASPECT
|
Public Member Functions | |
SphericalShell () | |
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 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 |
![]() | |
virtual | ~Interface ()=default |
virtual void | initialize () |
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 |
![]() | |
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 |
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 |
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::BlockVector & | get_current_linearization_point () const |
const LinearAlgebra::BlockVector & | get_solution () const |
const LinearAlgebra::BlockVector & | get_old_solution () const |
const LinearAlgebra::BlockVector & | get_old_old_solution () const |
const LinearAlgebra::BlockVector & | get_reaction_vector () const |
const LinearAlgebra::BlockVector & | get_mesh_velocity () const |
const DoFHandler< dim > & | get_dof_handler () const |
const FiniteElement< dim > & | get_fe () const |
const LinearAlgebra::BlockSparseMatrix & | get_system_matrix () const |
const LinearAlgebra::BlockSparseMatrix & | get_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 std::map< types::boundary_id, std::unique_ptr< BoundaryTraction::Interface< dim > > > & | get_boundary_traction () 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 |
const Particle::World< dim > & | get_particle_world () const |
Particle::World< dim > & | get_particle_world () |
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 void | declare_parameters (ParameterHandler &prm) |
![]() | |
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 |
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 |
const SphericalManifold< dim > | spherical_manifold |
bool | periodic |
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 47 of file spherical_shell.h.
|
private |
Specify the radial subdivision of the spherical shell mesh.
Enumerator | |
---|---|
none | |
list | |
slices |
Definition at line 236 of file spherical_shell.h.
aspect::GeometryModel::SphericalShell< dim >::SphericalShell | ( | ) |
Constructor.
|
overridevirtual |
Generate a coarse mesh for the geometry described by this class.
Implements aspect::GeometryModel::Interface< dim >.
|
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 >.
|
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 >.
|
overridevirtual |
Return the height of the given position relative to the outer radius of the shell.
Implements aspect::GeometryModel::Interface< dim >.
|
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 >.
|
overridevirtual |
Returns the maximal depth of this geometry.
Implements aspect::GeometryModel::Interface< dim >.
|
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 >.
|
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 >.
|
overridevirtual |
Return the set of periodic boundaries as described in the input file.
Reimplemented from aspect::GeometryModel::Interface< dim >.
|
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) in the same way. Adjusting connected_positions
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.
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 >.
|
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 >.
|
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 >.
|
overridevirtual |
Returns what the natural coordinate system for this geometry model is, which for a spherical shell is Spherical.
Implements aspect::GeometryModel::Interface< dim >.
|
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 >.
|
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 >.
|
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.
|
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::GeometryModel::Interface< dim >.
double aspect::GeometryModel::SphericalShell< dim >::inner_radius | ( | ) | const |
Return the inner radius of the shell.
double aspect::GeometryModel::SphericalShell< dim >::outer_radius | ( | ) | const |
Return the outer radius of the shell.
double aspect::GeometryModel::SphericalShell< dim >::opening_angle | ( | ) | const |
Return the opening angle of the shell sector.
|
overridevirtual |
Collects periodic boundaries constraints for the given geometry, which will be added to the existing constraints
.
Reimplemented from aspect::GeometryModel::Interface< dim >.
|
private |
Set the manifold ids on all cells (also boundaries) before refinement to generate well shaped cells.
|
private |
|
private |
Initial surface refinement for the custom mesh cases.
Definition at line 246 of file spherical_shell.h.
|
private |
Initial surface refinement for the custom mesh cases.
Definition at line 251 of file spherical_shell.h.
|
private |
List of radial values for the list custom mesh.
Definition at line 256 of file spherical_shell.h.
|
private |
Inner and outer radii of the spherical shell.
Definition at line 261 of file spherical_shell.h.
|
private |
Definition at line 261 of file spherical_shell.h.
|
private |
Opening angle of the section of the shell that we simulate.
Definition at line 266 of file spherical_shell.h.
|
private |
Number of tangential mesh cells in the initial, coarse mesh.
Definition at line 271 of file spherical_shell.h.
|
private |
The manifold that describes the geometry.
Definition at line 276 of file spherical_shell.h.
|
private |
Flag whether the 2D quarter shell is periodic in phi.
Definition at line 287 of file spherical_shell.h.