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

Public Member Functions

void initialize () override
 
void create_coarse_mesh (parallel::distributed::Triangulation< dim > &coarse_grid) 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
 
bool point_is_in_domain (const Point< dim > &point) 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
 
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
 
double get_radius (const Point< dim > &position) const
 
double get_semi_minor_axis_b () const
 
double get_semi_major_axis_a () const
 
double get_eccentricity () const
 
const std::vector< Point< 2 > > & get_corners () const
 
internal::EllipsoidalChunkGeometry< dim > get_manifold () const
 
- 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 bool has_curved_elements () 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_ids (parallel::distributed::Triangulation< dim > &coarse_grid) const
 

Private Attributes

std::vector< Point< 2 > > corners
 
double semi_major_axis_a
 
double eccentricity
 
double semi_minor_axis_b
 
double rot_para_to_para_angle
 
double para_to_rect_angle
 
double rotation_longitude
 
double rotation_latitude
 
double bottom_depth
 
double westLongitude
 
double eastLongitude
 
double northLatitude
 
double southLatitude
 
unsigned int EW_subdiv
 
unsigned int NS_subdiv
 
unsigned int depth_subdiv
 
std::unique_ptr< const internal::EllipsoidalChunkGeometry< dim > > manifold
 

Detailed Description

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

A class that describes a geometry for an ellipsoid such as the WGS84 model of the earth.

This geometry model implements a (3d) ellipsoidal chunk geometry where two of the axis have the same length. The ellipsoidal chunk can be a non-coordinate parallel part of the ellipsoid.

Author
This plugin is a joined effort of Menno Fraters, D. Sarah Stamps and Wolfgang Bangerth

Definition at line 140 of file ellipsoidal_chunk.h.

Member Function Documentation

§ initialize()

template<int dim>
void aspect::GeometryModel::EllipsoidalChunk< dim >::initialize ( )
overridevirtual

Initialize function

Reimplemented from aspect::Plugins::InterfaceBase.

§ create_coarse_mesh()

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

§ length_scale()

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

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

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

§ depth()

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

Placeholder for a function returning the height of the given position relative to the reference model surface.

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

§ representative_point()

template<int dim>
Point<dim> aspect::GeometryModel::EllipsoidalChunk< dim >::representative_point ( const double  depth) const
overridevirtual

Returns a point in the center of the domain.

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

§ point_is_in_domain()

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

§ maximal_depth()

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

Returns the bottom depth which was used to create the geometry and which is defined by the depth parameter.

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

§ get_used_boundary_indicators()

template<int dim>
std::set<types::boundary_id> aspect::GeometryModel::EllipsoidalChunk< 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::EllipsoidalChunk< dim >::get_symbolic_boundary_names_map ( ) const
overridevirtual

Set symbolic names for boundaries (mrtf)

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

§ natural_coordinate_system()

template<int dim>
aspect::Utilities::Coordinates::CoordinateSystem aspect::GeometryModel::EllipsoidalChunk< dim >::natural_coordinate_system ( ) const
overridevirtual

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

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

§ cartesian_to_natural_coordinates()

template<int dim>
std::array<double,dim> aspect::GeometryModel::EllipsoidalChunk< 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 ellipsoidal chunk this is (radius, longitude) in 2d and (radius, longitude, latitude) in 3d. Note that internally the coordinates are stored in longitude, latitude, depth.

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

§ natural_to_cartesian_coordinates()

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

Declare the parameters this class takes through input files.

§ parse_parameters()

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

Read the parameters this class declares from the parameter file.

Reimplemented from aspect::Plugins::InterfaceBase.

§ get_radius()

template<int dim>
double aspect::GeometryModel::EllipsoidalChunk< dim >::get_radius ( const Point< dim > &  position) const

Calculate radius at current position.

§ get_semi_minor_axis_b()

template<int dim>
double aspect::GeometryModel::EllipsoidalChunk< dim >::get_semi_minor_axis_b ( ) const

Retrieve the semi minor axis b value.

§ get_semi_major_axis_a()

template<int dim>
double aspect::GeometryModel::EllipsoidalChunk< dim >::get_semi_major_axis_a ( ) const

Retrieve the semi major axis a value.

§ get_eccentricity()

template<int dim>
double aspect::GeometryModel::EllipsoidalChunk< dim >::get_eccentricity ( ) const

