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

Public Member Functions

 StokesMatrixFreeHandlerImplementation (Simulator< dim > &, ParameterHandler &prm)
 
 ~StokesMatrixFreeHandlerImplementation () override=default
 
std::pair< double, double > solve (LinearAlgebra::BlockVector &solution_vector) override
 
void setup_dofs () override
 
void assemble () override
 
void build_preconditioner () override
 
const DoFHandler< dim > & get_dof_handler_v () const override
 
const DoFHandler< dim > & get_dof_handler_p () const override
 
const DoFHandler< dim > & get_dof_handler_projection () const override
 
const AffineConstraints< double > & get_constraints_v () const override
 
const AffineConstraints< double > & get_constraints_p () const override
 
const MGTransferMF< dim, GMGNumberType > & get_mg_transfer_A () const override
 
const MGTransferMF< dim, GMGNumberType > & get_mg_transfer_S () const override
 
std::size_t get_cell_data_memory_consumption () const override
 
- Public Member Functions inherited from aspect::StokesMatrixFreeHandler< dim >
virtual ~StokesMatrixFreeHandler ()=default
 

Static Public Member Functions

static void declare_parameters (ParameterHandler &prm)
 
- Static Public Member Functions inherited from aspect::StokesMatrixFreeHandler< dim >
static void declare_parameters (ParameterHandler &prm)
 

Private Types

using StokesMatrixType = MatrixFreeStokesOperators::StokesOperator< dim, velocity_degree, double >
 
using SchurComplementMatrixType = MatrixFreeStokesOperators::MassMatrixOperator< dim, velocity_degree-1, double >
 
using ABlockMatrixType = MatrixFreeStokesOperators::ABlockOperator< dim, velocity_degree, double >
 
using GMGSchurComplementMatrixType = MatrixFreeStokesOperators::MassMatrixOperator< dim, velocity_degree-1, GMGNumberType >
 
using GMGABlockMatrixType = MatrixFreeStokesOperators::ABlockOperator< dim, velocity_degree, GMGNumberType >
 

Private Member Functions

void parse_parameters (ParameterHandler &prm)
 
void evaluate_material_model ()
 
void correct_stokes_rhs ()
 

Private Attributes

Simulator< dim > & sim
 
bool print_details
 
bool do_timings
 
double minimum_viscosity
 
double maximum_viscosity
 
DoFHandler< dim > dof_handler_v
 
DoFHandler< dim > dof_handler_p
 
DoFHandler< dim > dof_handler_projection
 
FESystem< dim > fe_v
 
FESystem< dim > fe_p
 
FESystem< dim > fe_projection
 
MatrixFreeStokesOperators::OperatorCellData< dim, GMGNumberTypeactive_cell_data
 
MGLevelObject< MatrixFreeStokesOperators::OperatorCellData< dim, GMGNumberType > > level_cell_data
 
StokesMatrixType stokes_matrix
 
ABlockMatrixType A_block_matrix
 
SchurComplementMatrixType Schur_complement_block_matrix
 
AffineConstraints< double > constraints_v
 
AffineConstraints< double > constraints_p
 
MGLevelObject< GMGABlockMatrixTypemg_matrices_A_block
 
MGLevelObject< GMGSchurComplementMatrixTypemg_matrices_Schur_complement
 
MGConstrainedDoFs mg_constrained_dofs_A_block
 
MGConstrainedDoFs mg_constrained_dofs_Schur_complement
 
MGConstrainedDoFs mg_constrained_dofs_projection
 
MGTransferMF< dim, GMGNumberTypemg_transfer_A_block
 
MGTransferMF< dim, GMGNumberTypemg_transfer_Schur_complement
 
std::vector< std::shared_ptr< MatrixFree< dim, double > > > matrix_free_objects
 

Detailed Description

template<int dim, int velocity_degree>
class aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >

Main class of the Matrix-free method. Here are all the functions for setup, assembly and solving the Stokes system.

We need to derive from StokesMatrixFreeHandler to be able to introduce a second template argument for the degree of the Stokes finite element. This way, the main simulator does not need to know about the degree by using a pointer to the base class and we can pick the desired velocity degree at runtime.

Definition at line 90 of file simulator.h.

Member Typedef Documentation

§ StokesMatrixType

template<int dim, int velocity_degree>
using aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::StokesMatrixType = MatrixFreeStokesOperators::StokesOperator<dim,velocity_degree,double>
private

Definition at line 683 of file stokes_matrix_free.h.

§ SchurComplementMatrixType

template<int dim, int velocity_degree>
using aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::SchurComplementMatrixType = MatrixFreeStokesOperators::MassMatrixOperator<dim,velocity_degree-1,double>
private

Definition at line 684 of file stokes_matrix_free.h.

§ ABlockMatrixType

template<int dim, int velocity_degree>
using aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::ABlockMatrixType = MatrixFreeStokesOperators::ABlockOperator<dim,velocity_degree,double>
private

Definition at line 685 of file stokes_matrix_free.h.

§ GMGSchurComplementMatrixType

