ASPECT
Classes | Public Types | Public Member Functions | Static Public Member Functions | Private Attributes | Friends | List of all members
aspect::Simulator< dim > Class Template Reference

Classes

struct  AdvectionField
 
struct  IntermediaryConstructorAction
 

Public Types

using NonlinearSolver = typename Parameters< dim >::NonlinearSolver
 
using NullspaceRemoval = typename Parameters< dim >::NullspaceRemoval
 

Public Member Functions

 Simulator (const MPI_Comm mpi_communicator, ParameterHandler &prm)
 
 ~Simulator ()
 
void run ()
 
void write_plugin_graph (std::ostream &output_stream) const
 

Static Public Member Functions

static void declare_parameters (ParameterHandler &prm)
 

Private Member Functions

Top-level functions in the overall flow of the numerical

algorithm

void setup_dofs ()
 
void setup_introspection ()
 
void set_initial_temperature_and_compositional_fields ()
 
void compute_initial_pressure_field ()
 
void compute_initial_velocity_boundary_constraints (AffineConstraints< double > &constraints)
 
void compute_current_velocity_boundary_constraints (AffineConstraints< double > &constraints)
 
void compute_current_constraints ()
 
double compute_pressure_scaling_factor () const
 
void start_timestep ()
 
void solve_timestep ()
 
void solve_single_advection_single_stokes ()
 
void solve_no_advection_iterated_stokes ()
 
void solve_no_advection_single_stokes ()
 
void solve_first_timestep_only_single_stokes ()
 
void solve_iterated_advection_and_stokes ()
 
void solve_single_advection_iterated_stokes ()
 
void solve_no_advection_iterated_defect_correction_stokes ()
 
void solve_single_advection_iterated_defect_correction_stokes ()
 
void solve_iterated_advection_and_defect_correction_stokes ()
 
void solve_iterated_advection_and_newton_stokes ()
 
void solve_single_advection_and_iterated_newton_stokes ()
 
void solve_single_advection_no_stokes ()
 
void solve_no_advection_no_stokes ()
 
void build_stokes_preconditioner ()
 
void build_advection_preconditioner (const AdvectionField &advection_field, aspect::LinearAlgebra::PreconditionILU &preconditioner, const double diagonal_strengthening)
 
void assemble_stokes_system ()
 
double assemble_and_solve_temperature (const double &initial_residual=0, double *residual=nullptr)
 
std::vector< double > assemble_and_solve_composition (const std::vector< double > &initial_residual={}, std::vector< double > *residual=nullptr)
 
double assemble_and_solve_stokes (const double &initial_nonlinear_residual=0, double *nonlinear_residual=nullptr)
 
void do_one_defect_correction_Stokes_step (DefectCorrectionResiduals &dcr, const bool use_picard)
 
void assemble_advection_system (const AdvectionField &advection_field)
 
double solve_advection (const AdvectionField &advection_field)
 
DEAL_II_DEPRECATED void interpolate_particle_properties (const AdvectionField &advection_field)
 
void interpolate_particle_properties (const std::vector< AdvectionField > &advection_fields)
 
std::pair< double, double > solve_stokes ()
 
std::pair< double, double > solve_stokes_block_gmg ()
 
void postprocess ()
 
void refine_mesh (const unsigned int max_grid_level)
 
Functions used in saving the state of the program and

restarting from a saved state

void create_snapshot ()
 
void resume_from_snapshot ()
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 
Functions used in setting up linear systems
Table< 2, DoFTools::Coupling > setup_system_matrix_coupling () const
 
void setup_system_matrix (const std::vector< IndexSet > &system_partitioning)
 
void setup_system_preconditioner (const std::vector< IndexSet > &system_partitioning)
 
Helper functions
void make_pressure_rhs_compatible (LinearAlgebra::BlockVector &vector)
 
template<typename T >
void get_artificial_viscosity (Vector< T > &viscosity_per_cell, const AdvectionField &advection_field, const bool skip_interior_cells=false) const
 
void compute_Vs_anomaly (Vector< float > &values) const
 
void compute_Vp_anomaly (Vector< float > &values) const
 
double normalize_pressure (LinearAlgebra::BlockVector &vector) const
 
void denormalize_pressure (const double pressure_adjustment, LinearAlgebra::BlockVector &vector, const LinearAlgebra::BlockVector &relevant_vector) const
 
void apply_limiter_to_dg_solutions (const AdvectionField &advection_field)
 
void compute_reactions ()
 
void update_solution_vectors_with_reaction_results (const unsigned int block_index, const LinearAlgebra::BlockVector &distributed_vector, const LinearAlgebra::BlockVector &distributed_reaction_vector)
 
void initialize_current_linearization_point ()
 
void interpolate_material_output_into_advection_field (const AdvectionField &adv_field)
 
void interpolate_onto_velocity_system (const TensorFunction< 1, dim > &func, LinearAlgebra::Vector &vec)
 
void setup_nullspace_constraints (AffineConstraints< double > &constraints)
 
void remove_nullspace (LinearAlgebra::BlockVector &relevant_dst, LinearAlgebra::BlockVector &tmp_distributed_stokes)
 
RotationProperties< dim > compute_net_angular_momentum (const bool use_constant_density, const LinearAlgebra::BlockVector &solution, const bool limit_to_top_faces=false) const
 
void remove_net_angular_momentum (const bool use_constant_density, LinearAlgebra::BlockVector &relevant_dst, LinearAlgebra::BlockVector &tmp_distributed_stokes, const bool limit_to_top_faces=false)
 
void replace_outflow_boundary_ids (const unsigned int boundary_id_offset)
 
void restore_outflow_boundary_ids (const unsigned int boundary_id_offset)
 
void remove_net_linear_momentum (const bool use_constant_density, LinearAlgebra::BlockVector &relevant_dst, LinearAlgebra::BlockVector &tmp_distributed_stokes)
 
double get_maximal_velocity (const LinearAlgebra::BlockVector &solution) const
 
double get_entropy_variation (const double average_field, const AdvectionField &advection_field) const
 
std::pair< double, double > get_extrapolated_advection_field_range (const AdvectionField &advection_field) const
 
void exchange_refinement_flags ()
 
void maybe_write_timing_output () const
 
bool maybe_write_checkpoint (const time_t last_checkpoint_time, const bool force_writing_checkpoint)
 
bool maybe_do_initial_refinement (const unsigned int max_refinement_level)
 
void maybe_refine_mesh (const double new_time_step, unsigned int &max_refinement_level)
 
void advance_time (const double step_size)
 
double compute_viscosity (internal::Assembly::Scratch::AdvectionSystem< dim > &scratch, const double global_u_infty, const double global_field_variation, const double average_field, const double global_entropy_variation, const double cell_diameter, const AdvectionField &advection_field) const
 
void compute_advection_system_residual (internal::Assembly::Scratch::AdvectionSystem< dim > &scratch, const double average_field, const AdvectionField &advection_field, double &max_residual, double &max_velocity, double &max_density, double &max_specific_heat, double &conductivity) const
 
bool stokes_matrix_depends_on_solution () const
 
void check_consistency_of_formulation ()
 
void check_consistency_of_boundary_conditions () const
 
double compute_initial_newton_residual (const LinearAlgebra::BlockVector &linearized_stokes_initial_guess)
 
double compute_Eisenstat_Walker_linear_tolerance (const bool EisenstatWalkerChoiceOne, const double maximum_linear_stokes_solver_tolerance, const double linear_stokes_solver_tolerance, const double stokes_residual, const double newton_residual, const double newton_residual_old)
 
void output_statistics ()
 
double compute_initial_stokes_residual ()
 

Private Attributes

bool simulator_is_past_initialization
 
std::unique_ptr< MeshDeformation::MeshDeformationHandler< dim > > mesh_deformation
 
std::unique_ptr< StokesMatrixFreeHandler< dim > > stokes_matrix_free
 
Variables that describe the physical setup of the problem
const std::unique_ptr< InitialTopographyModel::Interface< dim > > initial_topography_model
 
const std::unique_ptr< GeometryModel::Interface< dim > > geometry_model
 
const IntermediaryConstructorAction post_geometry_model_creation_action
 
const std::unique_ptr< MaterialModel::Interface< dim > > material_model
 
const std::unique_ptr< GravityModel::Interface< dim > > gravity_model
 
BoundaryTemperature::Manager< dim > boundary_temperature_manager
 
BoundaryComposition::Manager< dim > boundary_composition_manager
 
const std::unique_ptr< PrescribedStokesSolution::Interface< dim > > prescribed_stokes_solution
 
std::shared_ptr< InitialTemperature::Manager< dim > > initial_temperature_manager
 
std::shared_ptr< InitialComposition::Manager< dim > > initial_composition_manager
 
const std::unique_ptr< AdiabaticConditions::Interface< dim > > adiabatic_conditions
 
BoundaryVelocity::Manager< dim > boundary_velocity_manager
 
BoundaryTraction::Manager< dim > boundary_traction_manager
 
const std::unique_ptr< BoundaryHeatFlux::Interface< dim > > boundary_heat_flux
 
std::unique_ptr< Particle::World< dim > > particle_world
 
::Particles::ParticleHandler< dim > particle_handler_copy
 
Variables that describe the time discretization
double time
 
double time_step
 
double old_time_step
 
unsigned int timestep_number
 
unsigned int pre_refinement_step
 
unsigned int nonlinear_iteration
 
Variables related to simulation time stepping
TimeStepping::Manager< dim > time_stepping_manager
 
Variables for doing lateral averaging
LateralAveraging< dim > lateral_averaging
 
Variables that describe the spatial discretization
parallel::distributed::Triangulation< dim > triangulation
 
double global_Omega_diameter
 
double global_volume
 
MeshRefinement::Manager< dim > mesh_refinement_manager
 
HeatingModel::Manager< dim > heating_model_manager
 
std::unique_ptr< Mapping< dim > > mapping
 
const FESystem< dim > finite_element
 
DoFHandler< dim > dof_handler
 
Postprocess::Manager< dim > postprocess_manager
 
AffineConstraints< double > constraints
 
AffineConstraints< double > current_constraints
 
double last_pressure_normalization_adjustment
 
double pressure_scaling
 
bool do_pressure_rhs_compatibility_modification
 
Variables that describe the linear systems and solution vectors
LinearAlgebra::BlockSparseMatrix system_matrix
 
LinearAlgebra::BlockSparseMatrix system_preconditioner_matrix
 
