ASPECT
|
Public Member Functions | |
virtual void | create_coarse_mesh (parallel::distributed::Triangulation< dim > &coarse_grid) const =0 |
virtual double | length_scale () const =0 |
virtual double | depth (const Point< dim > &position) const =0 |
virtual double | height_above_reference_surface (const Point< dim > &position) const =0 |
Utilities::NaturalCoordinate< dim > | cartesian_to_other_coordinates (const Point< dim > &position, const Utilities::Coordinates::CoordinateSystem &coordinate_system) const |
virtual aspect::Utilities::Coordinates::CoordinateSystem | natural_coordinate_system () const =0 |
virtual std::array< double, dim > | cartesian_to_natural_coordinates (const Point< dim > &position) const |
virtual Point< dim > | natural_to_cartesian_coordinates (const std::array< double, dim > &position) const |
virtual Point< dim > | representative_point (const double depth) const =0 |
virtual double | maximal_depth () const =0 |
virtual std::set< types::boundary_id > | get_used_boundary_indicators () const =0 |
virtual std::map< std::string, types::boundary_id > | get_symbolic_boundary_names_map () 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 bool | has_curved_elements () const |
virtual bool | point_is_in_domain (const Point< dim > &p) const =0 |
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 | initialize () |
virtual void | update () |
virtual void | parse_parameters (ParameterHandler &prm) |
Additional Inherited Members | |
Static Public Member Functions inherited from aspect::Plugins::InterfaceBase | |
static void | declare_parameters (ParameterHandler &prm) |
Base class for classes that describe particular geometries for the domain. These classes must also be able to create coarse meshes and describe what kinds of boundary conditions hold where.
Definition at line 41 of file parameters.h.
|
pure virtual |
Generate a coarse mesh for the geometry described by this class.
Implemented in aspect::GeometryModel::Chunk< dim >, aspect::GeometryModel::SphericalShell< dim >, aspect::GeometryModel::EllipsoidalChunk< dim >, aspect::GeometryModel::TwoMergedChunks< dim >, aspect::GeometryModel::Box< dim >, aspect::GeometryModel::TwoMergedBoxes< dim >, and aspect::GeometryModel::Sphere< dim >.
|
pure virtual |
Return the typical length scale one would expect of features in this geometry, assuming realistic parameters.
The result of this function is used in computing the scaling factor for the pressure in the Stokes equation. There, the scaling factor is chosen as the ratio of the reference viscosity divided by the length scale. As an example, in the step-32 tutorial program we have determined that a suitable length scale for scaling is 10km, in a domain that is 12,000km across. This length scale suitably matches the order of magnitude for the diameter of plumes in the earth.
Implemented in aspect::GeometryModel::Chunk< dim >, aspect::GeometryModel::SphericalShell< dim >, aspect::GeometryModel::EllipsoidalChunk< dim >, aspect::GeometryModel::TwoMergedChunks< dim >, aspect::GeometryModel::Box< dim >, aspect::GeometryModel::TwoMergedBoxes< dim >, and aspect::GeometryModel::Sphere< dim >.
|
pure virtual |
Return the depth that corresponds to the given position. The returned value is between 0 and maximal_depth(), where points with a zero depth correspond to the "surface" of the model.
Computing a depth requires a geometry model to define a "vertical" direction. For example, the Box model considers the \((0,1)^T\) vector in 2d (and the \((0,0,1)^T\) vector in 3d) as vertical and considers the "top" boundary as the "surface". Similarly, the spherical shell takes the radial vector as "vertical" and the outer boundary as "surface". In most cases, how geometry models define "vertical" and "surface" will be intuitive and obvious. In almost all cases one will use a gravity model that also matches these definitions.
Implemented in aspect::GeometryModel::Chunk< dim >, aspect::GeometryModel::SphericalShell< dim >, aspect::GeometryModel::EllipsoidalChunk< dim >, aspect::GeometryModel::TwoMergedChunks< dim >, aspect::GeometryModel::Box< dim >, aspect::GeometryModel::TwoMergedBoxes< dim >, and aspect::GeometryModel::Sphere< dim >.
|
pure virtual |
Return the height of the given position relative to the reference surface of the model. Positive returned value means that the point is above (i.e., farther from the center of the model) the reference surface, negative value means that the point is below the the reference surface.
Same limitations as for the depth function, apply here.
Implemented in aspect::GeometryModel::Chunk< dim >, aspect::GeometryModel::SphericalShell< dim >, aspect::GeometryModel::EllipsoidalChunk< dim >, aspect::GeometryModel::TwoMergedChunks< dim >, aspect::GeometryModel::Box< dim >, aspect::GeometryModel::TwoMergedBoxes< dim >, and aspect::GeometryModel::Sphere< dim >.
Utilities::NaturalCoordinate<dim> aspect::GeometryModel::Interface< dim >::cartesian_to_other_coordinates | ( | const Point< dim > & | position, |
const Utilities::Coordinates::CoordinateSystem & | coordinate_system | ||
) | const |
Converts a Cartesian Point into another coordinate system and returns it as a NaturalCoordinate.
|
pure virtual |
Returns what the natural coordinate system for this geometry model is.
Implemented in aspect::GeometryModel::Chunk< dim >, aspect::GeometryModel::SphericalShell< dim >, aspect::GeometryModel::EllipsoidalChunk< dim >, aspect::GeometryModel::TwoMergedChunks< dim >, aspect::GeometryModel::Box< dim >, aspect::GeometryModel::TwoMergedBoxes< dim >, and aspect::GeometryModel::Sphere< dim >.
|
virtual |
Takes the Cartesian points (x,z or x,y,z) and returns standardized coordinates which are most 'natural' to the geometry model. For a box this will be (x,z) in 2d or (x,y,z) in 3d, and for a spheroid geometry model it will be (radius, longitude) in 2d and (radius, longitude, latitude) in 3d.
Reimplemented in aspect::GeometryModel::Chunk< dim >, aspect::GeometryModel::SphericalShell< dim >, aspect::GeometryModel::EllipsoidalChunk< dim >, aspect::GeometryModel::TwoMergedChunks< dim >, aspect::GeometryModel::Box< dim >, aspect::GeometryModel::TwoMergedBoxes< dim >, and aspect::GeometryModel::Sphere< dim >.
|
virtual |
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 in aspect::GeometryModel::Chunk< dim >, aspect::GeometryModel::SphericalShell< dim >, aspect::GeometryModel::EllipsoidalChunk< dim >, aspect::GeometryModel::TwoMergedChunks< dim >, aspect::GeometryModel::Box< dim >, aspect::GeometryModel::TwoMergedBoxes< dim >, and aspect::GeometryModel::Sphere< dim >.
|
pure virtual |
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.
Implemented in aspect::GeometryModel::Chunk< dim >, aspect::GeometryModel::SphericalShell< dim >, aspect::GeometryModel::EllipsoidalChunk< dim >, aspect::GeometryModel::TwoMergedChunks< dim >, aspect::GeometryModel::Box< dim >, aspect::GeometryModel::TwoMergedBoxes< dim >, and aspect::GeometryModel::Sphere< dim >.
|
pure virtual |
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.
Implemented in aspect::GeometryModel::Chunk< dim >, aspect::GeometryModel::SphericalShell< dim >, aspect::GeometryModel::EllipsoidalChunk< dim >, aspect::GeometryModel::TwoMergedChunks< dim >, aspect::GeometryModel::Box< dim >, aspect::GeometryModel::TwoMergedBoxes< dim >, and aspect::GeometryModel::Sphere< dim >.
|
pure virtual |
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.
Implemented in aspect::GeometryModel::SphericalShell< dim >, aspect::GeometryModel::Chunk< dim >, aspect::GeometryModel::EllipsoidalChunk< dim >, aspect::GeometryModel::Box< dim >, aspect::GeometryModel::TwoMergedBoxes< dim >, aspect::GeometryModel::TwoMergedChunks< dim >, and aspect::GeometryModel::Sphere< dim >.
|
virtual |
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.
An example would be that the "box" geometry returns a map of the form {{"left"->0}, {"right"->1}, {"bottom"->2}, {"top"->3}}
in 2d.
The default implementation of this function returns an empty map. This still allows the use of a geometry model that does not implement this function but forces the user to identify parts of the boundary by their boundary indicator number, rather than using a symbolic name.
Reimplemented in aspect::GeometryModel::SphericalShell< dim >, aspect::GeometryModel::Chunk< dim >, aspect::GeometryModel::EllipsoidalChunk< dim >, aspect::GeometryModel::Box< dim >, aspect::GeometryModel::TwoMergedBoxes< dim >, aspect::GeometryModel::TwoMergedChunks< dim >, and aspect::GeometryModel::Sphere< dim >.
types::boundary_id aspect::GeometryModel::Interface< dim >::translate_symbolic_boundary_name_to_id | ( | const std::string & | name | ) | const |
For a given name of a boundary component, translate it to its numeric value – either by using one of the symbolic values in the mapping that derived classes provide through the get_symbolic_boundary_names_map() function, or by converting its string representation into a number.
name | A name or number (as string). Leading and trailing spaces are removed from this string. |
std::vector<types::boundary_id> aspect::GeometryModel::Interface< dim >::translate_symbolic_boundary_names_to_ids | ( | const std::vector< std::string > & | names | ) | const |
For each one of the given names of boundary components, translate it to its numeric value – either by using one of the symbolic values in the mapping that derived classes provide through the get_symbolic_boundary_names_map() function, or by converting its string representation into a number.
names | A list of names or numbers (as strings). Leading and trailing spaces are removed from the strings. |
std::string aspect::GeometryModel::Interface< dim >::translate_id_to_symbol_name | ( | const types::boundary_id | boundary_id | ) | const |
Given a boundary indicator, try and see whether this geometry model provides a symbolic name for it. If the geometry model does not provide a symbolic name, return the empty string.
boundary_id | The boundary indicator to be translated. |
|
virtual |
Returns a set of periodic boundary pairs. The elements of the set are a pair of boundary ids and a cartesian unit direction each. The base class returns an empty set, so this does nothing unless you specifically use a geometry model with periodic boundary conditions
Reimplemented in aspect::GeometryModel::SphericalShell< dim >, aspect::GeometryModel::Box< dim >, and aspect::GeometryModel::TwoMergedBoxes< dim >.
|
virtual |
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.
Reimplemented in aspect::GeometryModel::SphericalShell< dim >, aspect::GeometryModel::Box< dim >, and aspect::GeometryModel::TwoMergedBoxes< dim >.
|
virtual |
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
.
Reimplemented in aspect::GeometryModel::Chunk< dim >, aspect::GeometryModel::SphericalShell< dim >, aspect::GeometryModel::TwoMergedChunks< dim >, aspect::GeometryModel::Box< dim >, aspect::GeometryModel::TwoMergedBoxes< dim >, and aspect::GeometryModel::Sphere< dim >.
|
pure virtual |
If true, the queried point (in Cartesian coordinates) lies in the domain specified by the geometry.
Implemented in aspect::GeometryModel::Chunk< dim >, aspect::GeometryModel::SphericalShell< dim >, aspect::GeometryModel::TwoMergedChunks< dim >, aspect::GeometryModel::EllipsoidalChunk< dim >, aspect::GeometryModel::Box< dim >, aspect::GeometryModel::TwoMergedBoxes< dim >, and aspect::GeometryModel::Sphere< dim >.
|
virtual |
Collects periodic boundary constraints for the given geometry and dof_handler
, which will be added to the existing constraints
. The default implementation creates cartesian periodic boundary conditions for all periodic boundary indicators.
Reimplemented in aspect::GeometryModel::SphericalShell< dim >.