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

Public Member Functions

void initialize () override
 
void create_coarse_mesh (parallel::distributed::Triangulation< dim > &coarse_grid) 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
 
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
 
virtual double west_longitude () const
 
virtual double east_longitude () const
 
virtual double longitude_range () const
 
virtual double south_latitude () const
 
virtual double north_latitude () const
 
virtual double latitude_range () const
 
double maximal_depth () const override
 
virtual double inner_radius () const
 
virtual double outer_radius () const
 
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
 
- 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
 
virtual std::set< std::pair< std::pair< types::boundary_id, types::boundary_id >, unsigned int > > get_periodic_boundary_pairs () const
 
virtual void adjust_positions_for_periodicity (Point< dim > &position, const ArrayView< Point< dim >> &connected_positions={}, const ArrayView< Tensor< 1, dim >> &connected_velocities={}) const
 
virtual void make_periodicity_constraints (const DoFHandler< dim > &dof_handler, AffineConstraints< double > &constraints) 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 Member Functions

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

Private Attributes

bool use_merged_grids
 
Point< dim > point1
 
Point< dim > point2
 
Point< dim > point3
 
Point< dim > point4
 
std::array< unsigned int, dim > lower_repetitions
 
std::array< unsigned int, dim > upper_repetitions
 
std::unique_ptr< const internal::ChunkGeometry< dim > > manifold
 

Static Private Attributes

static constexpr types::manifold_id my_manifold_id = 15
 

Detailed Description

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

A geometry model class that describes a chunk of a spherical shell, but with two boundary indicators per side boundary. This allows for specifying two different types of boundary conditions on the side boundaries, for example prescribed plate motions along the edge of the lithosphere and an open boundary underneath. In 2D, this class simulates a sector with inner and outer radius, and minimum and maximum longitude. Longitude increases anticlockwise from the positive x-axis, as per the mathematical convention of phi. In 3D, the class simulates a chunk of a sphere, bounded by arbitrary lines of longitude, latitude and radius. Boundary indicator names are lowerwest, upperwest, lowereast, uppereast, lowersouth, uppersouth, lowernorth, uppernorth, bottom and top.

The parameters that describe this geometry and that are read from the input file are the inner, outer and middle boundary radius of the shell, the minimum and maximum longitude, minimum and maximum longitude, and the number of cells initialized in each dimension. In the radial direction, this number should be specified for both the lower and the upper part of the chunk (with their boundary lying at the middle boundary radius).

Initial topography can be added through a radial displacement of the mesh nodes.

Definition at line 65 of file two_merged_chunks.h.

Member Function Documentation

§ initialize()

template<int dim>
void aspect::GeometryModel::TwoMergedChunks< 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::TwoMergedChunks< 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 >.

§ get_used_boundary_indicators()

template<int dim>
std::set<types::boundary_id> aspect::GeometryModel::TwoMergedChunks< 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 chunk model uses boundary indicators zero through 2*dim+2*(dim-1)-1, with the first two being the faces perpendicular to the radius of the shell, the next four along lines of longitude and, in 3d, the next four along lines of latitude.

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

§ get_symbolic_boundary_names_map()

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

Return a mapping from symbolic names of each part of the boundary to the corresponding boundary indicator. This allows users to specify names, not just numbers in their input files when describing which parts of the boundary have to satisfy which boundary conditions.

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

§ length_scale()

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

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

As described in the first ASPECT paper, a length scale of 10km = 1e4m works well for the pressure scaling for earth sized spherical shells. use a length scale that yields this value for the R0,R1 corresponding to earth but otherwise scales like (R1-R0)

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

§ depth()

template<int dim>
double aspect::GeometryModel::TwoMergedChunks< 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 "top" 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::TwoMergedChunks< dim >::height_above_reference_surface ( const Point< dim > &  position) const
overridevirtual

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

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

§ representative_point()

template<int dim>
Point<dim> aspect::GeometryModel::TwoMergedChunks< 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 >.

§ west_longitude()

template<int dim>
virtual double aspect::GeometryModel::TwoMergedChunks< dim >::west_longitude ( ) const
virtual

Return the longitude at the western edge of the chunk measured in radians.

§ east_longitude()

template<int dim>
virtual double aspect::GeometryModel::TwoMergedChunks< dim >::east_longitude ( ) const
virtual

Return the longitude at the eastern edge of the chunk measured in radians.

§ longitude_range()