LinearAlgebra::BlockVector solution
 
LinearAlgebra::BlockVector old_solution
 
LinearAlgebra::BlockVector old_old_solution
 
LinearAlgebra::BlockVector system_rhs
 
LinearAlgebra::BlockVector current_linearization_point
 
LinearAlgebra::BlockVector pressure_shape_function_integrals
 
LinearAlgebra::BlockVector operator_split_reaction_vector
 
std::unique_ptr< LinearAlgebra::PreconditionAMGAmg_preconditioner
 
std::unique_ptr< LinearAlgebra::PreconditionBaseMp_preconditioner
 
bool rebuild_sparsity_and_matrices
 
bool rebuild_stokes_matrix
 
bool assemble_newton_stokes_matrix
 
bool assemble_newton_stokes_system
 
bool rebuild_stokes_preconditioner
 

Friends

class boost::serialization::access
 
class SimulatorAccess< dim >
 
class MeshDeformation::MeshDeformationHandler< dim >
 
class VolumeOfFluidHandler< dim >
 
class StokesMatrixFreeHandler< dim >
 
template<int dimension, int velocity_degree>
class StokesMatrixFreeHandlerImplementation
 
struct Parameters< dim >
 

Functions, classes, and variables used in the assembly of linear systems

std::unique_ptr< Assemblers::Manager< dim > > assemblers
 
void set_assemblers ()
 
void set_advection_assemblers ()
 
void set_stokes_assemblers ()
 
void assemble_stokes_preconditioner ()
 
void local_assemble_stokes_preconditioner (const typename DoFHandler< dim >::active_cell_iterator &cell, internal::Assembly::Scratch::StokesPreconditioner< dim > &scratch, internal::Assembly::CopyData::StokesPreconditioner< dim > &data)
 
void copy_local_to_global_stokes_preconditioner (const internal::Assembly::CopyData::StokesPreconditioner< dim > &data)
 
void local_assemble_stokes_system (const typename DoFHandler< dim >::active_cell_iterator &cell, internal::Assembly::Scratch::StokesSystem< dim > &scratch, internal::Assembly::CopyData::StokesSystem< dim > &data)
 
void copy_local_to_global_stokes_system (const internal::Assembly::CopyData::StokesSystem< dim > &data)
 
void local_assemble_advection_face_terms (const AdvectionField &advection_field, const typename DoFHandler< dim >::active_cell_iterator &cell, internal::Assembly::Scratch::AdvectionSystem< dim > &scratch, internal::Assembly::CopyData::AdvectionSystem< dim > &data)
 
void local_assemble_advection_system (const AdvectionField &advection_field, const Vector< double > &viscosity_per_cell, const typename DoFHandler< dim >::active_cell_iterator &cell, internal::Assembly::Scratch::AdvectionSystem< dim > &scratch, internal::Assembly::CopyData::AdvectionSystem< dim > &data)
 
void copy_local_to_global_advection_system (const AdvectionField &advection_field, const internal::Assembly::CopyData::AdvectionSystem< dim > &data)
 

Variables that have to do with input, output, parallel

communication and interfacing with other parts of the program.

using TeeDevice = boost::iostreams::tee_device< std::ostream, std::ofstream >
 
using TeeStream = boost::iostreams::stream< TeeDevice >
 
Parameters< dim > parameters
 
std::unique_ptr< MeltHandler< dim > > melt_handler
 
std::unique_ptr< NewtonHandler< dim > > newton_handler
 
SimulatorSignals< dim > signals
 
const IntermediaryConstructorAction post_signal_creation
 
std::unique_ptr< VolumeOfFluidHandler< dim > > volume_of_fluid_handler
 
Introspection< dim > introspection
 
MPI_Comm mpi_communicator
 
std::ofstream log_file_stream
 
TeeDevice iostream_tee_device
 
TeeStream iostream_tee_stream
 
ConditionalOStream pcout
 
TableHandler statistics
 
std::size_t statistics_last_write_size
 
std::size_t statistics_last_hash
 
TimerOutput computing_timer
 
Timer wall_timer
 
double total_walltime_until_last_snapshot
 
std::thread output_statistics_thread
 

Detailed Description

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

This is the main class of ASPECT. It implements the overall simulation algorithm using the numerical methods discussed in the papers and manuals that accompany ASPECT.

Definition at line 179 of file simulator.h.

Member Typedef Documentation

§ NonlinearSolver

template<int dim>
using aspect::Simulator< dim >::NonlinearSolver = typename Parameters<dim>::NonlinearSolver

Import Nonlinear Solver type.

Definition at line 245 of file simulator.h.

§ NullspaceRemoval

template<int dim>
using aspect::Simulator< dim >::NullspaceRemoval = typename Parameters<dim>::NullspaceRemoval

Import nullspace removal type.

Definition at line 250 of file simulator.h.

§ TeeDevice

template<int dim>
using aspect::Simulator< dim >::TeeDevice = boost::iostreams::tee_device<std::ostream, std::ofstream>
private

Definition at line 1831 of file simulator.h.

§ TeeStream

template<int dim>
using aspect::Simulator< dim >::TeeStream = boost::iostreams::stream<TeeDevice>
private

Definition at line 1832 of file simulator.h.

Constructor & Destructor Documentation

§ Simulator()

template<int dim>
aspect::Simulator< dim >::Simulator ( const MPI_Comm  mpi_communicator,
ParameterHandler &  prm 
)

Constructor.

Parameters
mpi_communicatorThe MPI communicator on which this class is to work. The class creates a clone of the actual communicator to make its communications private from the rest of the world.
prmThe run-time parameter object from which this class obtains its settings.

This function is implemented in source/simulator/core.cc.

§ ~Simulator()

template<int dim>
aspect::Simulator< dim >::~Simulator ( )

Destructor. Destroy what needs to be destroyed after waiting for all threads that may still be doing something in the background.

Member Function Documentation

§ declare_parameters()

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

Declare the run-time parameters this class takes, and call the respective declare_parameters functions of the namespaces that describe geometries, material models, etc.

Parameters
prmThe object in which the run-time parameters are to be declared.

This function is implemented in source/simulator/parameters.cc.

§ run()

template<int dim>
void aspect::Simulator< dim >::run ( )

The function that runs the overall algorithm. It contains the loop over all time steps as well as the logic of what to do when before the loop starts and within the time loop.

This function is implemented in source/simulator/core.cc.

§ write_plugin_graph()

template<int dim>
void aspect::Simulator< dim >::write_plugin_graph ( std::ostream &  output_stream) const

Write a connection graph of all of the plugins we know about, in the format that the programs dot and neato understand. This allows for a visualization of how all of the plugins that ASPECT knows about are interconnected, and connect to other parts of the ASPECT code.

This function is implemented in source/simulator/helper_functions.cc.

Parameters
output_streamThe stream to write the output to.

§ setup_dofs()

template<int dim>
void aspect::Simulator< dim >::setup_dofs ( )
private

The function that sets up the DoFHandler objects, It also sets up the various partitioners and computes those constraints on the Stokes variable and temperature that are the same between all time steps.

This function is implemented in source/simulator/core.cc.

§ setup_introspection()

template<int dim>
void aspect::Simulator< dim >::setup_introspection ( )
private

This function initializes the variables of the introspection object. It is called by setup_dofs() right after distributing degrees of freedom since this is when all of the information is available.

This function is implemented in source/simulator/core.cc.

§ set_initial_temperature_and_compositional_fields()

template<int dim>
void aspect::Simulator< dim >::set_initial_temperature_and_compositional_fields ( )
private

A function that is responsible for initializing the temperature/compositional field before the first time step. The temperature field then serves as the temperature from which the velocity is computed during the first time step, and is subsequently overwritten by the temperature field one gets by advancing by one time step.

This function is implemented in source/simulator/initial_conditions.cc.

§ compute_initial_pressure_field()

template<int dim>
void aspect::Simulator< dim >::compute_initial_pressure_field ( )
private

A function that initializes the pressure variable before the first time step. It does so by either interpolating (for continuous pressure finite elements) or projecting (for discontinuous elements) the adiabatic pressure computed from the material model.

Note that the pressure so set is overwritten by the pressure in fact computed during the first time step. We need this function, however, so that the evaluation of pressure-dependent coefficients (e.g. pressure dependent densities or thermal coefficients) during the first time step has some useful pressure to start with.

This function is implemented in source/simulator/initial_conditions.cc.

§ compute_initial_velocity_boundary_constraints()

template<int dim>
void aspect::Simulator< dim >::compute_initial_velocity_boundary_constraints ( AffineConstraints< double > &  constraints)
private

Fill the given constraints with constraints coming from the velocity boundary conditions that do not change over time. This function is used by setup_dofs();

§ compute_current_velocity_boundary_constraints()

template<int dim>
void aspect::Simulator< dim >::compute_current_velocity_boundary_constraints ( AffineConstraints< double > &  constraints)
private

Fill the given constraints with constraints coming from the velocity boundary conditions that do can change over time. This function is used by compute_current_constraints().

§ compute_current_constraints()

template<int dim>
void aspect::Simulator< dim >::compute_current_constraints ( )
private

Given the 'constraints' member that contains all constraints that are independent of the time (e.g., hanging node constraints, tangential flow constraints, etc), copy it over to 'current_constraints' and add to the latter all constraints that do depend on time such as temperature or velocity Dirichlet boundary conditions. This function is therefore called at the beginning of every time step in start_timestep(), but also when setting up the initial values.

This function is implemented in source/simulator/core.cc.

§ compute_pressure_scaling_factor()

template<int dim>
double aspect::Simulator< dim >::compute_pressure_scaling_factor ( ) const
private

Compute the factor by which we scale the second of the Stokes equations (the "pressure scaling factor"). We compute the factor by taking the logarithmic average of the viscosities we find on the cells in this domain and dividing this average by a reference length scale provided by the used geometry model.

This function returns the pressure scaling variable.

This function is implemented in source/simulator/helper_functions.cc.

§ start_timestep()

template<int dim>
void aspect::Simulator< dim >::start_timestep ( )
private

Do some housekeeping at the beginning of each time step. This includes generating some screen output, adding some information to the statistics file, and interpolating time-dependent boundary conditions specific to this particular time step (the time independent boundary conditions, for example for hanging nodes or for tangential flow, are computed only once per mesh in setup_dofs()).

This function is implemented in source/simulator/core.cc.

