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

Public Member Functions

void initialize () override
 
void topography (typename parallel::distributed::Triangulation< dim > &grid) const
 
Point< dim > add_topography (const Point< dim > &x_y_z) const
 
void create_coarse_mesh (parallel::distributed::Triangulation< dim > &coarse_grid) const override
 
virtual Point< dim > get_extents () const
 
const std::array< unsigned int, dim > & get_repetitions () const
 
virtual Point< dim > get_origin () const
 
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
 
- Public Member Functions inherited from aspect::GeometryModel::Interface< dim >
virtual ~Interface ()=default
 
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 void make_periodicity_constraints (const DoFHandler< dim > &dof_handler, AffineConstraints< double > &constraints) const
 
- 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
 
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 Public Member Functions inherited from aspect::GeometryModel::Interface< dim >
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 Attributes

InitialTopographyModel::Interface< dim > * topo_model
 
Point< dim > extents
 
Point< dim > box_origin
 
std::array< bool, dim > periodic
 
std::array< unsigned int, dim > repetitions
 

Detailed Description

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

A class that describes a box geometry of certain width, height, and depth (in 3d), and, possibly, topography.

Definition at line 40 of file box.h.

Member Function Documentation

§ initialize()

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

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

§ topography()

template<int dim>
void aspect::GeometryModel::Box< dim >::topography ( typename parallel::distributed::Triangulation< dim > &  grid) const

Add initial topography to the mesh.

§ add_topography()

template<int dim>
Point<dim> aspect::GeometryModel::Box< dim >::add_topography ( const Point< dim > &  x_y_z) const

Relocate the vertical coordinate of the given point based on the topography at the surface specified by the initial topography model.

§ create_coarse_mesh()

template<int dim>
void aspect::GeometryModel::Box< 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_extents()

template<int dim>
virtual Point<dim> aspect::GeometryModel::Box< dim >::get_extents ( ) const
virtual

Return a point that denotes the size of the box in each dimension of the domain.

§ get_repetitions()

template<int dim>
const std::array<unsigned int, dim>& aspect::GeometryModel::Box< dim >::get_repetitions ( ) const

Return an integer array that denotes the number of repetitions of the box's coarse mesh.

§ get_origin()

template<int dim>
virtual Point<dim> aspect::GeometryModel::Box< dim >::get_origin ( ) const
virtual

Return a point that denotes the lower left corner of the box domain.

§ length_scale()

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

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

We return 1/100th of the diameter of the box.

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

§ depth()

template<int dim>
double aspect::GeometryModel::Box< 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 \((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. In almost all cases one will use a gravity model that also matches these definitions.

Note that the depth is calculated with respect to the surface without initial topography.

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

§ height_above_reference_surface()

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

Return the height of the given position relative to the initial box height.

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

§ representative_point()

template<int dim>
Point<dim> aspect::GeometryModel::Box< 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::Box< 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::Box< 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 box model uses boundary indicators zero through 2*dim-1, with the first two being the faces perpendicular to the x-axis, the next two perpendicular to the y-axis, etc.

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

§ get_symbolic_boundary_names_map()

template<int dim>
std::map<std::string,types::boundary_id> aspect::GeometryModel::Box< 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.

This geometry returns the map {{"left"->0}, {"right"->1}, {"bottom"->2}, {"top"->3}} in 2d, and {{"left"->0}, {"right"->1}, {"front"->2}, {"back"->3}, {"bottom"->4}, {"top"->5}} in 3d.

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::Box< 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::Box< 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 translation 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::Box< 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 box has only straight boundaries and cells, so return false.

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

§ point_is_in_domain()

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

Returns what the natural coordinate system for this geometry model is, which for a box is Cartesian.

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

§ cartesian_to_natural_coordinates()

template<int dim>
std::array<double,dim> aspect::GeometryModel::Box< 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 box the results is unchanged and is (x,z) in 2d or (x,y,z) in 3d.

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

§ natural_to_cartesian_coordinates()

template<int dim>
Point<dim> aspect::GeometryModel::Box< 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::Box< dim >::declare_parameters ( ParameterHandler &  prm)
static

Declare the parameters this class takes through input files.

§ parse_parameters()

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

Read the parameters this class declares from the parameter file.

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

Member Data Documentation

§ topo_model

template<int dim>
InitialTopographyModel::Interface<dim>* aspect::GeometryModel::Box< dim >::topo_model
private

A pointer to the initial topography model.

Definition at line 229 of file box.h.

§ extents

template<int dim>
Point<dim> aspect::GeometryModel::Box< dim >::extents
private

Extent of the box in x-, y-, and z-direction (in 3d).

Definition at line 234 of file box.h.

§ box_origin

template<int dim>
Point<dim> aspect::GeometryModel::Box< dim >::box_origin
private

Origin of the box in x, y, and z (in 3d) coordinates.

Definition at line 239 of file box.h.

§ periodic

template<int dim>
std::array<bool, dim> aspect::GeometryModel::Box< dim >::periodic
private

Flag whether the box is periodic in the x-, y-, and z-direction.

Definition at line 244 of file box.h.

§ repetitions

template<int dim>
std::array<unsigned int, dim> aspect::GeometryModel::Box< dim >::repetitions
private

The number of cells in each coordinate direction.

Definition at line 249 of file box.h.


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