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

Public Member Functions

std::pair< std::string, std::string > execute (TableHandler &statistics) override
 
const LinearAlgebra::BlockVectortopography_vector () const
 
const Vector< float > & cellwise_topography () const
 
std::list< std::string > required_other_postprocessors () const override
 
void parse_parameters (ParameterHandler &prm) override
 
- Public Member Functions inherited from aspect::Postprocess::Interface< dim >
virtual ~Interface ()=default
 
virtual void initialize ()
 
virtual void update ()
 
virtual void save (std::map< std::string, std::string > &status_strings) const
 
virtual void load (const std::map< std::string, std::string > &status_strings)
 
- 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::Postprocess::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 output_to_file (const types::boundary_id boundary_id, const std::vector< std::pair< Point< dim >, double >> &position_and_topography)
 

Private Attributes

LinearAlgebra::BlockVector topo_vector
 
Vector< float > visualization_values
 
double density_above
 
double density_below
 
bool output_surface
 
bool output_bottom
 

Detailed Description

template<int dim>
class aspect::Postprocess::DynamicTopography< dim >

A postprocessor that computes dynamic topography at the top and bottom of the domain.

Definition at line 40 of file dynamic_topography.h.

Member Function Documentation

§ execute()

template<int dim>
std::pair<std::string,std::string> aspect::Postprocess::DynamicTopography< dim >::execute ( TableHandler &  statistics)
overridevirtual

Evaluate the solution for the dynamic topography.

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

§ topography_vector()

template<int dim>
const LinearAlgebra::BlockVector& aspect::Postprocess::DynamicTopography< dim >::topography_vector ( ) const

Return the topography vector as calculated by the consistent boundary flux (CBF) formulation. The velocity components store the surface-normal traction, and the temperature component stores the dynamic topography computed from that traction.

§ cellwise_topography()

template<int dim>
const Vector<float>& aspect::Postprocess::DynamicTopography< dim >::cellwise_topography ( ) const

Return the cell-wise topography vector as calculated by the CBF formulation, where indices of the vector correspond to cell indices. This vector is considerably smaller than the full topography vector returned by topography_vector(), and is useful for text output and visualization.

§ required_other_postprocessors()

template<int dim>
std::list<std::string> aspect::Postprocess::DynamicTopography< dim >::required_other_postprocessors ( ) const
overridevirtual

Register the other postprocessor that we need: BoundaryPressures

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

§ parse_parameters()

template<int dim>
void aspect::Postprocess::DynamicTopography< dim >::parse_parameters ( ParameterHandler &  prm)
overridevirtual

Parse the parameters for the postprocessor.

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

§ declare_parameters()

template<int dim>
static void aspect::Postprocess::DynamicTopography< dim >::declare_parameters ( ParameterHandler &  prm)
static

Declare the parameters for the postprocessor.

§ output_to_file()

template<int dim>
void aspect::Postprocess::DynamicTopography< dim >::output_to_file ( const types::boundary_id  boundary_id,
const std::vector< std::pair< Point< dim >, double >> &  position_and_topography 
)
private

Output the dynamic topography solution to a file for a given boundary id. All the values in position_and_topography are written to a file with a file name determined by boundary_id.

Member Data Documentation

§ topo_vector

template<int dim>
LinearAlgebra::BlockVector aspect::Postprocess::DynamicTopography< dim >::topo_vector
private

A vector which stores the surface stress values calculated by the postprocessor.

Definition at line 101 of file dynamic_topography.h.

§ visualization_values

template<int dim>
Vector<float> aspect::Postprocess::DynamicTopography< dim >::visualization_values
private

A vector which stores the surface stress values calculated at the midpoint of each surface cell face. This can be given to a visualization postprocessor to output dynamic topography.

Definition at line 108 of file dynamic_topography.h.

§ density_above

template<int dim>
double aspect::Postprocess::DynamicTopography< dim >::density_above
private

A parameter that allows users to set the density value above the top surface.

Definition at line 114 of file dynamic_topography.h.

§ density_below

template<int dim>
double aspect::Postprocess::DynamicTopography< dim >::density_below
private

A parameter that allows users to set the density value below the bottom surface.

Definition at line 120 of file dynamic_topography.h.

§ output_surface

template<int dim>
bool aspect::Postprocess::DynamicTopography< dim >::output_surface
private

Whether to output the surface topography.

Definition at line 125 of file dynamic_topography.h.

§ output_bottom

template<int dim>
bool aspect::Postprocess::DynamicTopography< dim >::output_bottom
private

Whether to output the bottom topography.

Definition at line 130 of file dynamic_topography.h.


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