Amandus: Simulations based on multilevel Schwarz methods
|
#include <amandus.h>
Public Types | |
typedef dealii::MeshWorker::IntegrationInfo< dim > | CellInfo |
Public Types inherited from AmandusApplicationSparse< dim > | |
typedef dealii::MeshWorker::IntegrationInfo< dim > | CellInfo |
Public Member Functions | |
AmandusApplication (dealii::Triangulation< dim > &triangulation, const dealii::FiniteElement< dim > &fe) | |
virtual void | parse_parameters (dealii::ParameterHandler ¶m) |
void | set_advection_directions (std::vector< dealii::Tensor< 1, dim >> &adv_dir) |
void | setup_system () |
void | setup_constraints () |
void | assemble_mg_matrix (const dealii::AnyData &in, const AmandusIntegrator< dim > &integrator) |
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) |
Public Member Functions inherited from AmandusApplicationSparse< dim > | |
AmandusApplicationSparse (dealii::Triangulation< dim > &triangulation, const dealii::FiniteElement< dim > &fe, bool use_umfpack=false) | |
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 |
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... | |
void | assemble_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) |
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 | |
bool | vertex_patches = true |
bool | boundary_patches = false |
Public Attributes inherited from AmandusApplicationSparse< dim > | |
dealii::ReductionControl | control |
dealii::SmartPointer< dealii::ParameterHandler > | param |
dealii::Triangulation< dim >::Signals & | signals |
Protected Attributes | |
dealii::MGConstrainedDoFs | mg_constraints |
dealii::MGLevelObject< dealii::SparsityPattern > | mg_sparsity |
dealii::MGLevelObject< dealii::SparsityPattern > | mg_sparsity_fluxes |
dealii::MGLevelObject< dealii::SparseMatrix< double > > | mg_matrix |
dealii::MGLevelObject< dealii::SparseMatrix< double > > | mg_matrix_down |
dealii::MGLevelObject< dealii::SparseMatrix< double > > | mg_matrix_up |
dealii::MGLevelObject< dealii::SparseMatrix< double > > | mg_matrix_flux_down |
dealii::MGLevelObject< dealii::SparseMatrix< double > > | mg_matrix_flux_up |
dealii::MGTransferPrebuilt< dealii::Vector< double > > | mg_transfer |
dealii::FullMatrix< double > | coarse_matrix |
dealii::MGCoarseGridSVD< double, dealii::Vector< double > > | mg_coarse |
dealii::MGLevelObject< typename RELAXATION::AdditionalData > | smoother_data |
dealii::mg::SmootherRelaxation< RELAXATION, dealii::Vector< double > > | mg_smoother |
std::vector< dealii::Tensor< 1, dim > > | advection_directions |
bool | log_smoother_statistics = false |
call dealii::PreconditionBlockBase::log_statistics on each level More... | |
bool | right_preconditioning = true |
refers to dealii::SolverGMRES::AdditionalData::right_preconditioning More... | |
bool | use_default_residual = true |
refers to dealii::SolverGMRES::AdditionalData::use_default_residual More... | |
double | smoother_relaxation = 1.0 |
relaxation parameter for multigrid smoother More... | |
Protected Attributes inherited from AmandusApplicationSparse< dim > | |
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 solvers based on local Schwarz smoothers without strong boundary conditions.
This class provides storage for the DoFHandler object and the matrices associated with a finite element discretization and its multilevel solver. 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 (first two are inherited from base class):
typedef dealii::MeshWorker::IntegrationInfo<dim> AmandusApplication< dim, RELAXATION >::CellInfo |
AmandusApplication< dim, RELAXATION >::AmandusApplication | ( | dealii::Triangulation< dim > & | triangulation, |
const dealii::FiniteElement< dim > & | fe | ||
) |
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().
|
virtual |
Solve the eigenvalue system stored in the first element of matrix with mass matrix in the second. Use iterative inverse of the first matrix for shifting.
Reimplemented from AmandusApplicationSparse< dim >.
|
virtual |
Use the integrator to build the matrix for the level meshes. This also automatically generates the transfer matrices needed for multigrid with local smoothing on locally refined meshes.
Reimplemented from AmandusApplicationSparse< dim >.
|
virtual |
Parse the paramaters from a handler
Reimplemented from AmandusApplicationSparse< dim >.
|
inline |
This function is used to set the advection directions in which sorting of vertex patches or cell patches should be done. This is useful while applying Gauss-Seidel type smoothers in MG when parameter is true. is vector of directions so that you can sort in more than one direction.
|
virtual |
Set up hanging node constraints for leaf mesh and for level meshes. Use redefinition in derived classes to add boundary constraints.
Reimplemented from AmandusApplicationSparse< dim >.
|
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 from AmandusApplicationSparse< dim >.
|
virtual |
Solve the linear system stored in matrix with the right hand side given. Uses the multigrid preconditioner.
Reimplemented from AmandusApplicationSparse< dim >.
|
protected |
used for gauss seidel like smoothers for determining the order in which the subproblems are solved (if parameter Sort is true)
bool AmandusApplication< dim, RELAXATION >::boundary_patches = false |
If using vertex patches, also add patches for vertices at the boundary. This is necessary for instance if only interior degrees of freedom of the patch are used for local smoothing operations and some boundary degrees of freedom are not constrained by essential boundary conditions.
|
protected |
|
protected |
call dealii::PreconditionBlockBase::log_statistics on each level
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
refers to dealii::SolverGMRES::AdditionalData::right_preconditioning
|
protected |
|
protected |
relaxation parameter for multigrid smoother
|
protected |
refers to dealii::SolverGMRES::AdditionalData::use_default_residual
bool AmandusApplication< dim, RELAXATION >::vertex_patches = true |
The multigrid smoother uses vertex patches instead of inverting cell matrices.