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

Public Member Functions

std::vector< std::vector< double > > get_averages (const unsigned int n_slices, const std::vector< std::string > &property_names) const
 
void get_temperature_averages (std::vector< double > &values) const
 
void get_composition_averages (const unsigned int composition_index, std::vector< double > &values) const
 
void get_viscosity_averages (std::vector< double > &values) const
 
void get_velocity_magnitude_averages (std::vector< double > &values) const
 
void get_sinking_velocity_averages (std::vector< double > &values) const
 
void get_Vs_averages (std::vector< double > &values) const
 
void get_Vp_averages (std::vector< double > &values) const
 
void get_vertical_heat_flux_averages (std::vector< double > &values) 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::unique_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_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 WorldBuilder::World & get_world_builder () const
 
const MeshDeformation::MeshDeformationHandler< dim > & get_mesh_deformation_handler () const
 
const LateralAveraging< dim > & get_lateral_averaging () const
 
const ConstraintMatrix & 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
 
TableHandlerget_statistics_object () const
 
template<typename PostprocessorType >
DEAL_II_DEPRECATED PostprocessorType * find_postprocessor () const
 
const Postprocess::Manager< dim > & get_postprocess_manager () const
 

Private Member Functions

std::vector< std::vector< double > > compute_lateral_averages (const unsigned int n_slices, std::vector< std::unique_ptr< internal::FunctorBase< dim > > > &functors) const
 

Additional Inherited Members

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

Detailed Description

template<int dim>
class aspect::LateralAveraging< dim >

LateralAveraging is a class that performs various averaging operations on the solution. The functions of this class take a vector as an argument. The model is divided into depth slices where the number of slices is the length of the output vector. Each function averages a specific quantity (as specified by their name), and that quantity is averaged laterally for each depth slice.

Plugins may access the LateralAveraging plugin through the SimulatorAccess function get_lateral_averaging(), and then query that for the desired averaged quantity.

Definition at line 110 of file lateral_averaging.h.

Member Function Documentation

§ get_averages()

template<int dim>
std::vector<std::vector<double> > aspect::LateralAveraging< dim >::get_averages ( const unsigned int  n_slices,
const std::vector< std::string > &  property_names 
) const

Fill the values with a set of lateral averages of the selected property_names. See the implementation of this function for a range of accepted names. This function is more efficient than calling multiple of the other functions that compute one property each.

Parameters
n_slicesThe number of depth slices to perform the averaging in.
Returns
The output vector of laterally averaged values. Each vector has the same size of n_slices, and there are as many vectors returned as names in property_names.
Parameters
property_namesNames of the available quantities to average. Check the implementation of this function for available names.

§ get_temperature_averages()

template<int dim>
void aspect::LateralAveraging< dim >::get_temperature_averages ( std::vector< double > &  values) const

Fill the argument with a set of lateral averages of the current temperature field. The function fills a vector that contains average field values over slices of the domain of same depth.

Parameters
valuesThe output vector of laterally averaged values. The function takes the pre-existing size of this vector as the number of depth slices.

§ get_composition_averages()

template<int dim>
void aspect::LateralAveraging< dim >::get_composition_averages ( const unsigned int  composition_index,
std::vector< double > &  values 
) const

Fill the argument with a set of lateral averages of the current compositional fields.

Parameters
composition_indexThe index of the compositional field whose matrix we want to assemble (0 <= composition_index < number of compositional fields in this problem).
valuesThe output vector of laterally averaged values. The function takes the pre-existing size of this vector as the number of depth slices.

§ get_viscosity_averages()

template<int dim>
void aspect::LateralAveraging< dim >::get_viscosity_averages ( std::vector< double > &  values) const

Compute a lateral average of the current viscosity.

Parameters
valuesThe output vector of laterally averaged values. The function takes the pre-existing size of this vector as the number of depth slices.

§ get_velocity_magnitude_averages()

template<int dim>
void aspect::LateralAveraging< dim >::get_velocity_magnitude_averages ( std::vector< double > &  values) const

Compute a lateral average of the current velocity magnitude.

Parameters
valuesThe output vector of laterally averaged values. The function takes the pre-existing size of this vector as the number of depth slices.

§ get_sinking_velocity_averages()

template<int dim>
void aspect::LateralAveraging< dim >::get_sinking_velocity_averages ( std::vector< double > &  values) const

Compute a lateral average of the current sinking velocity.

Parameters
valuesThe output vector of laterally averaged values. The function takes the pre-existing size of this vector as the number of depth slices.

§ get_Vs_averages()

template<int dim>
void aspect::LateralAveraging< dim >::get_Vs_averages ( std::vector< double > &  values) const

Compute a lateral average of the seismic shear wave speed: Vs.

Parameters
valuesThe output vector of laterally averaged values. The function takes the pre-existing size of this vector as the number of depth slices.

§ get_Vp_averages()

template<int dim>
void aspect::LateralAveraging< dim >::get_Vp_averages ( std::vector< double > &  values) const

Compute a lateral average of the seismic pressure wave speed: Vp.

Parameters
valuesThe output vector of laterally averaged values. The function takes the pre-existing size of this vector as the number of depth slices.

§ get_vertical_heat_flux_averages()

template<int dim>
void aspect::LateralAveraging< dim >::get_vertical_heat_flux_averages ( std::vector< double > &  values) const

Compute a lateral average of the heat flux, with the sign convention of positive heat flux when it flows upwards.

Parameters
valuesThe output vector of laterally averaged values. The function takes the pre-existing size of this vector as the number of depth slices.

§ compute_lateral_averages()

template<int dim>
std::vector<std::vector<double> > aspect::LateralAveraging< dim >::compute_lateral_averages ( const unsigned int  n_slices,
std::vector< std::unique_ptr< internal::FunctorBase< dim > > > &  functors 
) const
private

Internal routine to compute the depth averages of several quantities. All of the public functions that compute a single field also call this function. The vector of functors functors must contain one or more objects of classes that are derived from FunctorBase and are used to fill the values vectors.

Parameters
n_slicesNumber of depth slices to be computed.
functorsInstances of a class derived from FunctorBase that are used to compute the averaged properties.
Returns
The output vectors of depth averaged values. The function returns one vector of doubles per property and uses n_slices as the number of depth slices. Each returned vector has the same size.

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