![]() |
Reference documentation for deal.II version 8.1.0
|
#include <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) |
![]() | |
SolverBase (SolverControl &cn, const MPI_Comm &mpi_communicator) | |
virtual | ~SolverBase () |
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 |
DeclException1 (ExcPETScError, int,<< "An error with error number "<< arg1<< " occurred while calling a PETSc function") | |
Protected Member Functions | |
virtual void | set_solver_type (KSP &ksp) const |
Protected Attributes | |
const AdditionalData | additional_data |
![]() | |
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 | |
std_cxx1x::shared_ptr< SolverDataMUMPS > | solver_data |
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 1165 of file petsc_solver.h.
PETScWrappers::SparseDirectMUMPS::SparseDirectMUMPS | ( | SolverControl & | cn, |
const MPI_Comm & | mpi_communicator = PETSC_COMM_SELF , |
||
const AdditionalData & | data = AdditionalData() |
||
) |
Constructor
void PETScWrappers::SparseDirectMUMPS::solve | ( | const MatrixBase & | A, |
VectorBase & | x, | ||
const VectorBase & | b | ||
) |
The method to solve the linear system.
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 unstead of more expensive LU. The argument indicates whether the matrix is symmetric or not.
|
protectedvirtual |
Function that takes a Krylov Subspace Solver context object, and sets the type of solver that is requested by the derived class.
Implements PETScWrappers::SolverBase.
|
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.
|
protected |
Store a copy of flags for this particular solver.
Definition at line 1205 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 1254 of file petsc_solver.h.