Retrieve the value of the eccentricity.

§ get_corners()

template<int dim>
const std::vector<Point<2> >& aspect::GeometryModel::EllipsoidalChunk< dim >::get_corners ( ) const

Retrieve the corners used to create the ellipsoid. This variable contains four vectors with each two elements. Each set of two elements represents a longitude and latitude value. The four vectors represent respectively the point in the North-East, North-West, South-West and South-East.

§ get_manifold()

template<int dim>
internal::EllipsoidalChunkGeometry<dim> aspect::GeometryModel::EllipsoidalChunk< dim >::get_manifold ( ) const

Retrieve the manifold object.

§ set_boundary_ids()

template<int dim>
void aspect::GeometryModel::EllipsoidalChunk< dim >::set_boundary_ids ( parallel::distributed::Triangulation< dim > &  coarse_grid) const
private

Member Data Documentation

§ corners

template<int dim>
std::vector<Point<2> > aspect::GeometryModel::EllipsoidalChunk< dim >::corners
private

Definition at line 306 of file ellipsoidal_chunk.h.

§ semi_major_axis_a

template<int dim>
double aspect::GeometryModel::EllipsoidalChunk< dim >::semi_major_axis_a
private

Definition at line 307 of file ellipsoidal_chunk.h.

§ eccentricity

template<int dim>
double aspect::GeometryModel::EllipsoidalChunk< dim >::eccentricity
private

Definition at line 308 of file ellipsoidal_chunk.h.

§ semi_minor_axis_b

template<int dim>
double aspect::GeometryModel::EllipsoidalChunk< dim >::semi_minor_axis_b
private

Definition at line 309 of file ellipsoidal_chunk.h.

§ rot_para_to_para_angle

template<int dim>
double aspect::GeometryModel::EllipsoidalChunk< dim >::rot_para_to_para_angle
private

Definition at line 310 of file ellipsoidal_chunk.h.

§ para_to_rect_angle

template<int dim>
double aspect::GeometryModel::EllipsoidalChunk< dim >::para_to_rect_angle
private

Definition at line 311 of file ellipsoidal_chunk.h.

§ rotation_longitude

template<int dim>
double aspect::GeometryModel::EllipsoidalChunk< dim >::rotation_longitude
private

Definition at line 312 of file ellipsoidal_chunk.h.

§ rotation_latitude

template<int dim>
double aspect::GeometryModel::EllipsoidalChunk< dim >::rotation_latitude
private

Definition at line 313 of file ellipsoidal_chunk.h.

§ bottom_depth

template<int dim>
double aspect::GeometryModel::EllipsoidalChunk< dim >::bottom_depth
private

Definition at line 314 of file ellipsoidal_chunk.h.

§ westLongitude

template<int dim>
double aspect::GeometryModel::EllipsoidalChunk< dim >::westLongitude
private

Definition at line 315 of file ellipsoidal_chunk.h.

§ eastLongitude

template<int dim>
double aspect::GeometryModel::EllipsoidalChunk< dim >::eastLongitude
private

Definition at line 316 of file ellipsoidal_chunk.h.

§ northLatitude

template<int dim>
double aspect::GeometryModel::EllipsoidalChunk< dim >::northLatitude
private

Definition at line 317 of file ellipsoidal_chunk.h.

§ southLatitude

template<int dim>
double aspect::GeometryModel::EllipsoidalChunk< dim >::southLatitude
private

Definition at line 318 of file ellipsoidal_chunk.h.

§ EW_subdiv

template<int dim>
unsigned int aspect::GeometryModel::EllipsoidalChunk< dim >::EW_subdiv
private

Definition at line 320 of file ellipsoidal_chunk.h.

§ NS_subdiv

template<int dim>
unsigned int aspect::GeometryModel::EllipsoidalChunk< dim >::NS_subdiv
private

Definition at line 321 of file ellipsoidal_chunk.h.

§ depth_subdiv

template<int dim>
unsigned int aspect::GeometryModel::EllipsoidalChunk< dim >::depth_subdiv
private

Definition at line 322 of file ellipsoidal_chunk.h.

§ manifold

template<int dim>
std::unique_ptr<const internal::EllipsoidalChunkGeometry<dim> > aspect::GeometryModel::EllipsoidalChunk< 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 337 of file ellipsoidal_chunk.h.


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