Reference documentation for deal.II version 9.1.1
|
#include <deal.II/lac/petsc_solver.h>
Classes | |
struct | AdditionalData |
struct | SolverDataMUMPS |
Public Member Functions | |
SparseDirectMUMPS (SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData()) | |
void | solve (const MatrixBase &A, VectorBase &x, const VectorBase &b) |
void | set_symmetric_mode (const bool flag) |
Public Member Functions inherited from PETScWrappers::SolverBase | |
SolverBase (SolverControl &cn, const MPI_Comm &mpi_communicator) | |
virtual | ~SolverBase ()=default |
void | solve (const MatrixBase &A, VectorBase &x, const VectorBase &b, const PreconditionerBase &preconditioner) |
virtual void | reset () |
void | set_prefix (const std::string &prefix) |
SolverControl & | control () const |
void | initialize (const PreconditionerBase &preconditioner) |
Protected Member Functions | |
virtual void | set_solver_type (KSP &ksp) const override |
Protected Attributes | |
const AdditionalData | additional_data |
Protected Attributes inherited from PETScWrappers::SolverBase | |
SolverControl & | solver_control |
const MPI_Comm | mpi_communicator |
std::string | prefix_name |
Static Private Member Functions | |
static PetscErrorCode | convergence_test (KSP ksp, const PetscInt iteration, const PetscReal residual_norm, KSPConvergedReason *reason, void *solver_control) |
Private Attributes | |
bool | symmetric_mode |
An implementation of the solver interface using the sparse direct MUMPS solver through PETSc. This class has the usual interface of all other solver classes but it is of course different in that it doesn't implement an iterative solver. As a consequence, things like the SolverControl object have no particular meaning here.
MUMPS allows to make use of symmetry in this matrix. In this class this is made possible by the set_symmetric_mode() function. If your matrix is symmetric, you can use this class as follows:
Definition at line 934 of file petsc_solver.h.
PETScWrappers::SparseDirectMUMPS::SparseDirectMUMPS | ( | SolverControl & | cn, |
const MPI_Comm & | mpi_communicator = PETSC_COMM_SELF , |
||
const AdditionalData & | data = AdditionalData() |
||
) |
Constructor
Definition at line 634 of file petsc_solver.cc.
void PETScWrappers::SparseDirectMUMPS::solve | ( | const MatrixBase & | A, |
VectorBase & | x, | ||
const VectorBase & | b | ||
) |
The method to solve the linear system.
factorization matrix to be obtained from MUMPS
setting MUMPS integer control parameters ICNTL to be passed to MUMPS. Setting entry 7 of MUMPS ICNTL array (of size 40) to a value of 2. This sets use of Approximate Minimum Fill (AMF)
number of iterations to solution (should be 1) for a direct solver
norm of residual
creating a solver object if this is necessary
creates the default KSP context and puts it in the location solver_data->ksp
set the matrices involved. the last argument is irrelevant here, since we use the solver only once anyway
setting the solver type
getting the associated preconditioner context
build PETSc PC for particular PCLU or PCCHOLESKY preconditioner depending on whether the symmetric mode has been set
convergence monitor function that checks with the solver_control object for convergence
set the software that is to be used to perform the lu factorization here we start to see differences with the base class solve function
set up the package to call for the factorization
get the factored matrix F from the preconditioner context. This routine is valid only for LU, ILU, Cholesky, and incomplete Cholesky
Passing the control parameters to MUMPS
set the command line option prefix name
set the command line options provided by the user to override the defaults
solve the linear system
in case of failure throw exception
obtain convergence information. obtain the number of iterations and residual norm
Definition at line 671 of file petsc_solver.cc.
void PETScWrappers::SparseDirectMUMPS::set_symmetric_mode | ( | const bool | flag | ) |
The method allows to take advantage if the system matrix is symmetric by using LDL^T decomposition instead of more expensive LU. The argument indicates whether the matrix is symmetric or not.
Definition at line 884 of file petsc_solver.cc.
|
overrideprotectedvirtual |
Function that takes a Krylov Subspace Solver context object, and sets the type of solver that is requested by the derived class.
KSPPREONLY implements a stub method that applies only the preconditioner. Its use is due to SparseDirectMUMPS being a direct (rather than iterative) solver
The KSPPREONLY solver of PETSc never calls the convergence monitor, which leads to failure even when everything was ok. Therefore, the SolverControl status is set to some nice values, which guarantee a nice result at the end of the solution process.
Using a PREONLY solver with a nonzero initial guess leads PETSc to produce some error messages.
Implements PETScWrappers::SolverBase.
Definition at line 644 of file petsc_solver.cc.
|
staticprivate |
A function that is used in PETSc as a callback to check convergence. It takes the information provided from PETSc and checks it against deal.II's own SolverControl objects to see if convergence has been reached.
Definition at line 847 of file petsc_solver.cc.
|
protected |
Store a copy of flags for this particular solver.
Definition at line 967 of file petsc_solver.h.
|
private |
Flag specifies whether matrix being factorized is symmetric or not. It influences the type of the used preconditioner (PCLU or PCCHOLESKY)
Definition at line 1008 of file petsc_solver.h.