§ solve_timestep()

template<int dim>
void aspect::Simulator< dim >::solve_timestep ( )
private

Do the various steps necessary to assemble and solve the things necessary in each time step.

This function is implemented in source/simulator/core.cc.

§ solve_single_advection_single_stokes()

template<int dim>
void aspect::Simulator< dim >::solve_single_advection_single_stokes ( )
private

This function implements one scheme for the various steps necessary to assemble and solve the nonlinear problem.

If `single Advection, single Stokes' is selected as the nonlinear solver scheme, no nonlinear iterations are done, and the temperature, compositional fields and Stokes equations are solved exactly once per time step, one after the other.

This function is implemented in source/simulator/solver_schemes.cc.

§ solve_no_advection_iterated_stokes()

template<int dim>
void aspect::Simulator< dim >::solve_no_advection_iterated_stokes ( )
private

This function implements one scheme for the various steps necessary to assemble and solve the nonlinear problem.

The `no Advection, iterated Stokes' scheme only solves the Stokes system and ignores compositions and the temperature equation (careful, the material model must not depend on the temperature; mostly useful for Stokes benchmarks).

This function is implemented in source/simulator/solver_schemes.cc.

§ solve_no_advection_single_stokes()

template<int dim>
void aspect::Simulator< dim >::solve_no_advection_single_stokes ( )
private

This function implements one scheme for the various steps necessary to assemble and solve the nonlinear problem.

The `no Advection, single Stokes' scheme only solves the Stokes system and ignores compositions and the temperature equation.

This function is implemented in source/simulator/solver_schemes.cc.

§ solve_first_timestep_only_single_stokes()

template<int dim>
void aspect::Simulator< dim >::solve_first_timestep_only_single_stokes ( )
private

This function implements one scheme for the various steps necessary to assemble and solve the nonlinear problem.

The `first timestep only, single Stokes' scheme only solves the Stokes system, for the initial timestep. This results in a `steady state' velocity field for particle calculations.

This function is implemented in source/simulator/solver_schemes.cc.

§ solve_iterated_advection_and_stokes()

template<int dim>
void aspect::Simulator< dim >::solve_iterated_advection_and_stokes ( )
private

This function implements one scheme for the various steps necessary to assemble and solve the nonlinear problem.

The `iterated Advection and Stokes' scheme iterates by alternating the solution of the temperature, composition and Stokes systems. This is essentially a type of Picard iterations for the whole system of equations.

This function is implemented in source/simulator/solver_schemes.cc.

§ solve_single_advection_iterated_stokes()

template<int dim>
void aspect::Simulator< dim >::solve_single_advection_iterated_stokes ( )
private

This function implements one scheme for the various steps necessary to assemble and solve the nonlinear problem.

The `single Advection, iterated Stokes' scheme solves the temperature and composition equations once at the beginning of each time step and then iterates out the solution of the Stokes equation using Picard iterations.

This function is implemented in source/simulator/solver_schemes.cc.

§ solve_no_advection_iterated_defect_correction_stokes()

template<int dim>
void aspect::Simulator< dim >::solve_no_advection_iterated_defect_correction_stokes ( )
private

This function implements one scheme for the various steps necessary to assemble and solve the nonlinear problem.

The `no Advection, iterated defect correction Stokes' scheme does not solve the temperature and composition equations but only iterates out the solution of the Stokes equation using Defect Correction (DC) Picard iterations.

This function is implemented in source/simulator/solver_schemes.cc.

§ solve_single_advection_iterated_defect_correction_stokes()

template<int dim>
void aspect::Simulator< dim >::solve_single_advection_iterated_defect_correction_stokes ( )
private

This function implements one scheme for the various steps necessary to assemble and solve the nonlinear problem.

The `single Advection, iterated defect correction Stokes' scheme solves the temperature and composition equations once at the beginning of each time step and then iterates out the solution of the Stokes equation using Defect Correction (DC) Picard iterations.

This function is implemented in source/simulator/solver_schemes.cc.

§ solve_iterated_advection_and_defect_correction_stokes()

template<int dim>
void aspect::Simulator< dim >::solve_iterated_advection_and_defect_correction_stokes ( )
private

This function implements one scheme for the various steps necessary to assemble and solve the nonlinear problem.

The `iterated Advection and defect correction Stokes' scheme iterates over both the temperature and composition equations and the Stokes equations at the same time. The Stokes equation is iterated out using the Defect Correction (DC) form of the Picard iterations.

This function is implemented in source/simulator/solver_schemes.cc.

§ solve_iterated_advection_and_newton_stokes()

template<int dim>
void aspect::Simulator< dim >::solve_iterated_advection_and_newton_stokes ( )
private

This function implements one scheme for the various steps necessary to assemble and solve the nonlinear problem.

The `iterated Advection and Newton Stokes' scheme iterates over solving the temperature, composition, and Stokes equations just like `iterated Advection and Stokes', but for the Stokes system it is able to switch from a defect correction form of Picard iterations to Newton iterations after a certain tolerance or number of iterations is reached. This can greatly improve the convergence rate for particularly nonlinear viscosities.

This function is implemented in source/simulator/solver_schemes.cc.

§ solve_single_advection_and_iterated_newton_stokes()

template<int dim>
void aspect::Simulator< dim >::solve_single_advection_and_iterated_newton_stokes ( )
private

This function implements one scheme for the various steps necessary to assemble and solve the nonlinear problem.

The `single Advection, iterated Newton Stokes' scheme solves the temperature and composition equations once at the beginning of each time step and then iterates out the solution of the Stokes equation using Newton iterations. For the Stokes system it is able to switch from a defect correction form of Picard iterations to Newton iterations after a certain tolerance or number of iterations is reached. This can greatly improve the convergence rate for particularly nonlinear viscosities.

This function is implemented in source/simulator/solver_schemes.cc.

§ solve_single_advection_no_stokes()

template<int dim>
void aspect::Simulator< dim >::solve_single_advection_no_stokes ( )
private

This function implements one scheme for the various steps necessary to assemble and solve the nonlinear problem.

The `single Advection, no Stokes' scheme only solves the temperature and other advection systems and instead of solving for the Stokes system, a prescribed velocity and pressure is used."

This function is implemented in source/simulator/solver_schemes.cc.

§ solve_no_advection_no_stokes()

template<int dim>
void aspect::Simulator< dim >::solve_no_advection_no_stokes ( )
private

This function implements one scheme for the various steps necessary to assemble and solve the nonlinear problem.

