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

Public Member Functions

virtual ~Interface ()=default
 
virtual void initialize ()
 
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 parse_parameters (ParameterHandler &prm)
 
virtual void make_periodicity_constraints (const DoFHandler< dim > &dof_handler, AffineConstraints< double > &constraints) const
 

Static Public Member Functions

static void declare_parameters (ParameterHandler &prm)
 

Detailed Description

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

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.

Constructor & Destructor Documentation

§ ~Interface()

template<int dim>
virtual aspect::GeometryModel::Interface< dim >::~Interface ( )
virtualdefault

Destructor. Made virtual to enforce that derived classes also have virtual destructors.

Member Function Documentation

§ initialize()

template<int dim>
virtual void aspect::GeometryModel::Interface< dim >::initialize ( )
virtual

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.

Reimplemented in aspect::GeometryModel::Chunk< dim >, aspect::GeometryModel::SphericalShell< dim >, aspect::GeometryModel::EllipsoidalChunk< dim >, aspect::GeometryModel::TwoMergedChunks< dim >, and aspect::GeometryModel::Box< dim >.

§ create_coarse_mesh()

template<int dim>
virtual void aspect::GeometryModel::Interface< dim >::create_coarse_mesh ( parallel::distributed::Triangulation< dim > &  coarse_grid) const
pure virtual

§ length_scale()

template<int dim>
virtual double aspect::GeometryModel::Interface< dim >::length_scale ( ) const
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 >.

§ depth()

template<int dim>
virtual double aspect::GeometryModel::Interface< dim >::depth ( const Point< dim > &  position) const
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.

Note
Implementations of this function in derived classes can only compute the depth with regard to the reference configuration of the geometry, i.e., the geometry initially created. If you are using a dynamic topography in your models that changes in every time step, or if you apply initial topography to your model, then the actual depth of a point with regard to this dynamic topography will not match the value this function returns. This is so because computing the actual depth is difficult: In parallel computations, the processor on which you want to evaluate the depth of a point may know nothing about the displacement of the surface anywhere if it happens to store only interior cells. Furthermore, it is not even clear what "depth" one would compute in such situations: The distance to the closest surface point? The vertical distance to the surface point directly above? Or the length of the line from the given point to a surface point that is locally always parallel to the gravity vector? For all of these reasons, this function simply returns the geometric vertical depth of a point in the known and fixed reference geometry, to the reference surface that defines what zero depth 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 >.

§ height_above_reference_surface()

template<int dim>
virtual double aspect::GeometryModel::Interface< dim >::height_above_reference_surface ( const Point< dim > &  position) const
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 >.

§ cartesian_to_other_coordinates()

template<int 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.

§ natural_coordinate_system()

template<int dim>
virtual aspect::Utilities::Coordinates::CoordinateSystem aspect::GeometryModel::Interface< dim >::natural_coordinate_system ( ) const
pure virtual

§ cartesian_to_natural_coordinates()

template<int dim>
virtual std::array<double,dim> aspect::GeometryModel::Interface< dim >::cartesian_to_natural_coordinates ( const Point< dim > &  position) const
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 >.

§ natural_to_cartesian_coordinates()

template<int dim>
virtual Point<dim> aspect::GeometryModel::Interface< dim >::natural_to_cartesian_coordinates ( const std::array< double, dim > &  position) const
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 >.

§ representative_point()

template<int dim>
virtual Point<dim> aspect::GeometryModel::Interface< dim >::representative_point ( const double  depth) const
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 >.

§ maximal_depth()

template<int dim>
virtual double aspect::GeometryModel::Interface< dim >::maximal_depth ( ) const
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 >.

§ get_used_boundary_indicators()

template<int dim>
virtual std::set<types::boundary_id> aspect::GeometryModel::Interface< dim >::get_used_boundary_indicators ( ) const
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 >.

§ get_symbolic_boundary_names_map()

template<int dim>
virtual std::map<std::string,types::boundary_id> aspect::GeometryModel::Interface< dim >::get_symbolic_boundary_names_map ( ) const
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.

Note
Names may contain spaces and numbers, but they may not contain special characters and they should not equal the text representation of numbers (e.g., a name "10" is ill-advised).
Since in practice boundary indicators can be provided either via number or symbolic name, the mapping from something given in the input is not entirely trivial – in particular, because a function also has to do some error checking that a given string in fact matches any known boundary indicator. To this end, use GeometryModel::translate_boundary_indicator() and GeometryModel::translate_boundary_indicators().
Returns
A map from symbolic names to boundary indicators. The map should provide a symbolic name for each used boundary indicator as returned by get_user_boundary_indicators(). In the end, however, a geometry model may define multiple symbolic names for the same boundary or not define any.

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 >.

§ translate_symbolic_boundary_name_to_id()

template<int 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.

Parameters
nameA name or number (as string). Leading and trailing spaces are removed from this string.
Returns
A boundary indicator number corresponding to the given name. If the name does not represent either a symbolic name or a number, this function will throw an exception of type std::string that explains the error.

§ translate_symbolic_boundary_names_to_ids()

template<int dim>
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.

Parameters
namesA list of names or numbers (as strings). Leading and trailing spaces are removed from the strings.
Returns
A list of boundary indicator numbers corresponding to the given list of names. If one of the given names does not represent either a symbolic name or a number, this function will throw an exception of type std::string that explains the error.

§ translate_id_to_symbol_name()

template<int dim>
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.

Parameters
boundary_idThe boundary indicator to be translated.
Returns
A string representation for this boundary indicator, if one is available.

§ get_periodic_boundary_pairs()

template<int dim>
virtual std::set<std::pair<std::pair<types::boundary_id, types::boundary_id>, unsigned int> > aspect::GeometryModel::Interface< dim >::get_periodic_boundary_pairs ( ) const
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 >.

§ adjust_positions_for_periodicity()

template<int dim>
virtual void aspect::GeometryModel::Interface< dim >::adjust_positions_for_periodicity ( Point< dim > &  position,
const ArrayView< Point< dim >> &  connected_positions = {},
const ArrayView< Tensor< 1, dim >> &  connected_velocities = {} 
) const
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 >.

§ has_curved_elements()

template<int dim>
virtual bool aspect::GeometryModel::Interface< dim >::has_curved_elements ( ) const
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 >.

§ point_is_in_domain()

template<int dim>
virtual bool aspect::GeometryModel::Interface< dim >::point_is_in_domain ( const Point< dim > &  p) const
pure virtual

§ declare_parameters()

template<int dim>
static void aspect::GeometryModel::Interface< 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>
virtual void aspect::GeometryModel::Interface< dim >::parse_parameters ( ParameterHandler &  prm)
virtual

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 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 >.

§ make_periodicity_constraints()

template<int dim>
virtual void aspect::GeometryModel::Interface< dim >::make_periodicity_constraints ( const DoFHandler< dim > &  dof_handler,
AffineConstraints< double > &  constraints 
) const
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 >.


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