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

Public Member Functions

 GravityPointValues ()
 
void initialize () override
 
std::pair< std::string, std::string > execute (TableHandler &) override
 
void parse_parameters (ParameterHandler &prm) override
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 
void save (std::map< std::string, std::string > &status_strings) const override
 
void load (const std::map< std::string, std::string > &status_strings) override
 
- Public Member Functions inherited from aspect::Postprocess::Interface< dim >
virtual std::list< std::string > required_other_postprocessors () 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_worlds () const
 
const Particle::World< dim > & get_particle_world (const unsigned int particle_world_index) const
 
Particle::World< dim > & get_particle_world (const unsigned int particle_world_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 Types

enum  SamplingScheme { map, fibonacci_spiral, list_of_points }
 

Private Member Functions

void set_last_output_time (const double current_time)
 

Private Attributes

double output_interval
 
double last_output_time
 
unsigned int maximum_timesteps_between_outputs
 
unsigned int last_output_timestep
 
unsigned int output_file_number
 
double end_time
 
unsigned int precision
 
unsigned int quadrature_degree_increase
 
unsigned int n_points_spiral
 
unsigned int n_points_radius
 
unsigned int n_points_longitude
 
unsigned int n_points_latitude
 
double minimum_radius
 
double maximum_radius
 
double minimum_colongitude
 
double maximum_colongitude
 
double minimum_colatitude
 
double maximum_colatitude
 
double reference_density
 
enum aspect::Postprocess::GravityPointValues::SamplingScheme sampling_scheme
 
std::vector< double > radius_list
 
std::vector< double > longitude_list
 
std::vector< double > latitude_list
 
double model_outer_radius
 
double model_inner_radius
 
std::vector< std::array< double, dim > > satellite_positions_spherical
 
std::vector< Point< dim > > satellite_positions_cartesian
 

Detailed Description

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

A postprocessor that computes gravity, gravity anomalies, gravity potential and gravity gradients for a set of points (e.g. satellites) in or above the model surface for a user-defined range of latitudes, longitudes and radius, or a list of point coordinates. Spherical coordinates in the output file are radius, colatitude and colongitude. Gravity is here based on the density distribution from the material model (and non adiabatic). This means that the density may come directly from an ascii file. This postprocessor also computes theoretical gravity and its derivatives, which corresponds to the analytical solution of gravity in the same geometry but filled with a reference density. The reference density is also used to determine the density difference for computing gravity anomalies. Thus one must carefully evaluate the meaning of the gravity anomaly output, because the solution may not reflect the actual gravity anomaly (due to differences in the assumed reference density). One way to guarantee correct gravity anomalies is to subtract the gravity of a certain point from the average gravity on the map. Another way is to directly use density anomalies for this postprocessor. The average- minimum- and maximum gravity acceleration and potential are written into the statistics file.

Definition at line 57 of file gravity_point_values.h.

Member Enumeration Documentation

§ SamplingScheme

Specify the sampling scheme determining if gravity calculation is performed.

Enumerator
map 
fibonacci_spiral 
list_of_points 

Definition at line 246 of file gravity_point_values.h.

Constructor & Destructor Documentation

§ GravityPointValues()

Constructor.

Member Function Documentation

§ initialize()

template<int dim>
void aspect::Postprocess::GravityPointValues< 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.

The default implementation of this function does nothing, but plugins that derive from this class (via the Interface classes of their respective plugin systems) may overload it if they want something to happen upon startup of the Simulator object to which the plugin contributes.

Reimplemented from aspect::Plugins::InterfaceBase.

§ execute()

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

Specify the creation of output_gravity.txt.

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

§ declare_parameters()

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

Declare the parameters this class takes through input files.

§ parse_parameters()

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

Read the parameters this class declares from the parameter file.

Reimplemented from aspect::Plugins::InterfaceBase.

§ serialize()

template<int dim>
template<class Archive >
void aspect::Postprocess::GravityPointValues< dim >::serialize ( Archive &  ar,
const unsigned int  version 
)

Serialize the contents of this class as far as they are not read from input parameter files.

§ save()

template<int dim>
void aspect::Postprocess::GravityPointValues< dim >::save ( std::map< std::string, std::string > &  status_strings) const
overridevirtual

Save the state of this object.

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

§ load()

template<int dim>
void aspect::Postprocess::GravityPointValues< dim >::load ( const std::map< std::string, std::string > &  status_strings)
overridevirtual

Restore the state of the object.

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

§ set_last_output_time()

template<int dim>
void aspect::Postprocess::GravityPointValues< dim >::set_last_output_time ( const double  current_time)
private

Set the time output was supposed to be written. In the simplest case, this is the previous last output time plus the interval, but in general we'd like to ensure that it is the largest supposed output time, which is smaller than the current time, to avoid falling behind with last_output_time and having to catch up once the time step becomes larger. This is done after every output.

Member Data Documentation

§ output_interval

template<int dim>
double aspect::Postprocess::GravityPointValues< dim >::output_interval
private

Interval between the generation of gravity output.

Definition at line 109 of file gravity_point_values.h.

§ last_output_time

template<int dim>
double aspect::Postprocess::GravityPointValues< dim >::last_output_time
private

A time (in seconds) at which the last graphical output was supposed to be produced. Used to check for the next necessary output time.

Definition at line 115 of file gravity_point_values.h.

§ maximum_timesteps_between_outputs

template<int dim>
unsigned int aspect::Postprocess::GravityPointValues< dim >::maximum_timesteps_between_outputs
private

Maximum number of steps between the generation of gravity output.

Definition at line 120 of file gravity_point_values.h.

§ last_output_timestep

template<int dim>
unsigned int aspect::Postprocess::GravityPointValues< dim >::last_output_timestep
private

Timestep at which the last gravity output was produced Used to check for the next necessary output time.

Definition at line 126 of file gravity_point_values.h.

§ output_file_number

template<int dim>
unsigned int aspect::Postprocess::GravityPointValues< dim >::output_file_number
private

Consecutively counted number indicating the how-manyth time we will create output the next time we get to it.

Definition at line 132 of file gravity_point_values.h.

§ end_time

template<int dim>
double aspect::Postprocess::GravityPointValues< dim >::end_time
private

end_time is taken from parameter file. It is used to tell the postprocessor to write gravity output at the end time.

Definition at line 138 of file gravity_point_values.h.

§ precision

template<int dim>
unsigned int aspect::Postprocess::GravityPointValues< dim >::precision
private

Set the precision of the gravity acceleration, potential and gradients in the gravity output and statistics file.

Definition at line 154 of file gravity_point_values.h.

§ quadrature_degree_increase

template<int dim>
unsigned int aspect::Postprocess::GravityPointValues< dim >::quadrature_degree_increase
private

Quadrature degree increase over the velocity element degree may be required when gravity is calculated near the surface or inside the model. An increase in the quadrature element adds accuracy to the gravity solution from noise due to the model grid.

Definition at line 162 of file gravity_point_values.h.

§ n_points_spiral

template<int dim>
unsigned int aspect::Postprocess::GravityPointValues< dim >::n_points_spiral
private

Parameter for the fibonacci spiral sampling scheme:

Definition at line 167 of file gravity_point_values.h.

§ n_points_radius

template<int dim>
unsigned int aspect::Postprocess::GravityPointValues< dim >::n_points_radius
private

Parameter for the map and fibonacci spiral sampling scheme: Gravity may be calculated for a sets of points along the radius (e.g. depth profile) between a minimum and maximum radius. Number of points along the radius is specified with n_points_radius.

Definition at line 175 of file gravity_point_values.h.

§ n_points_longitude

template<int dim>
unsigned int aspect::Postprocess::GravityPointValues< dim >::n_points_longitude
private

Parameter for the map sampling scheme: Gravity may be calculated for a sets of points along the longitude (e.g. satellite mapping) between a minimum and maximum longitude. Number of points along the longitude is specified with n_points_longitude.

Definition at line 183 of file gravity_point_values.h.

§ n_points_latitude

template<int dim>
unsigned int aspect::Postprocess::GravityPointValues< dim >::n_points_latitude
private

Parameter for the map sampling scheme: Gravity may be calculated for a sets of points along the latitude (e.g. satellite mapping) between a minimum and maximum latitude. Number of points along the latitude is specified with n_points_latitude.

Definition at line 191 of file gravity_point_values.h.

§ minimum_radius

template<int dim>
double aspect::Postprocess::GravityPointValues< dim >::minimum_radius
private

Parameter for the map and fibonacci spiral sampling scheme: Prescribe a minimum radius for a sampling coverage at a specific height. May be set in- or outside the model domain.

Definition at line 198 of file gravity_point_values.h.

§ maximum_radius

template<int dim>
double aspect::Postprocess::GravityPointValues< dim >::maximum_radius
private

Parameter for the map and fibonacci spiral sampling scheme: Maximum radius for the radius range. May be set in- or outside the model domain. No need to specify maximum_radius if n_points_radius is 1.

Definition at line 206 of file gravity_point_values.h.

§ minimum_colongitude

template<int dim>
double aspect::Postprocess::GravityPointValues< dim >::minimum_colongitude
private

Parameter for the map sampling scheme: Minimum longitude for longitude range.

Definition at line 212 of file gravity_point_values.h.

§ maximum_colongitude

template<int dim>
double aspect::Postprocess::GravityPointValues< dim >::maximum_colongitude
private

Parameter for the map sampling scheme: Maximum longitude for the longitude range. No need to specify maximum_longitude if n_points_longitude is 1.

Definition at line 219 of file gravity_point_values.h.

§ minimum_colatitude

template<int dim>
double aspect::Postprocess::GravityPointValues< dim >::minimum_colatitude
private

Parameter for the map sampling scheme: Minimum latitude for the latitude range.

Definition at line 225 of file gravity_point_values.h.

§ maximum_colatitude

template<int dim>
double aspect::Postprocess::GravityPointValues< dim >::maximum_colatitude
private

Parameter for the map sampling scheme: Maximum latitude for the latitude range. No need to specify maximum_latitude if n_points_latitude is 1.

Definition at line 232 of file gravity_point_values.h.

§ reference_density

template<int dim>
double aspect::Postprocess::GravityPointValues< dim >::reference_density
private

A reference density is required for two purposes: 1) benchmark the gravity postprocessor with computing analytically gravity acceleration and gradients in the domain filled with a constant density. 2) calculate the density difference of the density given by the material model and a constant density, in order to compute gravity anomalies.