The `no Advection, no Stokes' scheme skips solving the temperature, composition and Stokes equations, which permits to go directly to postprocessing after setting up the initial condition.

This function is implemented in source/simulator/solver_schemes.cc.

§ build_stokes_preconditioner()

template<int dim>
void aspect::Simulator< dim >::build_stokes_preconditioner ( )
private

Initiate the assembly of the Stokes preconditioner matrix via assemble_stokes_preconditioner(), then set up the data structures to actually build a preconditioner from this matrix.

This function is implemented in source/simulator/assembly.cc.

§ build_advection_preconditioner()

template<int dim>
void aspect::Simulator< dim >::build_advection_preconditioner ( const AdvectionField advection_field,
aspect::LinearAlgebra::PreconditionILU preconditioner,
const double  diagonal_strengthening 
)
private

Initialize the preconditioner for the advection equation of field index.

This function is implemented in source/simulator/assembly.cc.

§ assemble_stokes_system()

template<int dim>
void aspect::Simulator< dim >::assemble_stokes_system ( )
private

Initiate the assembly of the Stokes matrix and right hand side.

This function is implemented in source/simulator/assembly.cc.

§ assemble_and_solve_temperature()

template<int dim>
double aspect::Simulator< dim >::assemble_and_solve_temperature ( const double &  initial_residual = 0,
double *  residual = nullptr 
)
private

Assemble and solve the temperature equation. This function returns the residual after solving.

If the residual argument is not a nullptr, the function computes the residual and puts it into this variable. The function returns the current residual divided by the initial residual given as the first argument. The two arguments may point to the same variable, in which case the function first computes the residual and at the end scales that residual by itself, thus returning 1.0.

This function is implemented in source/simulator/solver_schemes.cc.

§ assemble_and_solve_composition()

template<int dim>
std::vector<double> aspect::Simulator< dim >::assemble_and_solve_composition ( const std::vector< double > &  initial_residual = {},
std::vector< double > *  residual = nullptr 
)
private

Solve the composition equations with whatever method is selected (fields or particles). This function returns the residuals for all fields after solving.

If the residual argument is not a nullptr, the function computes the residual and puts it into this variable. The function returns the current residual divided by the initial residual given as the first argument. The two arguments may point to the same variable, in which case the function first computes the residual and at the end scales that residual by itself, thus returning 1.0.

This function is implemented in source/simulator/solver_schemes.cc.

§ assemble_and_solve_stokes()

template<int dim>
double aspect::Simulator< dim >::assemble_and_solve_stokes ( const double &  initial_nonlinear_residual = 0,
double *  nonlinear_residual = nullptr 
)
private

Assemble and solve the Stokes equation. This function returns the nonlinear residual after solving.

If the residual argument is not a nullptr, the function computes the residual and puts it into this variable. The function returns the current residual divided by the initial residual given as the first argument. The two arguments may point to the same variable, in which case the function first computes the residual and at the end scales that residual by itself, thus returning 1.0.

This function is implemented in source/simulator/solver_schemes.cc.

§ do_one_defect_correction_Stokes_step()

template<int dim>
void aspect::Simulator< dim >::do_one_defect_correction_Stokes_step ( DefectCorrectionResiduals dcr,
const bool  use_picard 
)
private

Do one step of the defect correction form of the Stokes equation; i.e., assemble and solve the defect correction equations. This function takes a structure of DefectCorrectionResiduals which contains information about different residuals. The information in this structure is updated by this function. The parameter use_picard forces the use of the defect correction Picard iteration, if set to 'true' no Newton derivatives are added to the matrix.

This function is implemented in source/simulator/solver_schemes.cc.

§ assemble_advection_system()

template<int dim>
void aspect::Simulator< dim >::assemble_advection_system ( const AdvectionField advection_field)
private

Initiate the assembly of one advection matrix and right hand side and build a preconditioner for the matrix.

This function is implemented in source/simulator/assembly.cc.

§ solve_advection()

template<int dim>
double aspect::Simulator< dim >::solve_advection ( const AdvectionField advection_field)
private

Solve one block of the temperature/composition linear system. Return the initial nonlinear residual, i.e., if the linear system to be solved is \(Ax=b\), then return \(\|Ax_0-b\|\) where \(x_0\) is the initial guess for the solution variable and is taken from the current_linearization_point member variable.

This function is implemented in source/simulator/solver.cc.

§ interpolate_particle_properties() [1/2]

template<int dim>
DEAL_II_DEPRECATED void aspect::Simulator< dim >::interpolate_particle_properties ( const AdvectionField advection_field)
private

Interpolate a particular particle property to the solution field.

Deprecated:
: Use interpolate_particle_property_vector() instead.

§ interpolate_particle_properties() [2/2]

template<int dim>
void aspect::Simulator< dim >::interpolate_particle_properties ( const std::vector< AdvectionField > &  advection_fields)
private

Interpolate the corresponding particle properties into the given advection_fields solution fields.

§ solve_stokes()

template<int dim>
std::pair<double,double> aspect::Simulator< dim >::solve_stokes ( )
private

Solve the Stokes linear system.

The function returns two pieces of information as a pair of doubles:

  • The initial nonlinear residual, i.e., if the linear system to be solved is \(Ax_{k+1}=b\), then we use \(\|Ax_k-b\|\) where \(x_k\) is the initial guess for the solution variable and is taken from the current_linearization_point member variable. For the purpose of this function, this residual is computed only from the velocity and pressure equations (i.e., for the 2x2 block system involving the velocity and pressure variables). A rationale for why this number is computed is given below.
  • The final linear residual, i.e., if the linear system to be solved is \(Ax_{k+1}=b\), then we use \(\|Ax_{k+1}-b\|\) where \(x_{k+1}\) is the solution just computed. If we use a direct solver to compute the solution of the linear system, then this linear residual is of course zero (or at least quite close to it) and the function just returns a zero value without even attempting to compute the actual value. On the other hand, if the function uses an iterative solver, then the value of the final linear residual is related to the tolerance with which we solve the linear system and generally indicates how accurately or inaccurately the linear system has been solved.

The two values are used in nonlinear solver schemes to assess how accurate the solution was before the current solve (for the first element of the returned pair) and how accurately the next iteration will have to be solved (for the second element of the pair) when using the Eisenstat-Walker method.

Note
If this function is called from a nonlinear solver – e.g., the `single Advection, iterated Stokes', or the `iterated Advection and Stokes' solvers schemes –, then the current_linearization_point is the solution of the previous iteration (or the solution extrapolated from the previous time steps, if this is the first nonlinear iteration). Let us call this solution \(x_k\) as above, where \(x\) is a two-component block vector that consists of velocity and pressure. This function then assumes that we have already built the system matrix \(A_k=A(x_k)\) and \(F_k=F(x_k)\), both linearized around the previous solution. The function solves the linear system \(A_k x_{k+1} = F_k\) for the solution \(x_{k+1}\). If the linear system were solved exactly, then that would imply that the linear residual \(\|A_k x_{k+1} - F_k\|\) were zero, or at least small. In other words, its size does not tell us anything about how accurately we have solved the nonlinear system. On the other hand, the nonlinear residual \(\|A_k x_k - F_k\|\) tells us something about how accurately the previous guess \(x_k\) already solved the nonlinear system. Consequently, this is what this function returns. (In some sense, this is not really what we are interested in: it tells us how accurate the solution already was, and if it was already pretty accurate, then we may not want to actually solve for \(x_{k+1}\). But, this would require that this function receives a tolerance so that it can bail out early without actually solving for \(x_{k+1}\) if the tolerance is already reached. This function does not actually do that – in some sense, one may argue that if we have already built the matrix and right hand side, we may as well solve with them, whether or not the solution was already good. If it happens to have been good already, then it will be even better after the solve. If it was not good enough yet, then we have to solve anyway.) In contrast to all of this, if we are using a Newton solver, then \(x_{k+1}\) is actually the Newton update vector, for which we have no initial guess other than the zero vector. In this case, the function simply returns \(\|F_k\|\) as the first element of the pair, where \(F_k=F(x_k)\) is the residual vector for the previous solution \(x_k\).

This function is implemented in source/simulator/solver.cc.

§ solve_stokes_block_gmg()

template<int dim>
std::pair<double,double> aspect::Simulator< dim >::solve_stokes_block_gmg ( )
private

Solve the Stokes system using a block preconditioner and GMG.

§ postprocess()

template<int dim>
void aspect::Simulator< dim >::postprocess ( )
private

This function is called at the end of every time step. It runs all the postprocessors that have been listed in the input parameter file (see the manual) in turn. In particular, this usually includes generating graphical output every few time steps.

The function also updates the statistics output file at the end of each time step.

This function is implemented in source/simulator/core.cc.

§ refine_mesh()

template<int dim>
void aspect::Simulator< dim >::refine_mesh ( const unsigned int  max_grid_level)
private

Refine the mesh according to error indicators calculated by compute_refinement_criterion(), set up all necessary data structures on this new mesh, and interpolate the old solutions onto the new mesh.

Parameters
[in]max_grid_levelThe maximum refinement level of the mesh. This is the sum of the initial global refinement and the initial adaptive refinement (as provided by the user in the input file) and in addition it gets increased by one at each additional refinement time.

This function is implemented in source/simulator/core.cc.

§ create_snapshot()

template<int dim>
void aspect::Simulator< dim >::create_snapshot ( )
private

Save the state of this program to a set of files in the output directory. In reality, however, only some variables are stored (in particular the mesh, the solution vectors, etc) whereas others can either be re-generated (matrices, DoFHandler objects, etc) or are read from the input parameter file. See the manual for more information.

This function is implemented in source/simulator/checkpoint_restart.cc.

§ resume_from_snapshot()

template<int dim>
void aspect::Simulator< dim >::resume_from_snapshot ( )
private

Restore the state of this program from a set of files in the output directory. In reality, however, only some variables are stored (in particular the mesh, the solution vectors, etc) whereas others can either be re-generated (matrices, DoFHandler objects, etc) or are read from the input parameter file. See the manual for more information. This function only restores those variables that can neither be re-generated from other information nor are read from the input parameter file.

This function is implemented in source/simulator/checkpoint_restart.cc.

§ serialize()

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

Save a number of variables using BOOST serialization mechanism.

This function is implemented in source/simulator/checkpoint_restart.cc.

§ setup_system_matrix_coupling()

template<int dim>
Table<2,DoFTools::Coupling> aspect::Simulator< dim >::setup_system_matrix_coupling ( ) const
private

Determine which of the components of our finite-element system couple to each other. Depending on which equations are solved and which solver is used this varies widely.

This function is implemented in source/simulator/core.cc.

§ setup_system_matrix()

template<int dim>
void aspect::Simulator< dim >::setup_system_matrix ( const std::vector< IndexSet > &  system_partitioning)
private

Set up the size and structure of the matrix used to store the elements of the linear system.

This function is implemented in source/simulator/core.cc.

§ setup_system_preconditioner()

template<int dim>
void aspect::Simulator< dim >::setup_system_preconditioner ( const std::vector< IndexSet > &  system_partitioning)
private

Set up the size and structure of the matrix used to store the elements of the matrix that is used to build the preconditioner for the system. This matrix is only used for the Stokes system, so while it has the size of the whole system, it only has entries in the velocity and pressure blocks.

This function is implemented in source/simulator/core.cc.

§ set_assemblers()

template<int dim>
void aspect::Simulator< dim >::set_assemblers ( )
private

Determine, based on the run-time parameters of the current simulation, which functions need to be called in order to assemble linear systems, matrices, and right hand side vectors.

This function is implemented in source/simulator/assembly.cc.

§ set_advection_assemblers()

template<int dim>
void aspect::Simulator< dim >::set_advection_assemblers ( )
private

Determine, based on the run-time parameters of the current simulation, which functions need to be called in order to assemble linear systems, matrices, and right hand side vectors for the advection. This function is used by both the default full Stokes solver and the Newton solvers, but not by the two-phase flow solver.

This function is implemented in source/simulator/assembly.cc.

§ set_stokes_assemblers()

template<int dim>
void aspect::Simulator< dim >::set_stokes_assemblers ( )
private

Determine, based on the run-time parameters of the current simulation, which functions need to be called in order to assemble linear systems, matrices, and right hand side vectors for the default full Stokes solver i.e. without considering two-phase flow, or Newton solvers.

This function is implemented in source/simulator/assembly.cc.

§ assemble_stokes_preconditioner()

template<int dim>
void aspect::Simulator< dim >::assemble_stokes_preconditioner ( )
private

Initiate the assembly of the preconditioner for the Stokes system.

This function is implemented in source/simulator/assembly.cc.

§ local_assemble_stokes_preconditioner()

template<int dim>
void aspect::Simulator< dim >::local_assemble_stokes_preconditioner ( const typename DoFHandler< dim >::active_cell_iterator &  cell,
internal::Assembly::Scratch::StokesPreconditioner< dim > &  scratch,
internal::Assembly::CopyData::StokesPreconditioner< dim > &  data 
)
private

Compute the integrals for the preconditioner for the Stokes system on a single cell.

This function is implemented in source/simulator/assembly.cc.

§ copy_local_to_global_stokes_preconditioner()

template<int dim>
void aspect::Simulator< dim >::copy_local_to_global_stokes_preconditioner ( const internal::Assembly::CopyData::StokesPreconditioner< dim > &  data)
private

Copy the contribution to the preconditioner for the Stokes system from a single cell into the global matrix that stores these elements.

This function is implemented in source/simulator/assembly.cc.

§ local_assemble_stokes_system()

template<int dim>
void aspect::Simulator< dim >::local_assemble_stokes_system ( const typename DoFHandler< dim >::active_cell_iterator &  cell,
internal::Assembly::Scratch::StokesSystem< dim > &  scratch,
internal::Assembly::CopyData::StokesSystem< dim > &  data 
)
private

Compute the integrals for the Stokes matrix and right hand side on a single cell.

This function is implemented in source/simulator/assembly.cc.

§ copy_local_to_global_stokes_system()

template<int dim>
void aspect::Simulator< dim >::copy_local_to_global_stokes_system ( const internal::Assembly::CopyData::StokesSystem< dim > &  data)
private

Copy the contribution to the Stokes system from a single cell into the global matrix that stores these elements.

This function is implemented in source/simulator/assembly.cc.

§ local_assemble_advection_face_terms()