template<int dim>
virtual double aspect::GeometryModel::TwoMergedChunks< dim >::longitude_range ( ) const
virtual

Return the longitude range of the chunk measured in radians.

§ south_latitude()

template<int dim>
virtual double aspect::GeometryModel::TwoMergedChunks< dim >::south_latitude ( ) const
virtual

Return the latitude at the southern edge of the chunk measured in radians from the equator.

§ north_latitude()

template<int dim>
virtual double aspect::GeometryModel::TwoMergedChunks< dim >::north_latitude ( ) const
virtual

Return the latitude at the northern edge of the chunk measured in radians from the equator.

§ latitude_range()

template<int dim>
virtual double aspect::GeometryModel::TwoMergedChunks< dim >::latitude_range ( ) const
virtual

Return the latitude range of the chunk measured in radians

§ maximal_depth()

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

Return the maximum depth from the surface of the model measured in meters.

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

§ inner_radius()

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

Return the inner radius of the chunk measured in meters.

§ outer_radius()

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

Return the outer radius of the chunk measured in meters.

§ has_curved_elements()

template<int dim>
bool aspect::GeometryModel::TwoMergedChunks< 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.

A chunk has curved boundaries and cells, so return true.

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

§ point_is_in_domain()

template<int dim>
bool aspect::GeometryModel::TwoMergedChunks< 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::TwoMergedChunks< dim >::natural_coordinate_system ( ) const
overridevirtual

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

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

§ cartesian_to_natural_coordinates()

template<int dim>
std::array<double,dim> aspect::GeometryModel::TwoMergedChunks< 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 chunk 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::TwoMergedChunks< 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::TwoMergedChunks< dim >::declare_parameters ( ParameterHandler &  prm)
static

Declare the parameters this class takes through input files.

§ parse_parameters()

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

Read the parameters this class declares from the parameter file.

Reimplemented from aspect::Plugins::InterfaceBase.

§ set_boundary_indicators()

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

Bind boundary indicators to child cells after each mesh refinement round.

Member Data Documentation

§ use_merged_grids

template<int dim>
bool aspect::GeometryModel::TwoMergedChunks< dim >::use_merged_grids
private

Whether to make the grid by gluing together two chunks, or just use one chunk to make the grid. Using two grids glued together is a safer option, since it forces the boundary conditions to be always applied to the same depth, but one unified grid allows for a more flexible usage of the adaptive refinement.

Definition at line 261 of file two_merged_chunks.h.

§ point1

template<int dim>
Point<dim> aspect::GeometryModel::TwoMergedChunks< dim >::point1
private

Minimum longitude-depth (2D) or longitude-latitude-depth (3D) point of the entire merged chunk.

Definition at line 268 of file two_merged_chunks.h.

§ point2

template<int dim>
Point<dim> aspect::GeometryModel::TwoMergedChunks< dim >::point2
private

Maximum longitude-depth (2D) or longitude-latitude-depth (3D) point of the entire merged chunk.

Definition at line 275 of file two_merged_chunks.h.

§ point3

template<int dim>
Point<dim> aspect::GeometryModel::TwoMergedChunks< dim >::point3
private

Minimum longitude-depth (2D) or longitude-latitude-depth (3D) point for the upper chunk.

Definition at line 282 of file two_merged_chunks.h.

§ point4

template<int dim>
Point<dim> aspect::GeometryModel::TwoMergedChunks< dim >::point4
private

Maximum longitude-depth (2D) or longitude-latitude-depth (3D) point for the lower chunk.

Definition at line 289 of file two_merged_chunks.h.

§ lower_repetitions

template<int dim>
std::array<unsigned int, dim> aspect::GeometryModel::TwoMergedChunks< dim >::lower_repetitions
private

The number of cells in each coordinate direction for the lower and upper chunk.

Definition at line 295 of file two_merged_chunks.h.

§ upper_repetitions

template<int dim>
std::array<unsigned int, dim> aspect::GeometryModel::TwoMergedChunks< dim >::upper_repetitions
private

Definition at line 296 of file two_merged_chunks.h.

§ manifold

template<int dim>
std::unique_ptr<const internal::ChunkGeometry<dim> > aspect::GeometryModel::TwoMergedChunks< 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 309 of file two_merged_chunks.h.

§ my_manifold_id

template<int dim>
constexpr types::manifold_id aspect::GeometryModel::TwoMergedChunks< dim >::my_manifold_id = 15
staticprivate

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

Definition at line 314 of file two_merged_chunks.h.


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