Definition at line 241 of file gravity_point_values.h.

§ sampling_scheme

§ radius_list

template<int dim>
std::vector<double> aspect::Postprocess::GravityPointValues< dim >::radius_list
private

Parameter for the list of points sampling scheme: List of radius coordinates for the list of points sampling scheme. Must follow the same order as the lists of longitude and latitude.

Definition at line 258 of file gravity_point_values.h.

§ longitude_list

template<int dim>
std::vector<double> aspect::Postprocess::GravityPointValues< dim >::longitude_list
private

Parameter for the list of points sampling scheme: List of longitude coordinates for the list of points sampling scheme. Must follow the same order as the lists of radius and latitude.

Definition at line 265 of file gravity_point_values.h.

§ latitude_list

template<int dim>
std::vector<double> aspect::Postprocess::GravityPointValues< dim >::latitude_list
private

Parameter for the list of points sampling scheme: List of latitude coordinates for the list of points sampling scheme. Must follow the same order as the lists of longitude and longitude.

Definition at line 272 of file gravity_point_values.h.

§ model_outer_radius

template<int dim>
double aspect::Postprocess::GravityPointValues< dim >::model_outer_radius
private

Definition at line 277 of file gravity_point_values.h.

§ model_inner_radius

template<int dim>
double aspect::Postprocess::GravityPointValues< dim >::model_inner_radius
private

Definition at line 278 of file gravity_point_values.h.

§ satellite_positions_spherical

template<int dim>
std::vector<std::array<double,dim> > aspect::Postprocess::GravityPointValues< dim >::satellite_positions_spherical
private

The positions of all satellite positions in spherical and Cartesian coordinate systems.

Definition at line 284 of file gravity_point_values.h.

§ satellite_positions_cartesian

template<int dim>
std::vector<Point<dim> > aspect::Postprocess::GravityPointValues< dim >::satellite_positions_cartesian
private

Definition at line 285 of file gravity_point_values.h.


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