template<int dim>
void aspect::Simulator< dim >::local_assemble_advection_face_terms ( const AdvectionField advection_field,
const typename DoFHandler< dim >::active_cell_iterator &  cell,
internal::Assembly::Scratch::AdvectionSystem< dim > &  scratch,
internal::Assembly::CopyData::AdvectionSystem< dim > &  data 
)
private

Compute the integrals for one advection matrix and right hand side on the faces of a single cell.

This function is implemented in source/simulator/assembly.cc.

§ local_assemble_advection_system()

template<int dim>
void aspect::Simulator< dim >::local_assemble_advection_system ( const AdvectionField advection_field,
const Vector< double > &  viscosity_per_cell,
const typename DoFHandler< dim >::active_cell_iterator &  cell,
internal::Assembly::Scratch::AdvectionSystem< dim > &  scratch,
internal::Assembly::CopyData::AdvectionSystem< dim > &  data 
)
private

Compute the integrals for one advection matrix and right hand side on a single cell.

This function is implemented in source/simulator/assembly.cc.

§ copy_local_to_global_advection_system()

template<int dim>
void aspect::Simulator< dim >::copy_local_to_global_advection_system ( const AdvectionField advection_field,
const internal::Assembly::CopyData::AdvectionSystem< dim > &  data 
)
private

Copy the contribution to the advection system from a single cell into the global matrix that stores these elements.

This function is implemented in source/simulator/assembly.cc.

§ make_pressure_rhs_compatible()

template<int dim>
void aspect::Simulator< dim >::make_pressure_rhs_compatible ( LinearAlgebra::BlockVector vector)
private