template<int dim, int velocity_degree>
using aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::GMGSchurComplementMatrixType = MatrixFreeStokesOperators::MassMatrixOperator<dim,velocity_degree-1,GMGNumberType>
private

Definition at line 687 of file stokes_matrix_free.h.

§ GMGABlockMatrixType

template<int dim, int velocity_degree>
using aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::GMGABlockMatrixType = MatrixFreeStokesOperators::ABlockOperator<dim,velocity_degree,GMGNumberType>
private

Definition at line 688 of file stokes_matrix_free.h.

Constructor & Destructor Documentation

§ StokesMatrixFreeHandlerImplementation()

template<int dim, int velocity_degree>
aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::StokesMatrixFreeHandlerImplementation ( Simulator< dim > &  ,
ParameterHandler &  prm 
)

Initialize this class, allowing it to read in relevant parameters as well as giving it a reference to the Simulator that owns it, since it needs to make fairly extensive changes to the internals of the simulator.

§ ~StokesMatrixFreeHandlerImplementation()

template<int dim, int velocity_degree>
aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::~StokesMatrixFreeHandlerImplementation ( )
overridedefault

Destructor.

Member Function Documentation

§ solve()

template<int dim, int velocity_degree>
std::pair<double,double> aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::solve ( LinearAlgebra::BlockVector solution_vector)
overridevirtual

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.

Implements aspect::StokesMatrixFreeHandler< dim >.

§ setup_dofs()

template<int dim, int velocity_degree>
void aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::setup_dofs ( )
overridevirtual

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

Implements aspect::StokesMatrixFreeHandler< dim >.

§ assemble()

template<int dim, int velocity_degree>
void aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::assemble ( )
overridevirtual

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

Implements aspect::StokesMatrixFreeHandler< dim >.

§ build_preconditioner()

template<int dim, int velocity_degree>
void aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::build_preconditioner ( )
overridevirtual

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.

Implements aspect::StokesMatrixFreeHandler< dim >.

§ declare_parameters()

template<int dim, int velocity_degree>
static void aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::declare_parameters ( ParameterHandler &  prm)
static

Declare parameters. (No actual parameters at the moment).

§ get_dof_handler_v()

template<int dim, int velocity_degree>
const DoFHandler<dim>& aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::get_dof_handler_v ( ) const
overridevirtual

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

Implements aspect::StokesMatrixFreeHandler< dim >.

§ get_dof_handler_p()

template<int dim, int velocity_degree>
const DoFHandler<dim>& aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::get_dof_handler_p ( ) const
overridevirtual

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

Implements aspect::StokesMatrixFreeHandler< dim >.

§ get_dof_handler_projection()

template<int dim, int velocity_degree>
const DoFHandler<dim>& aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::get_dof_handler_projection ( ) const
overridevirtual

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

Implements aspect::StokesMatrixFreeHandler< dim >.

§ get_constraints_v()

template<int dim, int velocity_degree>
const AffineConstraints<double>& aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::get_constraints_v ( ) const
overridevirtual

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

Implements aspect::StokesMatrixFreeHandler< dim >.

§ get_constraints_p()

template<int dim, int velocity_degree>
const AffineConstraints<double>& aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::get_constraints_p ( ) const
overridevirtual

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

Implements aspect::StokesMatrixFreeHandler< dim >.

§ get_mg_transfer_A()

template<int dim, int velocity_degree>
const MGTransferMF<dim,GMGNumberType>& aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::get_mg_transfer_A ( ) const
overridevirtual

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

Implements aspect::StokesMatrixFreeHandler< dim >.

§ get_mg_transfer_S()

template<int dim, int velocity_degree>
const MGTransferMF<dim,GMGNumberType>& aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::get_mg_transfer_S ( ) const
overridevirtual

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

Implements aspect::StokesMatrixFreeHandler< dim >.

§ get_cell_data_memory_consumption()

template<int dim, int velocity_degree>
std::size_t aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::get_cell_data_memory_consumption ( ) const
overridevirtual

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

Implements aspect::StokesMatrixFreeHandler< dim >.

§ parse_parameters()

template<int dim, int velocity_degree>
void aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::parse_parameters ( ParameterHandler &  prm)
private

Parse parameters.

§ evaluate_material_model()

template<int dim, int velocity_degree>
void aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::evaluate_material_model ( )
private

Evaluate the MaterialModel to query information like the viscosity and project this viscosity to the multigrid hierarchy. Also queries other parameters like pressure scaling.

§ correct_stokes_rhs()

template<int dim, int velocity_degree>
void aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::correct_stokes_rhs ( )
private

Add correction to system RHS for non-zero boundary condition. See description in StokesMatrixFreeHandler::correct_stokes_rhs() for more information.

Member Data Documentation

§ sim

template<int dim, int velocity_degree>
Simulator<dim>& aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::sim
private

Definition at line 649 of file stokes_matrix_free.h.

§ print_details

template<int dim, int velocity_degree>
bool aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::print_details
private

Definition at line 651 of file stokes_matrix_free.h.

§ do_timings

