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

Classes

class  ChunkGeometry
 

Public Member Functions

void initialize ()
 
virtual void create_coarse_mesh (parallel::distributed::Triangulation< dim > &coarse_grid) const
 
void connect_to_signal (SimulatorSignals< dim > &signals)
 
virtual std::set< types::boundary_idget_used_boundary_indicators () const
 
virtual std::map< std::string, types::boundary_idget_symbolic_boundary_names_map () const
 
virtual double length_scale () const
 
virtual double depth (const Point< dim > &position) const
 
virtual double height_above_reference_surface (const Point< dim > &position) const
 
virtual Point< dim > representative_point (const double depth) const
 
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
 
virtual double maximal_depth () const
 
virtual double inner_radius () const
 
virtual double outer_radius () const
 
virtual bool has_curved_elements () const
 
virtual bool point_is_in_domain (const Point< dim > &p) const
 
virtual aspect::Utilities::Coordinates::CoordinateSystem natural_coordinate_system () const
 
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 void parse_parameters (ParameterHandler &prm)
 
- Public Member Functions inherited from aspect::GeometryModel::Interface< dim >
virtual ~Interface ()
 
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_idtranslate_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
 
- Public Member Functions inherited from aspect::SimulatorAccess< dim >
 SimulatorAccess ()
 
 SimulatorAccess (const Simulator< dim > &simulator_object)
 
virtual ~SimulatorAccess ()
 
virtual void initialize_simulator (const Simulator< dim > &simulator_object)
 
template<typename PostprocessorType >
PostprocessorType * find_postprocessor () const
 
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
 
TimerOutputget_computing_timer () const
 
const ConditionalOStreamget_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::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
 
void compute_material_model_input_values (const LinearAlgebra::BlockVector &input_solution, const FEValuesBase< dim, dim > &input_finite_element_values, const typename DoFHandler< dim >::active_cell_iterator &cell, const bool compute_strainrate, MaterialModel::MaterialModelInputs< dim > &material_model_inputs) 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
 
DEAL_II_DEPRECATED const BoundaryTemperature::Interface< dim > & get_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
 
DEAL_II_DEPRECATED const BoundaryComposition::Interface< dim > & get_boundary_composition () const
 
const BoundaryComposition::Manager< dim > & get_boundary_composition_manager () const
 
const std::map< types::boundary_id, std::shared_ptr< BoundaryTraction::Interface< dim > > > & get_boundary_traction () const
 
DEAL_II_DEPRECATED const InitialTemperature::Interface< dim > & get_initial_temperature () const
 
const InitialTemperature::Manager< dim > & get_initial_temperature_manager () const
 
DEAL_II_DEPRECATED const InitialComposition::Interface< dim > & get_initial_composition () 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_free_surface_boundary_indicators () const
 
DEAL_II_DEPRECATED const std::map< types::boundary_id, std::shared_ptr< BoundaryVelocity::Interface< dim > > > get_prescribed_boundary_velocity () 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 NewtonHandler< dim > & get_newton_handler () const
 
const WorldBuilder::World & get_world_builder () const
 
const FreeSurfaceHandler< dim > & get_free_surface_handler () const
 
const LateralAveraging< dim > & get_lateral_averaging () const
 
const ConstraintMatrix & get_current_constraints () const
 
bool simulator_is_initialized () const
 
double get_pressure_scaling () const
 
bool pressure_rhs_needs_compatibility_modification () const
 
bool model_has_prescribed_stokes_solution () const
 
TableHandlerget_statistics_object () const
 
template<typename PostprocessorType >
DEAL_II_DEPRECATED PostprocessorType * find_postprocessor () const
 
const Postprocess::Manager< dim > & get_postprocess_manager () 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 Member Functions

void set_manifold_ids (typename parallel::distributed::Triangulation< dim > &triangulation)
 
void clear_manifold_ids (typename parallel::distributed::Triangulation< dim > &triangulation)
 

Private Attributes

Point< dim > point1
 
Point< dim > point2
 
unsigned int repetitions [dim]
 
ChunkGeometry manifold
 

Detailed Description

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

A geometry model class that describes a chunk of a spherical shell. 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 west, east, south, north, inner and outer.