This routine adjusts the second block of the right hand side of a Stokes system (containing the term that comes from compressibility, so that the system becomes compatible: \(0=\int div u = \int g\). The vector to adjust is given as the argument of this function. This function makes use of the helper vector pressure_shape_function_integrals that contains \(h_i=(q_i,1)\) with the pressure functions \(q_i\) and we adjust the right hand side \(g\) by \(h_i \int g / |\Omega|\).

The purpose of this function is described in the second paper on the numerical methods in ASPECT.

This function is implemented in source/simulator/helper_functions.cc.

§ get_artificial_viscosity()

template<int dim>
template<typename T >
void aspect::Simulator< dim >::get_artificial_viscosity ( Vector< T > &  viscosity_per_cell,
const AdvectionField advection_field,
const bool  skip_interior_cells = false 
) const
private

Fills a vector with the artificial viscosity for the temperature or composition on each local cell.

Parameters
viscosity_per_cellOutput vector
advection_fieldDetermines whether this variable should select the temperature field or a compositional field.
skip_interior_cellsA boolean flag. If set to true the function will only compute the artificial viscosity in cells at boundaries.

§ compute_Vs_anomaly()

template<int dim>
void aspect::Simulator< dim >::compute_Vs_anomaly ( Vector< float > &  values) const
private

Compute the seismic shear wave speed, Vs anomaly per element. we compute the anomaly by computing a smoothed (over 200 km or so) laterally averaged temperature profile and associated seismic velocity that is then subtracted from the seismic velocity at the current pressure temperature conditions

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

§ compute_Vp_anomaly()

template<int dim>
void aspect::Simulator< dim >::compute_Vp_anomaly ( Vector< float > &  values) const
private

Compute the seismic pressure wave speed, Vp anomaly per element. we compute the anomaly by computing a smoothed (over 200 km or so) laterally averaged temperature profile and associated seismic velocity that is then subtracted from the seismic velocity at the current pressure temperature conditions

This function is implemented in source/simulator/helper_functions.cc.

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

§ normalize_pressure()

template<int dim>
double aspect::Simulator< dim >::normalize_pressure ( LinearAlgebra::BlockVector vector) const
private

Adjust the pressure variable (which is only determined up to a constant by the equations, though its value may enter traction boundary conditions) by adding a constant to it in such a way that the pressure on the surface or within the entire volume has a known average value. The point of this function is that the pressure that results from solving the linear system may not coincide with what we think of as the "physical pressure"; in particular, we typically think of the pressure as zero (on average) along the surface because it is the sum of the hydrostatic and dynamic pressure, where the former is thought of as zero along the surface. This function therefore converts from the "mathematical" pressure to the "physical" pressure so that all following postprocessing steps can use the latter.

In the case of the surface average, whether a face is part of the surface is determined by asking whether its depth of its midpoint (as determined by the geometry model) is less than 1/3*1/sqrt(dim-1)*diameter of the face. For reasonably curved boundaries, this rules out side faces that are perpendicular to the surface boundary but includes those faces that are along the boundary even if the real boundary is curved.

Whether the pressure should be normalized based on the surface or volume average is decided by a parameter in the input file.

Note
This function is called after setting the initial pressure field in compute_initial_pressure_field() and at the end of solve_stokes(). This makes sense because these are exactly the places where the pressure is modified or re-computed.
Returns
This function returns the pressure adjustment by value. This is so that its negative can later be used again in denormalize_pressure().

This function is implemented in source/simulator/helper_functions.cc.

§ denormalize_pressure()

template<int dim>
void aspect::Simulator< dim >::denormalize_pressure ( const double  pressure_adjustment,
LinearAlgebra::BlockVector vector,
const LinearAlgebra::BlockVector relevant_vector 
) const
private

Invert the action of the normalize_pressure() function above. This means that we move from a pressure that satisfies the pressure normalization (e.g., has a zero average pressure, or a zero average surface pressure) to one that does not actually satisfy this normalization, and this doesn't seem to make sense because we are not interested in such a pressure.

Indeed, this function is only called at the very beginning of solve_stokes() before we compute the initial (linear) residual of the linear system \(Ax=b\) that corresponds to the Stokes system, where \(x\) is the variable for which the pressure is adjusted back to the "wrong" form. Because stokes_system() calls normalize_pressure() at the end of its operations, no such "wrong" pressure ever escapes the realm of solve_stokes(). The "wrong" pressure is then used for two purposes in that function: (i) To compute the initial Stokes residual, which makes sense because it can only be zero (if we solved the same linear system twice) if we re-use the exact same pressure as we got from the previous solve – i.e., before we called normalize_pressure() at the end of the solve. (ii) To initialize the solution vector before calling the GMRES solver, which also makes sense because the best guess vector for GMRES is the one that had previously come out of GMRES, namely the one on which we later called normalize_pressure().

This function modifies vector in-place. In some cases, we need locally_relevant values of the pressure. To avoid creating a new vector and transferring data, this function uses a second vector with relevant dofs (relevant_vector) for accessing these pressure values. Both vector and relevant_vector are expected to already contain the correct pressure values.

Note
The adjustment made in this function is done using the negative of the pressure_adjustment function argument that would typically have been computed and returned by the normalize_pressure() function. This value is typically stored in the member variable last_pressure_normalization_adjustment, but the current function doesn't read this variable but instead gets the adjustment variable from the given argument.

This function is implemented in source/simulator/helper_functions.cc.

§ apply_limiter_to_dg_solutions()

template<int dim>
void aspect::Simulator< dim >::apply_limiter_to_dg_solutions ( const AdvectionField advection_field)
private

Apply the bound preserving limiter to the discontinuous Galerkin solutions: i.e., given two fixed upper and lower bound [min, max], after applying the limiter, the discontinuous Galerkin solution will stay in the prescribed bounds.

This function is implemented in source/simulator/helper_functions.cc.

§ compute_reactions()

template<int dim>
void aspect::Simulator< dim >::compute_reactions ( )
private

Compute the reactions in case of operator splitting: Using the current solution vector, this function makes a number of time steps determined by the size of the reaction time step, and solves a system of coupled ordinary differential equations for the reactions between compositional fields and temperature in each of them. To do that, is uses the reaction rates outputs from the material and heating models used in the computation. The solution vector is then updated with the new values of temperature and composition after the reactions.

As the ordinary differential equation in any given point is independent from the solution at all other points, we do not have to assemble a matrix, but just need to loop over all node locations for the temperature and compositional fields and compute the update to the solution.

The function also updates the old solution vectors with the reaction update so that the advection time stepping scheme will have the correct field terms for the right-hand side when assembling the advection system.

This function is implemented in source/simulator/helper_functions.cc.

§ update_solution_vectors_with_reaction_results()

template<int dim>
void aspect::Simulator< dim >::update_solution_vectors_with_reaction_results ( const unsigned int  block_index,
const LinearAlgebra::BlockVector distributed_vector,
const LinearAlgebra::BlockVector distributed_reaction_vector 
)
private

Update the indicated block of the solution vector with the corresponding block of the handed over distributed_vector. Also update reaction_vector with the corresponding block of distributed_reaction_vector.

This function is implemented in source/simulator/helper_functions.cc.

§ initialize_current_linearization_point()

template<int dim>
void aspect::Simulator< dim >::initialize_current_linearization_point ( )
private

Initialize the current linearization point vector from the old solution vector(s). Depending on the time of the call this can be simply a copy of the solution of the last timestep, or a linear extrapolation of old and old_old timestep to the new timestep.

This function is implemented in source/simulator/helper_functions.cc.

§ interpolate_material_output_into_advection_field()

template<int dim>
void aspect::Simulator< dim >::interpolate_material_output_into_advection_field ( const AdvectionField adv_field)
private

Interpolate material model outputs onto an advection field (temperature or composition). For the field identified by the AdvectionField adv_field, this function asks the material model to fill a MaterialModel::MaterialModelOutputs object that has an attached MaterialModel::PrescribedFieldOutputs "additional outputs" object (for composition) or an attached MaterialModel::PrescribedTemperatureOutputs (for temperature). The MaterialModel::MaterialModelInputs object passed to the material model then contains the support points of the advection field, thereby allowing the outputs to be interpolated to the finite element space and consequently into the solution vector. This is useful for advection fields whose advection mode is set to Parameters::AdvectionFieldMethod::prescribed_field or Parameters::AdvectionFieldMethod::prescribed_field_with_diffusion.

This function is implemented in source/simulator/helper_functions.cc.

§ interpolate_onto_velocity_system()

template<int dim>
void aspect::Simulator< dim >::interpolate_onto_velocity_system ( const TensorFunction< 1, dim > &  func,
LinearAlgebra::Vector vec 
)
private

Interpolate the given function onto the velocity FE space and write it into the given vector.

This function is implemented in source/simulator/helper_functions.cc.

§ setup_nullspace_constraints()

template<int dim>
void aspect::Simulator< dim >::setup_nullspace_constraints ( AffineConstraints< double > &  constraints)
private

Add constraints to the given constraints object that are required for unique solvability of the velocity block based on the nullspace removal settings.

This method will add a zero Dirichlet constraint for the first velocity unknown in the domain for each velocity component, which is later being processed for translational or linear momentum removal. This avoids breakdowns of the linear solvers that otherwise occurred in some instances.

Note
: Rotational modes are currently not handled and don't appear to require constraints so far.

§ remove_nullspace()

template<int dim>
void aspect::Simulator< dim >::remove_nullspace ( LinearAlgebra::BlockVector relevant_dst,
LinearAlgebra::BlockVector tmp_distributed_stokes 
)
private

Eliminate the nullspace of the velocity in the given vector. Both vectors are expected to contain the up to date data.

Parameters
relevant_dstlocally relevant vector for the whole FE, will be filled at the end.
tmp_distributed_stokesonly contains velocity and pressure.

This function is implemented in source/simulator/nullspace.cc.

§ compute_net_angular_momentum()

template<int dim>
RotationProperties<dim> aspect::Simulator< dim >::compute_net_angular_momentum ( const bool  use_constant_density,
const LinearAlgebra::BlockVector solution,
const bool  limit_to_top_faces = false 
) const
private

Compute the angular momentum and other rotation properties of the velocities in the given solution vector.

Parameters
use_constant_densitydetermines whether to use a constant density (which corresponds to computing a net rotation instead of net angular momentum).
solutionSolution vector to compute the properties for.
limit_to_top_facesallows to only compute the net angular momentum (or net rotation) of the top surface.

This function is implemented in source/simulator/nullspace.cc.

§ remove_net_angular_momentum()

template<int dim>
void aspect::Simulator< dim >::remove_net_angular_momentum ( const bool  use_constant_density,
LinearAlgebra::BlockVector relevant_dst,
LinearAlgebra::BlockVector tmp_distributed_stokes,
const bool  limit_to_top_faces = false 
)
private

Remove the angular momentum of the given vector

Parameters
use_constant_densitydetermines whether to use a constant density (which corresponds to removing a net rotation instead of net angular momentum).
relevant_dstlocally relevant vector for the whole FE, will be filled at the end.
tmp_distributed_stokesonly contains velocity and pressure.
limit_to_top_facesallows to only remove the net angular momentum (or net rotation) of the top surface. This can be useful to compare surface motions against plate reconstructions in no net rotation reference frames.

This function is implemented in source/simulator/nullspace.cc.

§ replace_outflow_boundary_ids()

template<int dim>
void aspect::Simulator< dim >::replace_outflow_boundary_ids ( const unsigned int  boundary_id_offset)
private

Offset the boundary id of all faces located on an outflow boundary by a fixed value given by the input parameter boundary_id_offset.

This function is implemented in source/simulator/helper_functions.cc.

§ restore_outflow_boundary_ids()

template<int dim>
void aspect::Simulator< dim >::restore_outflow_boundary_ids ( const unsigned int  boundary_id_offset)
private

Undo the offset of the boundary ids done in replace_outflow_boundary_ids by resetting all boundary ids to their original value.

This function is implemented in source/simulator/helper_functions.cc.

§ remove_net_linear_momentum()

template<int dim>
void aspect::Simulator< dim >::remove_net_linear_momentum ( const bool  use_constant_density,
LinearAlgebra::BlockVector relevant_dst,
LinearAlgebra::BlockVector tmp_distributed_stokes 
)
private

Remove the linear momentum of the given vector

Parameters
use_constant_densitydetermines whether to use a constant density (which corresponds to removing a net translation instead of net linear momentum).
relevant_dstlocally relevant vector for the whole FE, will be filled at the end.
tmp_distributed_stokesonly contains velocity and pressure.

This function is implemented in source/simulator/nullspace.cc.

§ get_maximal_velocity()

template<int dim>
double aspect::Simulator< dim >::get_maximal_velocity ( const LinearAlgebra::BlockVector solution) const
private

Compute the maximal velocity throughout the domain. This is needed to compute the size of the time step.

This function is implemented in source/simulator/helper_functions.cc.

§ get_entropy_variation()

template<int dim>
double aspect::Simulator< dim >::get_entropy_variation ( const double  average_field,
const AdvectionField advection_field 
) const
private

Compute the variation (i.e., the difference between maximal and minimal value) of the entropy \((T-\bar T)^2\) where \(\bar T\) is the average temperature throughout the domain given as argument to this function.

This function is used in computing the artificial diffusion stabilization term.

This function is implemented in source/simulator/assembly.cc.

§ get_extrapolated_advection_field_range()

template<int dim>
std::pair<double,double> aspect::Simulator< dim >::get_extrapolated_advection_field_range ( const AdvectionField advection_field) const
private

Compute the minimal and maximal temperature throughout the domain from a solution vector extrapolated from the previous time steps. This is needed to compute the artificial diffusion stabilization terms.

This function is implemented in source/simulator/helper_functions.cc.

§ exchange_refinement_flags()

template<int dim>
void aspect::Simulator< dim >::exchange_refinement_flags ( )
private

Exchange coarsen/refinement flags set between processors so that we have the correct settings on all ghost cells.

This function is implemented in source/simulator/helper_functions.cc.

§ maybe_write_timing_output()

template<int dim>
void aspect::Simulator< dim >::maybe_write_timing_output ( ) const
private

Check if timing output should be written in this timestep, and if so write it.

This function is implemented in source/simulator/helper_functions.cc.

§ maybe_write_checkpoint()

template<int dim>
bool aspect::Simulator< dim >::maybe_write_checkpoint ( const time_t  last_checkpoint_time,
const bool  force_writing_checkpoint 
)
private

Check if a checkpoint should be written in this timestep. If so create one. Returns whether a checkpoint was written.

This function is implemented in source/simulator/helper_functions.cc.

§ maybe_do_initial_refinement()

template<int dim>
bool aspect::Simulator< dim >::maybe_do_initial_refinement ( const unsigned int  max_refinement_level)
private

Check if we should do an initial refinement cycle in this timestep. This will only be checked in timestep 0, afterwards the variable pre_refinement_step variable is invalidated, and this function will return without doing refinement. An initial refinement cycle is different from a regular one, because time is not increased. Instead the same timestep is solved using the new mesh. Therefore, only output timing information and postprocessor output if required in the input file. But always output statistics (to have a history of the number of cells in the statistics file). This function returns whether an initial refinement was done.

This function is implemented in source/simulator/helper_functions.cc.

§ maybe_refine_mesh()

template<int dim>
void aspect::Simulator< dim >::maybe_refine_mesh ( const double  new_time_step,
unsigned int &  max_refinement_level 
)
private

Check if refinement is requested in this timestep. If so: Refine mesh. The max_refinement_level might be increased from this time on if this is an additional refinement cycle.

This function is implemented in source/simulator/helper_functions.cc.

§ advance_time()

template<int dim>
void aspect::Simulator< dim >::advance_time ( const double  step_size)
private

Advance the current time by the given step_size and update the solution vectors as needed.

This function is implemented in source/simulator/helper_functions.cc.

§ compute_viscosity()

template<int dim>
double aspect::Simulator< dim >::compute_viscosity ( internal::Assembly::Scratch::AdvectionSystem< dim > &  scratch,
const double  global_u_infty,
const double  global_field_variation,
const double  average_field,
const double  global_entropy_variation,
const double  cell_diameter,
const AdvectionField advection_field 
) const
private

Compute the artificial diffusion coefficient value on a cell given the values and gradients of the solution passed as arguments.

This function is implemented in source/simulator/assembly.cc.

§ compute_advection_system_residual()

template<int dim>
void aspect::Simulator< dim >::compute_advection_system_residual ( internal::Assembly::Scratch::AdvectionSystem< dim > &  scratch,
const double  average_field,
const AdvectionField advection_field,
double &  max_residual,
double &  max_velocity,
double &  max_density,
double &  max_specific_heat,
double &  conductivity 
) const
private

Compute the residual of one advection equation to be used for the artificial diffusion coefficient value on a cell given the values and gradients of the solution passed as arguments.

This function is implemented in source/simulator/assembly.cc.

§ stokes_matrix_depends_on_solution()

template<int dim>
bool aspect::Simulator< dim >::stokes_matrix_depends_on_solution ( ) const
private

Return whether the Stokes matrix depends on the values of the solution at the previous time step. This is the case is the coefficients that appear in the matrix (i.e., the viscosity and, in the case of a compressible model, the density) depend on the solution.

This function exists to ensure that the Stokes matrix is rebuilt in time steps where it may have changed, while we want to save the effort of rebuilding it whenever we don't need to.

This function is implemented in source/simulator/helper_functions.cc.

§ check_consistency_of_formulation()

template<int dim>
void aspect::Simulator< dim >::check_consistency_of_formulation ( )
private

This function checks that the user-selected formulations of the equations are consistent with the other inputs. If an incorrect selection is detected it throws an exception. It for example assures that correct heating terms are selected, and the material model supports the selection of the mass conservation formulation (e.g. incompressible)). If the parameter 'parameters.formulation' is set to 'custom' it only ensures very basic consistency.

This function is implemented in source/simulator/helper_functions.cc.

§ check_consistency_of_boundary_conditions()

template<int dim>
void aspect::Simulator< dim >::check_consistency_of_boundary_conditions ( ) const
private

This function checks that the user-selected boundary conditions do not contain contradictions. If an incorrect selection is detected it throws an exception. This for example assures that not both velocity and traction boundary conditions are prescribed at the same boundary, and that no boundary temperatures are prescribed at a periodic boundary.

This function is implemented in source/simulator/helper_functions.cc.

§ compute_initial_newton_residual()

template<int dim>
double aspect::Simulator< dim >::compute_initial_newton_residual ( const LinearAlgebra::BlockVector linearized_stokes_initial_guess)
private

Computes the initial Newton residual.

§ compute_Eisenstat_Walker_linear_tolerance()

template<int dim>
double aspect::Simulator< dim >::compute_Eisenstat_Walker_linear_tolerance ( const bool  EisenstatWalkerChoiceOne,
const double  maximum_linear_stokes_solver_tolerance,
const double  linear_stokes_solver_tolerance,
const double  stokes_residual,
const double  newton_residual,
const double  newton_residual_old 
)
private

