ASPECT
Public Member Functions | Static Public Member Functions | List of all members
aspect::StokesMatrixFreeHandler< dim > Class Template Referenceabstract
Inheritance diagram for aspect::StokesMatrixFreeHandler< dim >:
Inheritance graph
[legend]

Public Member Functions

virtual ~StokesMatrixFreeHandler ()=default
 
virtual std::pair< double, double > solve (LinearAlgebra::BlockVector &solution_vector)=0
 
virtual void setup_dofs ()=0
 
virtual void assemble ()=0
 
virtual void build_preconditioner ()=0
 
virtual const DoFHandler< dim > & get_dof_handler_v () const =0
 
virtual const DoFHandler< dim > & get_dof_handler_p () const =0
 
virtual const DoFHandler< dim > & get_dof_handler_projection () const =0
 
virtual const AffineConstraints< double > & get_constraints_v () const =0
 
virtual const AffineConstraints< double > & get_constraints_p () const =0
 
virtual const MGTransferMF< dim, GMGNumberType > & get_mg_transfer_A () const =0
 
virtual const MGTransferMF< dim, GMGNumberType > & get_mg_transfer_S () const =0
 
virtual std::size_t get_cell_data_memory_consumption () const =0
 

Static Public Member Functions

static void declare_parameters (ParameterHandler &prm)
 

Detailed Description

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

Base class for the matrix free GMG solver for the Stokes system. The actual implementation is found inside StokesMatrixFreeHandlerImplementation below.

Definition at line 87 of file simulator.h.

Constructor & Destructor Documentation

§ ~StokesMatrixFreeHandler()

template<int dim>
virtual aspect::StokesMatrixFreeHandler< dim >::~StokesMatrixFreeHandler ( )
virtualdefault

virtual Destructor.

Member Function Documentation

§ solve()

template<int dim>
virtual std::pair<double,double> aspect::StokesMatrixFreeHandler< dim >::solve ( LinearAlgebra::BlockVector solution_vector)
pure virtual

Solves the Stokes linear system using the matrix-free solver.

Parameters
solution_vectorThe existing solution vector that will be updated with the new solution. This vector is expected to have the block structure of the full solution vector, and its velocity and pressure blocks will be updated with the new solution.

Implemented in aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >.

§ setup_dofs()

template<int dim>
virtual void aspect::StokesMatrixFreeHandler< dim >::setup_dofs ( )
pure virtual

Allocates and sets up the members of the StokesMatrixFreeHandler. This is called by Simulator<dim>::setup_dofs()

Implemented in aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >.

§ assemble()

template<int dim>
virtual void aspect::StokesMatrixFreeHandler< dim >::assemble ( )
pure virtual

Perform various tasks to update the linear system to solve for. Note that we are not assembling a matrix (as this is a matrix-free algorithm), but we are evaluating the material model and storing the information necessary for a later call to solve().

Implemented in aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >.

§ build_preconditioner()

template<int dim>
virtual void aspect::StokesMatrixFreeHandler< dim >::build_preconditioner ( )
pure virtual

Computes and sets the diagonal for both the mass matrix operator and the A-block operators on each level for the purpose of smoothing inside the multigrid v-cycle.

Implemented in aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >.

§ declare_parameters()

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

Declare parameters.

§ get_dof_handler_v()

template<int dim>
virtual const DoFHandler<dim>& aspect::StokesMatrixFreeHandler< dim >::get_dof_handler_v ( ) const
pure virtual

Return a reference to the DoFHandler that is used for velocity in the block GMG solver.

Implemented in aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >.

§ get_dof_handler_p()

template<int dim>
virtual const DoFHandler<dim>& aspect::StokesMatrixFreeHandler< dim >::get_dof_handler_p ( ) const
pure virtual

Return a reference to the DoFHandler that is used for pressure in the block GMG solver.

Implemented in aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >.

§ get_dof_handler_projection()

template<int dim>
virtual const DoFHandler<dim>& aspect::StokesMatrixFreeHandler< dim >::get_dof_handler_projection ( ) const
pure virtual

Return a reference to the DoFHandler that is used for the coefficient projection in the block GMG solver.

Implemented in aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >.

§ get_constraints_v()

template<int dim>
virtual const AffineConstraints<double>& aspect::StokesMatrixFreeHandler< dim >::get_constraints_v ( ) const
pure virtual

Return a pointer to the object that describes the velocity DoF constraints for the block GMG Stokes solver.

Implemented in aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >.

§ get_constraints_p()

template<int dim>
virtual const AffineConstraints<double>& aspect::StokesMatrixFreeHandler< dim >::get_constraints_p ( ) const
pure virtual

Return a pointer to the object that describes the pressure DoF constraints for the block GMG Stokes solver.

Implemented in aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >.

§ get_mg_transfer_A()

template<int dim>
virtual const MGTransferMF<dim,GMGNumberType>& aspect::StokesMatrixFreeHandler< dim >::get_mg_transfer_A ( ) const
pure virtual

Return a pointer to the MGTransfer object used for the A block of the block GMG Stokes solver.

Implemented in aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >.

§ get_mg_transfer_S()

template<int dim>
virtual const MGTransferMF<dim,GMGNumberType>& aspect::StokesMatrixFreeHandler< dim >::get_mg_transfer_S ( ) const
pure virtual

Return a pointer to the MGTransfer object used for the Schur complement block of the block GMG Stokes solver.

Implemented in aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >.

§ get_cell_data_memory_consumption()

template<int dim>
virtual std::size_t aspect::StokesMatrixFreeHandler< dim >::get_cell_data_memory_consumption ( ) const
pure virtual

Return the memory consumption in bytes that are used to store equation data like viscosity to be able to apply the operators.

Implemented in aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >.


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