The parameters that describe this geometry and that are read from the input file are the inner and outer radii of the shell, the minimum and maximum longitude, minimum and maximum longitude, and the number of cells initialised in each dimension.

Definition at line 54 of file chunk.h.

Member Function Documentation

§ initialize()

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

§ create_coarse_mesh()

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

Generate a coarse mesh for the geometry described by this class.

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

§ connect_to_signal()

template<int dim>
void aspect::GeometryModel::Chunk< dim >::connect_to_signal ( SimulatorSignals< dim > &  signals)

Connect the set/clear manifold id functions to the post/pre_compute_no_normal_flux signals.

§ get_used_boundary_indicators()

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

The chunk model uses boundary indicators zero through 2*dim-1, with the first two being the faces perpendicular to the radius of the shell, the next two along lines of longitude and, in 3d, the next two along lines of latitude.

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

§ get_symbolic_boundary_names_map()

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

This geometry returns the map {{"bottom"->0}, {"top"->1}, {"west"->2}, {"east"->3}} in 2d, and {{"bottom"->0}, {"top"->1}, {"west"->2}, {"east"->3}, {"south"->4}, {"north"->5}} in 3d.

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

§ length_scale()

template<int dim>
virtual double aspect::GeometryModel::Chunk< dim >::length_scale ( ) const
virtual

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

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

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

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

§ representative_point()

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

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

§ west_longitude()

template<int dim>
virtual double aspect::GeometryModel::Chunk< 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::Chunk< 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::Chunk< dim >::longitude_range ( ) const
virtual

Return the longitude range of the chunk measured in radians.

§ south_latitude()

template<int dim>
virtual double aspect::GeometryModel::Chunk< 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::Chunk< 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::Chunk< dim >::latitude_range ( ) const
virtual

Return the latitude range of the chunk Measured in radians

§ maximal_depth()

template<int dim>
virtual double aspect::GeometryModel::Chunk< dim >::maximal_depth ( ) const
virtual

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::Chunk< dim >::inner_radius ( ) const
virtual

Return the inner radius of the chunk measured in meters.

§ outer_radius()

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

Return the outer radius of the chunk measured in meters.

§ has_curved_elements()

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

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

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

§ point_is_in_domain()

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

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>
virtual aspect::Utilities::Coordinates::CoordinateSystem aspect::GeometryModel::Chunk< dim >::natural_coordinate_system ( ) const
virtual

Returns what the natural coordinate system for this geometry model is.

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

§ cartesian_to_natural_coordinates()

template<int dim>
virtual std::array<double,dim> aspect::GeometryModel::Chunk< 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 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>
virtual Point<dim> aspect::GeometryModel::Chunk< 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 from aspect::GeometryModel::Interface< dim >.

§ declare_parameters()

template<int dim>
static void aspect::GeometryModel::Chunk< dim >::declare_parameters ( ParameterHandler prm)
static

Declare the parameters this class takes through input files.

§ parse_parameters()

template<int dim>
virtual void aspect::GeometryModel::Chunk< dim >::parse_parameters ( ParameterHandler prm)
virtual

Read the parameters this class declares from the parameter file.

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

§ set_manifold_ids()

template<int dim>
void aspect::GeometryModel::Chunk< dim >::set_manifold_ids ( typename parallel::distributed::Triangulation< dim > &  triangulation)
private

§ clear_manifold_ids()

template<int dim>
void aspect::GeometryModel::Chunk< dim >::clear_manifold_ids ( typename parallel::distributed::Triangulation< dim > &  triangulation)
private

Member Data Documentation

§ point1

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

Minimum depth, longitude-depth or longitude-latitude-depth point

Definition at line 266 of file chunk.h.

§ point2

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

Maximum depth, longitude-depth or longitude-latitude-depth point

Definition at line 272 of file chunk.h.

§ repetitions

template<int dim>
unsigned int aspect::GeometryModel::Chunk< dim >::repetitions[dim]
private

The number of cells in each coordinate direction

Definition at line 277 of file chunk.h.

§ manifold

template<int dim>
ChunkGeometry aspect::GeometryModel::Chunk< dim >::manifold
private

An object that describes the geometry.

Definition at line 339 of file chunk.h.


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