This function computes the Eisenstat Walker linear tolerance used for the Newton iterations in the `iterated Advection and Newton Stokes' and `single Advection, iterated Newton Stokes' solver schemes. The Eisenstat and Walker (1996) method is used for determining the linear tolerance of the iteration after the first iteration. The paper gives two preferred choices of computing this tolerance. Both choices are implemented here with the suggested parameter values and safeguards.

§ output_statistics()

template<int dim>
void aspect::Simulator< dim >::output_statistics ( )
private

This function is called at the end of each time step and writes the statistics object that contains data like the current time, the number of linear solver iterations, and whatever the postprocessors have generated, to disk.

This function is implemented in source/simulator/helper_functions.cc.

§ compute_initial_stokes_residual()

template<int dim>
double aspect::Simulator< dim >::compute_initial_stokes_residual ( )
private

This routine computes the initial (nonlinear) Stokes residual that is needed as a convergence criterion in models with solver schemes that do nonlinear iterations. We calculate it in the same way as the tolerance for the linear solver, using the norm of the pressure RHS for the pressure part and a residual with zero velocity for the velocity part to get the part of the RHS not balanced by the static pressure.

This function is implemented in source/simulator/helper_functions.cc.

Friends And Related Function Documentation

§ boost::serialization::access

template<int dim>
friend class boost::serialization::access
friend

Definition at line 2127 of file simulator.h.

§ SimulatorAccess< dim >

template<int dim>
friend class SimulatorAccess< dim >
friend

Definition at line 2128 of file simulator.h.

§ MeshDeformation::MeshDeformationHandler< dim >

template<int dim>
friend class MeshDeformation::MeshDeformationHandler< dim >
friend

Definition at line 2129 of file simulator.h.

§ VolumeOfFluidHandler< dim >

template<int dim>
friend class VolumeOfFluidHandler< dim >
friend

Definition at line 2130 of file simulator.h.

§ StokesMatrixFreeHandler< dim >

template<int dim>
friend class StokesMatrixFreeHandler< dim >
friend

Definition at line 2131 of file simulator.h.

§ StokesMatrixFreeHandlerImplementation

template<int dim>
template<int dimension, int velocity_degree>
friend class StokesMatrixFreeHandlerImplementation
friend

Definition at line 2133 of file simulator.h.

§ Parameters< dim >

template<int dim>
friend struct Parameters< dim >
friend

Definition at line 2134 of file simulator.h.

Member Data Documentation

§ simulator_is_past_initialization

template<int dim>
bool aspect::Simulator< dim >::simulator_is_past_initialization
private

A member variable that tracks whether we are completely done with initialization and have started the time loop. This variable needs to be the first one that is initialized in the constructor, so that its value correctly tracks the status of the overall object. As a consequence, it has to be the first one declared in this class.

The variable is set to true just before we start the time stepping loop, but may be temporarily reset to false if, for example, we are during the initial mesh refinement steps where we start the time loop, but then go back to initialization steps (mesh refinement, interpolation of initial conditions, etc.) before re-starting the time loop.

This variable is queried by SimulatorAccess::simulator_is_past_initialization().

Definition at line 396 of file simulator.h.

§ assemblers

template<int dim>
std::unique_ptr<Assemblers::Manager<dim> > aspect::Simulator< dim >::assemblers
private

A member variable that stores, for the current simulation, what functions need to be called in order to assemble linear systems, matrices, and right hand side vectors.

One would probably want this variable to just be a member of type Assemblers::Manager<dim>, but this requires that this type is declared in the current scope, and that would require including <simulator/assemblers/interface.h> which we don't want because it's big. Consequently, we just store a pointer to such an object, and create the object pointed to at the top of set_assemblers().

Definition at line 1070 of file simulator.h.

§ parameters

template<int dim>
Parameters<dim> aspect::Simulator< dim >::parameters
private

Definition at line 1791 of file simulator.h.

§ melt_handler

template<int dim>
std::unique_ptr<MeltHandler<dim> > aspect::Simulator< dim >::melt_handler
private

Unique pointer for an instance of the MeltHandler. This way, if we do not need the machinery for doing melt stuff, we do not even allocate it.

Definition at line 1798 of file simulator.h.

§ newton_handler

template<int dim>
std::unique_ptr<NewtonHandler<dim> > aspect::Simulator< dim >::newton_handler
private

Unique pointer for an instance of the NewtonHandler. This way, if we do not need the machinery for doing Newton stuff, we do not even allocate it.

Definition at line 1805 of file simulator.h.

§ signals

template<int dim>
SimulatorSignals<dim> aspect::Simulator< dim >::signals
private

Definition at line 1807 of file simulator.h.

§ post_signal_creation

template<int dim>
const IntermediaryConstructorAction aspect::Simulator< dim >::post_signal_creation
private

Definition at line 1809 of file simulator.h.

§ volume_of_fluid_handler

template<int dim>
std::unique_ptr<VolumeOfFluidHandler<dim> > aspect::Simulator< dim >::volume_of_fluid_handler
private

Unique pointer for an instance of the VolumeOfFluidHandler. This way, if we do not need the machinery for doing volume_of_fluid stuff, we do not even allocate it.

Located here due to needing signals access

Definition at line 1818 of file simulator.h.

§ introspection

template<int dim>
Introspection<dim> aspect::Simulator< dim >::introspection
private

Definition at line 1820 of file simulator.h.

§ mpi_communicator

template<int dim>
MPI_Comm aspect::Simulator< dim >::mpi_communicator
private

Definition at line 1823 of file simulator.h.

§ log_file_stream

template<int dim>
std::ofstream aspect::Simulator< dim >::log_file_stream
private

This stream will log into the file output/log.txt (used automatically by pcout).

Definition at line 1829 of file simulator.h.

§ iostream_tee_device

template<int dim>
TeeDevice aspect::Simulator< dim >::iostream_tee_device
private

Definition at line 1834 of file simulator.h.

§ iostream_tee_stream

template<int dim>
TeeStream aspect::Simulator< dim >::iostream_tee_stream
private

Definition at line 1835 of file simulator.h.

§ pcout

template<int dim>
ConditionalOStream aspect::Simulator< dim >::pcout
private

Output stream for logging information. Will only output on processor 0.

Definition at line 1841 of file simulator.h.

§ statistics

template<int dim>
TableHandler aspect::Simulator< dim >::statistics
private

An object that stores a bunch of statistics such as the number of linear solver iterations, the time corresponding to each time step, etc, as well as whatever the various postprocessors want to put into it.

This variable is written to disk after every time step, by the Simulator::output_statistics() function.

Definition at line 1852 of file simulator.h.

§ statistics_last_write_size

template<int dim>
std::size_t aspect::Simulator< dim >::statistics_last_write_size
private

The following two variables keep track which parts of the statistics object have already been written. This is because the TableHandler class has no way to keep track what it has already written, and so we can not just append the last row of the table to the output file. Rather, we keep track how many bytes we already wrote, and a hash of what they contained, and if these so-many bytes have not changed between the previous and current write operation, then we only open the file in 'append' mode to add the new bytes from the last row. If what we would write now has changed from what we wrote back then in the first so-many bytes (e.g., because column widths of the table have changed), then we just replace the previous file by the current table contents in their entirety.

Definition at line 1868 of file simulator.h.

§ statistics_last_hash

template<int dim>
std::size_t aspect::Simulator< dim >::statistics_last_hash
private

Definition at line 1869 of file simulator.h.

§ computing_timer

template<int dim>
TimerOutput aspect::Simulator< dim >::computing_timer
mutableprivate

Definition at line 1871 of file simulator.h.

§ wall_timer

template<int dim>
Timer aspect::Simulator< dim >::wall_timer
private

A timer used to track the current wall time since the last snapshot (or since the program started).

Definition at line 1877 of file simulator.h.

§ total_walltime_until_last_snapshot

template<int dim>
double aspect::Simulator< dim >::total_walltime_until_last_snapshot
private

The total wall time that has elapsed up to the last snapshot that was created.

Definition at line 1883 of file simulator.h.

§ output_statistics_thread

template<int dim>
std::thread aspect::Simulator< dim >::output_statistics_thread
private

In output_statistics(), where we output the statistics object above, we do the actual writing on a separate thread. This variable is the handle we get for this thread so that we can wait for it to finish, either if we want to write the statistics object for the next thread, or if we want to terminate altogether.

Definition at line 1892 of file simulator.h.

§ initial_topography_model

template<int dim>
const std::unique_ptr<InitialTopographyModel::Interface<dim> > aspect::Simulator< dim >::initial_topography_model
private

Definition at line 1902 of file simulator.h.

§ geometry_model

template<int dim>
const std::unique_ptr<GeometryModel::Interface<dim> > aspect::Simulator< dim >::geometry_model
private

Definition at line 1903 of file simulator.h.

§ post_geometry_model_creation_action

template<int dim>
const IntermediaryConstructorAction aspect::Simulator< dim >::post_geometry_model_creation_action
private

Definition at line 1904 of file simulator.h.

§ material_model

template<int dim>
const std::unique_ptr<MaterialModel::Interface<dim> > aspect::Simulator< dim >::material_model
private

Definition at line 1905 of file simulator.h.

§ gravity_model

template<int dim>
const std::unique_ptr<GravityModel::Interface<dim> > aspect::Simulator< dim >::gravity_model
private

Definition at line 1906 of file simulator.h.

§ boundary_temperature_manager

template<int dim>
BoundaryTemperature::Manager<dim> aspect::Simulator< dim >::boundary_temperature_manager
private

Definition at line 1907 of file simulator.h.

§ boundary_composition_manager

template<int dim>
BoundaryComposition::Manager<dim> aspect::Simulator< dim >::boundary_composition_manager
private

Definition at line 1908 of file simulator.h.

§ prescribed_stokes_solution

template<int dim>
const std::unique_ptr<PrescribedStokesSolution::Interface<dim> > aspect::Simulator< dim >::prescribed_stokes_solution
private

Definition at line 1909 of file simulator.h.

§ initial_temperature_manager

template<int dim>
std::shared_ptr<InitialTemperature::Manager<dim> > aspect::Simulator< dim >::initial_temperature_manager
private

The following two variables are pointers to objects that describe the initial temperature and composition values. The Simulator class itself releases these pointers once they are no longer needed, somewhere during the first time step once it is known that they are no longer needed. However, plugins can have their own shared pointers to these objects, and the lifetime of the objects pointed to is then until the last of these plugins gets deleted.

