Amandus: Simulations based on multilevel Schwarz methods
|
#include <amandus.h>
Public Types | |
typedef dealii::MeshWorker::IntegrationInfo< dim > | CellInfo |
Public Member Functions | |
AmandusApplicationSparse (dealii::Triangulation< dim > &triangulation, const dealii::FiniteElement< dim > &fe, bool use_umfpack=false) | |
virtual void | parse_parameters (dealii::ParameterHandler ¶m) |
void | set_number_of_matrices (unsigned int n) |
void | set_boundary (dealii::types::boundary_id index, dealii::ComponentMask mask=dealii::ComponentMask()) |
void | set_meanvalue (dealii::ComponentMask mask=dealii::ComponentMask()) |
virtual void | setup_vector (dealii::Vector< double > &v) const |
virtual void | update_vector_inhom_boundary (dealii::Vector< double > &v, const dealii::Function< dim > &inhom_boundary, bool projection=false) const |
virtual void | setup_system () |
void | assemble_right_hand_side (dealii::AnyData &out, const dealii::AnyData &in, const AmandusIntegrator< dim > &integrator) const |
void | refine_mesh (const bool global=false) |
const dealii::DoFHandler< dim > & | dofs () const |
The object describing the finite element space. More... | |
const dealii::ConstraintMatrix & | constraints () const |
The object describing the constraints. More... | |
const dealii::ConstraintMatrix & | hanging_nodes () const |
The object describing the constraints for hanging nodes, not for the boundary. More... | |
virtual void | setup_constraints () |
void | assemble_matrix (const dealii::AnyData &in, const AmandusIntegrator< dim > &integrator) |
virtual void | assemble_mg_matrix (const dealii::AnyData &in, const AmandusIntegrator< dim > &integrator) |
const dealii::Vector< double > & | indicators () const |
The error indicators. More... | |
double | estimate (const dealii::AnyData &in, AmandusIntegrator< dim > &integrator) |
void | error (dealii::BlockVector< double > &out, const dealii::AnyData &in, const AmandusIntegrator< dim > &integrator) |
void | error (dealii::BlockVector< double > &out, const dealii::AnyData &in, const ErrorIntegrator< dim > &integrator) |
void | error (const dealii::AnyData &in, const AmandusIntegrator< dim > &integrator, unsigned int num_errs) |
virtual void | solve (dealii::Vector< double > &sol, const dealii::Vector< double > &rhs) |
virtual void | arpack_solve (std::vector< std::complex< double >> &eigenvalues, std::vector< dealii::Vector< double >> &eigenvectors) |
void | output_results (unsigned int refinement_cycle, const dealii::AnyData *data=0) const |
void | verify_residual (dealii::AnyData &out, const dealii::AnyData &in, const AmandusIntegrator< dim > &integrator) const |
Public Attributes | |
dealii::ReductionControl | control |
dealii::SmartPointer< dealii::ParameterHandler > | param |
dealii::Triangulation< dim >::Signals & | signals |
Protected Attributes | |
dealii::SmartPointer< dealii::Triangulation< dim >, AmandusApplicationSparse< dim > > | triangulation |
The mesh. More... | |
const dealii::MappingQ1< dim > | mapping |
The default mapping. More... | |
dealii::SmartPointer< const dealii::FiniteElement< dim >, AmandusApplicationSparse< dim > > | fe |
The finite element constructed from the string. More... | |
dealii::DoFHandler< dim > | dof_handler |
The object handling the degrees of freedom. More... | |
std::map< dealii::types::boundary_id, dealii::ComponentMask > | boundary_masks |
The masks used to set boundary conditions, indexed by the boundary indicator. More... | |
dealii::ComponentMask | meanvalue_mask |
dealii::ConstraintMatrix | constraint_matrix |
The object holding the constraints for the active mesh. More... | |
dealii::ConstraintMatrix | hanging_node_constraints |
The object holding the hanging node constraints for the active mesh. More... | |
dealii::SparsityPattern | sparsity |
std::vector< dealii::SparseMatrix< double > > | matrix |
const bool | use_umfpack |
dealii::SparseDirectUMFPACK | inverse |
dealii::BlockVector< double > | estimates |
std::vector< dealii::DataComponentInterpretation::DataComponentInterpretation > | output_data_types |
An application class with a plain GMRES solver and optional UMFPack as a preconditioner (see constructor arguments). This is mostly to test discretizations and it sould be improved by a derived class with a multigrid solver like AmandusApplication.
This class provides storage for the DoFHandler object and the matrix associated with a finite element discretization. The basic structures only depend on the Triangulation and the FiniteElement, which are provided to the constructor.
The purpose of this class is not so much implementing an application program, but providing the data structures used by all application programs with the following characteristics:
typedef dealii::MeshWorker::IntegrationInfo<dim> AmandusApplicationSparse< dim >::CellInfo |
AmandusApplicationSparse< dim >::AmandusApplicationSparse | ( | dealii::Triangulation< dim > & | triangulation, |
const dealii::FiniteElement< dim > & | fe, | ||
bool | use_umfpack = false |
||
) |
Constructor, setting the finite element and the triangulation. This constructor does not distribute the degrees of freedom or initialize sparsity patterns, which has to be achieved by calling setup_system().
triangulation | The mesh used to build the application |
fe | The finite element space used for discretization |
use_umfpack | if true implies that assemble_matrix() not only generates a matrix, but also the inverse, using SparseDirectUMFPACK. |
|
virtual |
Solve the eigenvalue system stored in the first element of matrix with mass matrix in the second. Use the UMFPack inverse of the first matrix aslo for shifting.
Reimplemented in AmandusApplication< dim, RELAXATION >.
void AmandusApplicationSparse< dim >::assemble_matrix | ( | const dealii::AnyData & | in, |
const AmandusIntegrator< dim > & | integrator | ||
) |
Use the integrator to build the matrix for the leaf mesh.
|
virtual |
Empty function here, but it is reimplemented in AmandusApplicationMultigrid.
Reimplemented in AmandusApplication< dim, RELAXATION >.
void AmandusApplicationSparse< dim >::assemble_right_hand_side | ( | dealii::AnyData & | out, |
const dealii::AnyData & | in, | ||
const AmandusIntegrator< dim > & | integrator | ||
) | const |
Apply the local operator integrator
to the vectors in in
and write the result to the first vector in out
.
|
inline |
The object describing the constraints.
|
inline |
The object describing the finite element space.
void AmandusApplicationSparse< dim >::error | ( | dealii::BlockVector< double > & | out, |
const dealii::AnyData & | in, | ||
const AmandusIntegrator< dim > & | integrator | ||
) |
Compute several error values using the integrator and return them in a BlockVector.
The number of errors computed is the number of blocks in the vector. The size of each block is adjusted inside this function to match the number of cells. Then, the error contribution of each cell is stored in this vector.
void AmandusApplicationSparse< dim >::error | ( | dealii::BlockVector< double > & | out, |
const dealii::AnyData & | in, | ||
const ErrorIntegrator< dim > & | integrator | ||
) |
Compute several error values using the error integrator and return them in a BlockVector.
The BlockVector will be resized to match the number of errors calulcated by the integrator.
void AmandusApplicationSparse< dim >::error | ( | const dealii::AnyData & | in, |
const AmandusIntegrator< dim > & | integrator, | ||
unsigned int | num_errs | ||
) |
Compute errors and print them to #deallog
double AmandusApplicationSparse< dim >::estimate | ( | const dealii::AnyData & | in, |
AmandusIntegrator< dim > & | integrator | ||
) |
Currently disabled.
|
inline |
The object describing the constraints for hanging nodes, not for the boundary.
|
inline |
The error indicators.
void AmandusApplicationSparse< dim >::output_results | ( | unsigned int | refinement_cycle, |
const dealii::AnyData * | data = 0 |
||
) | const |
|
virtual |
Parse the paramaters from a handler
Reimplemented in AmandusApplication< dim, RELAXATION >.
void AmandusApplicationSparse< dim >::refine_mesh | ( | const bool | global = false | ) |
Refine the mesh globally. For more sophisticated (adaptive) refinement strategies, use a Remesher.
void AmandusApplicationSparse< dim >::set_boundary | ( | dealii::types::boundary_id | index, |
dealii::ComponentMask | mask = dealii::ComponentMask() |
||
) |
Set the boundary components that should be constrained if the boundary indicator is equal to index.
The constraints for these boundary values are set in setup_constraints(), and they are always constrained to zero. Inhomogeneous boundary conditions are obtained by setting the boundary values of the start vector of a dealii::Newton or a dealii::ThetaTimestepping method.
index | the boundary indicator for which the boundary constraints re being set. |
mask | the object selecting the blocks of an dealii::FESystem to which the constraints are to be applied. |
void AmandusApplicationSparse< dim >::set_meanvalue | ( | dealii::ComponentMask | mask = dealii::ComponentMask() | ) |
Constrain solution to be mean value free.
mask | the object selecting the blocks of an dealii::FESystem to which the constraints are to be applied. |
|
inline |
Change the number of matrices assembled. Default is one, but for instance a second matrix (the mass matrix) is needed for eigenvelue problems.
|
virtual |
Set up hanging node constraints for leaf mesh and for level meshes. Use redefinition in derived classes to add boundary constraints.
Reimplemented in AmandusApplication< dim, RELAXATION >.
|
virtual |
Initialize the finite element system on the current mesh. This involves distributing degrees of freedom for the leaf mesh and the levels as well as setting up sparsity patterns for the matrices.
Reimplemented in AmandusApplication< dim, RELAXATION >.
|
virtual |
Initialize the vector v
to the size matching the DoFHandler. This requires that setup_system() is called before.
|
virtual |
Solve the linear system stored in matrix with the right hand side given. Uses the multigrid preconditioner.
Reimplemented in AmandusApplication< dim, RELAXATION >.
|
virtual |
Sets degrees on the boundary to their inhomogeneous Dirichlet constraints with an option of interpolation or projection. This requires that setup_system() and setup_vector() are called before.
void AmandusApplicationSparse< dim >::verify_residual | ( | dealii::AnyData & | out, |
const dealii::AnyData & | in, | ||
const AmandusIntegrator< dim > & | integrator | ||
) | const |
For testing purposes compare the residual operator applied to to the vectors in in
with the result of the matrix vector multiplication.
|
protected |
The masks used to set boundary conditions, indexed by the boundary indicator.
|
protected |
The object holding the constraints for the active mesh.
dealii::ReductionControl AmandusApplicationSparse< dim >::control |
The object controlling the iteration in solve().
|
protected |
The object handling the degrees of freedom.
|
protected |
|
protected |
The finite element constructed from the string.
|
protected |
The object holding the hanging node constraints for the active mesh.
|
protected |
|
protected |
The default mapping.
|
protected |
|
protected |
|
protected |
dealii::SmartPointer<dealii::ParameterHandler> AmandusApplicationSparse< dim >::param |
Reference to parameters read by parse_parameters().
dealii::Triangulation<dim>::Signals& AmandusApplicationSparse< dim >::signals |
|
protected |
|
protected |
The mesh.
|
protected |