template<int dim, int velocity_degree>
bool aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::do_timings
private

If true, it will time the key components of this matrix-free implementation, such as vmult of different matrices, solver IDR with the cheap preconditioner, etc.

Definition at line 657 of file stokes_matrix_free.h.

§ minimum_viscosity

template<int dim, int velocity_degree>
double aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::minimum_viscosity
private

The max/min of the evaluated viscosities.

Definition at line 662 of file stokes_matrix_free.h.

§ maximum_viscosity

template<int dim, int velocity_degree>
double aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::maximum_viscosity
private

Definition at line 663 of file stokes_matrix_free.h.

§ dof_handler_v

template<int dim, int velocity_degree>
DoFHandler<dim> aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::dof_handler_v
private

Definition at line 665 of file stokes_matrix_free.h.

§ dof_handler_p

template<int dim, int velocity_degree>
DoFHandler<dim> aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::dof_handler_p
private

Definition at line 666 of file stokes_matrix_free.h.

§ dof_handler_projection

template<int dim, int velocity_degree>
DoFHandler<dim> aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::dof_handler_projection
private

Definition at line 667 of file stokes_matrix_free.h.

§ fe_v

template<int dim, int velocity_degree>
FESystem<dim> aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::fe_v
private

Definition at line 669 of file stokes_matrix_free.h.

§ fe_p

template<int dim, int velocity_degree>
FESystem<dim> aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::fe_p
private

Definition at line 670 of file stokes_matrix_free.h.

§ fe_projection

template<int dim, int velocity_degree>
FESystem<dim> aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::fe_projection
private

Definition at line 671 of file stokes_matrix_free.h.

§ active_cell_data

template<int dim, int velocity_degree>
MatrixFreeStokesOperators::OperatorCellData<dim, GMGNumberType> aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::active_cell_data
private

Store the data for the Stokes operator (viscosity, etc.) for the active cells.

Definition at line 676 of file stokes_matrix_free.h.

§ level_cell_data

template<int dim, int velocity_degree>
MGLevelObject<MatrixFreeStokesOperators::OperatorCellData<dim, GMGNumberType> > aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::level_cell_data
private

Store the data for the Stokes operator (viscosity, etc.) for each multigrid level.

Definition at line 681 of file stokes_matrix_free.h.

§ stokes_matrix

template<int dim, int velocity_degree>
StokesMatrixType aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::stokes_matrix
private

Definition at line 690 of file stokes_matrix_free.h.

§ A_block_matrix

template<int dim, int velocity_degree>
ABlockMatrixType aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::A_block_matrix
private

Definition at line 691 of file stokes_matrix_free.h.

§ Schur_complement_block_matrix

template<int dim, int velocity_degree>
SchurComplementMatrixType aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::Schur_complement_block_matrix
private

Definition at line 692 of file stokes_matrix_free.h.

§ constraints_v

template<int dim, int velocity_degree>
AffineConstraints<double> aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::constraints_v
private

Definition at line 694 of file stokes_matrix_free.h.

§ constraints_p

template<int dim, int velocity_degree>
AffineConstraints<double> aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::constraints_p
private

Definition at line 695 of file stokes_matrix_free.h.

§ mg_matrices_A_block

template<int dim, int velocity_degree>
MGLevelObject<GMGABlockMatrixType> aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::mg_matrices_A_block
private

Definition at line 697 of file stokes_matrix_free.h.

§ mg_matrices_Schur_complement

template<int dim, int velocity_degree>
MGLevelObject<GMGSchurComplementMatrixType> aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::mg_matrices_Schur_complement
private

Definition at line 698 of file stokes_matrix_free.h.

§ mg_constrained_dofs_A_block

template<int dim, int velocity_degree>
MGConstrainedDoFs aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::mg_constrained_dofs_A_block
private

Definition at line 700 of file stokes_matrix_free.h.

§ mg_constrained_dofs_Schur_complement

template<int dim, int velocity_degree>
MGConstrainedDoFs aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::mg_constrained_dofs_Schur_complement
private

Definition at line 701 of file stokes_matrix_free.h.

§ mg_constrained_dofs_projection

template<int dim, int velocity_degree>
MGConstrainedDoFs aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::mg_constrained_dofs_projection
private

Definition at line 702 of file stokes_matrix_free.h.

§ mg_transfer_A_block

template<int dim, int velocity_degree>
MGTransferMF<dim,GMGNumberType> aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::mg_transfer_A_block
private

Definition at line 704 of file stokes_matrix_free.h.

§ mg_transfer_Schur_complement

template<int dim, int velocity_degree>
MGTransferMF<dim,GMGNumberType> aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::mg_transfer_Schur_complement
private

Definition at line 705 of file stokes_matrix_free.h.

§ matrix_free_objects

template<int dim, int velocity_degree>
std::vector<std::shared_ptr<MatrixFree<dim,double> > > aspect::StokesMatrixFreeHandlerImplementation< dim, velocity_degree >::matrix_free_objects
private

Definition at line 707 of file stokes_matrix_free.h.


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