Definition at line 1921 of file simulator.h.

§ initial_composition_manager

template<int dim>
std::shared_ptr<InitialComposition::Manager<dim> > aspect::Simulator< dim >::initial_composition_manager
private

Definition at line 1922 of file simulator.h.

§ adiabatic_conditions

template<int dim>
const std::unique_ptr<AdiabaticConditions::Interface<dim> > aspect::Simulator< dim >::adiabatic_conditions
private

Definition at line 1924 of file simulator.h.

§ boundary_velocity_manager

template<int dim>
BoundaryVelocity::Manager<dim> aspect::Simulator< dim >::boundary_velocity_manager
private

Definition at line 1939 of file simulator.h.

§ boundary_traction_manager

template<int dim>
BoundaryTraction::Manager<dim> aspect::Simulator< dim >::boundary_traction_manager
private

Definition at line 1940 of file simulator.h.

§ boundary_heat_flux

template<int dim>
const std::unique_ptr<BoundaryHeatFlux::Interface<dim> > aspect::Simulator< dim >::boundary_heat_flux
private

Definition at line 1941 of file simulator.h.

§ particle_world

template<int dim>
std::unique_ptr<Particle::World<dim> > aspect::Simulator< dim >::particle_world
private

The world holding the particles

Definition at line 1946 of file simulator.h.

§ particle_handler_copy

template<int dim>
::Particles::ParticleHandler<dim> aspect::Simulator< dim >::particle_handler_copy
private

A copy of the particle handler to reset the particles when repeating a time step.

Definition at line 1952 of file simulator.h.

§ time

template<int dim>
double aspect::Simulator< dim >::time
private

Definition at line 1961 of file simulator.h.

§ time_step

template<int dim>
double aspect::Simulator< dim >::time_step
private

Definition at line 1962 of file simulator.h.

§ old_time_step

template<int dim>
double aspect::Simulator< dim >::old_time_step
private

Definition at line 1963 of file simulator.h.

§ timestep_number

template<int dim>
unsigned int aspect::Simulator< dim >::timestep_number
private

Definition at line 1964 of file simulator.h.

§ pre_refinement_step

template<int dim>
unsigned int aspect::Simulator< dim >::pre_refinement_step
private

Definition at line 1965 of file simulator.h.

§ nonlinear_iteration

template<int dim>
unsigned int aspect::Simulator< dim >::nonlinear_iteration
private

Definition at line 1966 of file simulator.h.

§ time_stepping_manager

template<int dim>
TimeStepping::Manager<dim> aspect::Simulator< dim >::time_stepping_manager
private

Definition at line 1975 of file simulator.h.

§ lateral_averaging

template<int dim>
LateralAveraging<dim> aspect::Simulator< dim >::lateral_averaging
private

Definition at line 1984 of file simulator.h.

§ triangulation

template<int dim>
parallel::distributed::Triangulation<dim> aspect::Simulator< dim >::triangulation
private

Definition at line 1993 of file simulator.h.

§ global_Omega_diameter

template<int dim>
double aspect::Simulator< dim >::global_Omega_diameter
private

Definition at line 1994 of file simulator.h.

§ global_volume

template<int dim>
double aspect::Simulator< dim >::global_volume
private

Definition at line 1995 of file simulator.h.

§ mesh_refinement_manager

template<int dim>
MeshRefinement::Manager<dim> aspect::Simulator< dim >::mesh_refinement_manager
private

Definition at line 1997 of file simulator.h.

§ heating_model_manager

template<int dim>
HeatingModel::Manager<dim> aspect::Simulator< dim >::heating_model_manager
private

Definition at line 1998 of file simulator.h.

§ mapping

template<int dim>
std::unique_ptr<Mapping<dim> > aspect::Simulator< dim >::mapping
private

Pointer to the Mapping object used by the finite elements when going from the reference cell to the cell in the computational domain. We use a pointer since different mapping objects may be useful. In particular, when the mesh is deformable we use a MappingQ1Eulerian object to describe the mesh deformation, swapping it in for the original MappingQ or MappingCartesian object.

Definition at line 2008 of file simulator.h.

§ finite_element

template<int dim>
const FESystem<dim> aspect::Simulator< dim >::finite_element
private

Definition at line 2010 of file simulator.h.

§ dof_handler

template<int dim>
DoFHandler<dim> aspect::Simulator< dim >::dof_handler
private

Definition at line 2012 of file simulator.h.

§ postprocess_manager

template<int dim>
Postprocess::Manager<dim> aspect::Simulator< dim >::postprocess_manager
private

Definition at line 2014 of file simulator.h.

§ constraints

template<int dim>
AffineConstraints<double> aspect::Simulator< dim >::constraints
private

Constraint objects. The first of these describes all constraints that are not time dependent (e.g., hanging nodes, no-normal-flux constraints), whereas the second one is initialized at the top of every time step by copying from the first and then adding to it constraints that are time dependent (e.g., time dependent velocity or temperature boundary conditions).

'constraints' is computed in setup_dofs(), 'current_constraints' is done in compute_current_constraints().

Definition at line 2027 of file simulator.h.

§ current_constraints

template<int dim>
AffineConstraints<double> aspect::Simulator< dim >::current_constraints
private

Definition at line 2028 of file simulator.h.

§ last_pressure_normalization_adjustment

template<int dim>
double aspect::Simulator< dim >::last_pressure_normalization_adjustment
private

A place to store the latest correction computed by normalize_pressure(). We store this so we can undo the correction in denormalize_pressure().

Definition at line 2034 of file simulator.h.

§ pressure_scaling

template<int dim>
double aspect::Simulator< dim >::pressure_scaling
private

Scaling factor for the pressure as explained in the Kronbichler/Heister/Bangerth paper to ensure that the linear system that results from the Stokes equations is well conditioned.

Definition at line 2041 of file simulator.h.

§ do_pressure_rhs_compatibility_modification

template<int dim>
bool aspect::Simulator< dim >::do_pressure_rhs_compatibility_modification
private

A variable that determines whether we need to do the correction of the Stokes right hand side vector to ensure that the average divergence is zero. This is necessary for compressible models, but only if there are no in/outflow boundaries.

Definition at line 2049 of file simulator.h.

§ system_matrix

template<int dim>
LinearAlgebra::BlockSparseMatrix aspect::Simulator< dim >::system_matrix
private

An object that contains the entries of the system matrix. It has a size equal to the total number of degrees of freedom, but since we typically do not solve for all variables at once, the content of the matrix at any given time is only appropriate for the part of the system we are currently solving.

Definition at line 2069 of file simulator.h.

§ system_preconditioner_matrix

template<int dim>
LinearAlgebra::BlockSparseMatrix aspect::Simulator< dim >::system_preconditioner_matrix
private

An object that contains the entries of preconditioner matrices for the system matrix. It has a size equal to the total number of degrees of freedom, but is only used for the Stokes system (that's the only part of the system where we use a matrix for preconditioning that is different from the matrix we solve). Consequently, the blocks in rows and columns corresponding to temperature or compositional fields are left empty when building the sparsity pattern of this matrix in the Simulator::setup_system_preconditioner() function.

Definition at line 2083 of file simulator.h.

§ solution

template<int dim>
LinearAlgebra::BlockVector aspect::Simulator< dim >::solution
private

Definition at line 2085 of file simulator.h.

§ old_solution

template<int dim>
LinearAlgebra::BlockVector aspect::Simulator< dim >::old_solution
private

Definition at line 2086 of file simulator.h.

§ old_old_solution

template<int dim>
LinearAlgebra::BlockVector aspect::Simulator< dim >::old_old_solution
private

Definition at line 2087 of file simulator.h.

§ system_rhs

template<int dim>
LinearAlgebra::BlockVector aspect::Simulator< dim >::system_rhs
private

Definition at line 2088 of file simulator.h.

§ current_linearization_point

template<int dim>
LinearAlgebra::BlockVector aspect::Simulator< dim >::current_linearization_point
private

Definition at line 2090 of file simulator.h.

§ pressure_shape_function_integrals

template<int dim>
LinearAlgebra::BlockVector aspect::Simulator< dim >::pressure_shape_function_integrals
private

Definition at line 2093 of file simulator.h.

§ operator_split_reaction_vector

template<int dim>
LinearAlgebra::BlockVector aspect::Simulator< dim >::operator_split_reaction_vector
private

Definition at line 2096 of file simulator.h.

§ Amg_preconditioner

template<int dim>
std::unique_ptr<LinearAlgebra::PreconditionAMG> aspect::Simulator< dim >::Amg_preconditioner
private

Definition at line 2100 of file simulator.h.

§ Mp_preconditioner

template<int dim>
std::unique_ptr<LinearAlgebra::PreconditionBase> aspect::Simulator< dim >::Mp_preconditioner
private

Definition at line 2101 of file simulator.h.

§ rebuild_sparsity_and_matrices

template<int dim>
bool aspect::Simulator< dim >::rebuild_sparsity_and_matrices
private

Definition at line 2103 of file simulator.h.

§ rebuild_stokes_matrix

template<int dim>
bool aspect::Simulator< dim >::rebuild_stokes_matrix
private

Definition at line 2104 of file simulator.h.

§ assemble_newton_stokes_matrix

template<int dim>
bool aspect::Simulator< dim >::assemble_newton_stokes_matrix
private

Definition at line 2105 of file simulator.h.

§ assemble_newton_stokes_system

template<int dim>
bool aspect::Simulator< dim >::assemble_newton_stokes_system
private

Definition at line 2106 of file simulator.h.

§ rebuild_stokes_preconditioner

template<int dim>
bool aspect::Simulator< dim >::rebuild_stokes_preconditioner
private

Definition at line 2107 of file simulator.h.

§ mesh_deformation

template<int dim>
std::unique_ptr<MeshDeformation::MeshDeformationHandler<dim> > aspect::Simulator< dim >::mesh_deformation
private

Unique pointer for an instance of the MeshDeformationHandler. this way, if we do not need the machinery for doing mesh deformation stuff, we do not even allocate it.

Definition at line 2120 of file simulator.h.

§ stokes_matrix_free

template<int dim>
std::unique_ptr<StokesMatrixFreeHandler<dim> > aspect::Simulator< dim >::stokes_matrix_free
private

Unique pointer for the matrix-free Stokes solver

Definition at line 2125 of file simulator.h.


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