Reference documentation for deal.II version GIT relicensing-218-g1f873f3929 2024-03-28 15:00:02+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Classes | Public Member Functions | Protected Member Functions | Protected Attributes | Static Private Member Functions | List of all members
PETScWrappers::SparseDirectMUMPS Class Reference

#include <deal.II/lac/petsc_solver.h>

Inheritance diagram for PETScWrappers::SparseDirectMUMPS:
Inheritance graph
[legend]

Classes

struct  AdditionalData
 

Public Member Functions

 SparseDirectMUMPS (SolverControl &cn, const AdditionalData &data=AdditionalData())
 
 SparseDirectMUMPS (SolverControl &cn, const MPI_Comm mpi_communicator, const AdditionalData &data=AdditionalData())
 
void solve (const MatrixBase &A, VectorBase &x, const VectorBase &b)
 
void set_symmetric_mode (const bool matrix_is_symmetric)
 
void solve (const MatrixBase &A, VectorBase &x, const VectorBase &b, const PreconditionBase &preconditioner)
 
virtual void reset ()
 
void set_prefix (const std::string &prefix)
 
SolverControlcontrol () const
 
void initialize (const PreconditionBase &preconditioner)
 
KSP petsc_ksp ()
 
 operator KSP () const
 

Protected Member Functions

virtual void set_solver_type (KSP &ksp) const override
 
void initialize_ksp_with_comm (const MPI_Comm comm)
 
void perhaps_set_convergence_test () const
 

Protected Attributes

const AdditionalData additional_data
 
bool symmetric_mode
 
KSP ksp
 
SmartPointer< SolverControl, SolverBasesolver_control
 
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)
 

Detailed Description

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:

PETScWrappers::SparseDirectMUMPS solver(cn, mpi_communicator);
solver.set_symmetric_mode(true);
solver.solve(system_matrix, solution, system_rhs);
Note
The class internally calls KSPSetFromOptions thus you are able to use all the PETSc parameters for MATSOLVERMUMPS package. See http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MATSOLVERMUMPS.html

Definition at line 903 of file petsc_solver.h.

Constructor & Destructor Documentation

◆ SparseDirectMUMPS() [1/2]

PETScWrappers::SparseDirectMUMPS::SparseDirectMUMPS ( SolverControl cn,
const AdditionalData data = AdditionalData() 
)

Constructor.

Definition at line 641 of file petsc_solver.cc.

◆ SparseDirectMUMPS() [2/2]

PETScWrappers::SparseDirectMUMPS::SparseDirectMUMPS ( SolverControl cn,
const MPI_Comm  mpi_communicator,
const AdditionalData data = AdditionalData() 
)

Constructor. This constructor is deprecated and ignores the MPI communicator argument. Use the other constructor instead.

Deprecated:

Definition at line 650 of file petsc_solver.cc.

Member Function Documentation

◆ solve() [1/2]

void PETScWrappers::SparseDirectMUMPS::solve ( const MatrixBase A,
VectorBase x,
const VectorBase b 
)

The method to solve the linear system.

Definition at line 684 of file petsc_solver.cc.

◆ set_symmetric_mode()

void PETScWrappers::SparseDirectMUMPS::set_symmetric_mode ( const bool  matrix_is_symmetric)

If called with true as argument, tell the direct solver to assume that the system matrix is symmetric. It does so by computing the LDL^T decomposition (in effect, a Cholesky decomposition) instead of more expensive LU decomposition. The argument indicates whether the matrix can be assumed to be symmetric or not.

Note that most finite element matrices are "structurally symmetric", i.e., the sparsity pattern is symmetric, even though the matrix is not. An example of a matrix that is structurally symmetric but not symmetric is the matrix you obtain by discretizing the advection equation \(\nabla \cdot (\vec\beta u) = f\) (see, for example step-12). Because the operator here is not symmetric, the matrix is not symmetric either; however, if matrix entry \(A_{ij}\) is nonzero, then matrix entry \(A_{ji}\) is generally not zero either (and in any case, DoFTools::make_sparsity_pattern() will create a symmetric sparsity pattern). That said, the current function is not meant to indicate whether the sparsity pattern is symmetric, but whether the matrix itself is symmetric, and this typically requires that the differential operator you are considering is symmetric.

Definition at line 814 of file petsc_solver.cc.

◆ set_solver_type()

void PETScWrappers::SparseDirectMUMPS::set_solver_type ( KSP &  ksp) const
overrideprotectedvirtual

Function that takes a Krylov Subspace Solver context object, and sets the type of solver that is requested by the derived class.

Reimplemented from PETScWrappers::SolverBase.

Definition at line 659 of file petsc_solver.cc.

◆ solve() [2/2]

void SolverBase< VectorType >::solve ( const MatrixBase A,
VectorBase x,
const VectorBase b,
const PreconditionBase preconditioner 
)
inherited

Solve the linear system Ax=b. Depending on the information provided by derived classes and the object passed as a preconditioner, one of the linear solvers and preconditioners of PETSc is chosen. Repeated calls to solve() do not reconstruct the preconditioner for performance reasons. See class Documentation.

Definition at line 83 of file petsc_solver.cc.

◆ reset()

void SolverBase< VectorType >::reset ( )
virtualinherited

Resets the contained preconditioner and solver object. See class description for more details.

Definition at line 148 of file petsc_solver.cc.

◆ set_prefix()

void SolverBase< VectorType >::set_prefix ( const std::string &  prefix)
inherited

Sets a prefix name for the solver object. Useful when customizing the PETSc KSP object with command-line options.

Definition at line 141 of file petsc_solver.cc.

◆ control()

SolverControl & SolverBase< VectorType >::control ( ) const
inherited

Access to object that controls convergence.

Definition at line 155 of file petsc_solver.cc.

◆ initialize()

void SolverBase< VectorType >::initialize ( const PreconditionBase preconditioner)
inherited

initialize the solver with the preconditioner. This function is intended for use with SLEPc spectral transformation class.

Definition at line 231 of file petsc_solver.cc.

◆ petsc_ksp()

KSP SolverBase< VectorType >::petsc_ksp ( )
inherited

Return the PETSc KSP object.

Definition at line 68 of file petsc_solver.cc.

◆ operator KSP()

SolverBase< VectorType >::operator KSP ( ) const
inherited

Conversion operator to gain access to the underlying PETSc type. If you do this, you cut this class off some information it may need, so this conversion operator should only be used if you know what you do.

Definition at line 75 of file petsc_solver.cc.

◆ initialize_ksp_with_comm()

void SolverBase< VectorType >::initialize_ksp_with_comm ( const MPI_Comm  comm)
protectedinherited

Utility to create the KSP object and attach convergence test.

Definition at line 206 of file petsc_solver.cc.

◆ perhaps_set_convergence_test()

void SolverBase< VectorType >::perhaps_set_convergence_test ( ) const
protectedinherited

Utility to use deal.II convergence testing.

This call changes the convergence criterion when the instance of the class has a SolverControl object associated.

Definition at line 222 of file petsc_solver.cc.

◆ convergence_test()

PetscErrorCode SolverBase< VectorType >::convergence_test ( KSP  ksp,
const PetscInt  iteration,
const PetscReal  residual_norm,
KSPConvergedReason *  reason,
void *  solver_control 
)
staticprivateinherited

A function that is used in PETSc as a callback to check on 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 166 of file petsc_solver.cc.

Member Data Documentation

◆ additional_data

const AdditionalData PETScWrappers::SparseDirectMUMPS::additional_data
protected

Store a copy of flags for this particular solver.

Definition at line 965 of file petsc_solver.h.

◆ symmetric_mode

bool PETScWrappers::SparseDirectMUMPS::symmetric_mode
protected

Flag specifies whether matrix being factorized is symmetric or not. It influences the type of the used preconditioner (PCLU or PCCHOLESKY)

Definition at line 974 of file petsc_solver.h.

◆ ksp

KSP PETScWrappers::SolverBase::ksp
protectedinherited

The PETSc KSP object.

Definition at line 167 of file petsc_solver.h.

◆ solver_control

SmartPointer<SolverControl, SolverBase> PETScWrappers::SolverBase::solver_control
protectedinherited

Reference to the object that controls convergence of the iterative solver. In fact, for these PETSc wrappers, PETSc does so itself, but we copy the data from this object before starting the solution process, and copy the data back into it afterwards.

Definition at line 175 of file petsc_solver.h.

◆ prefix_name

std::string PETScWrappers::SolverBase::prefix_name
protectedinherited

Solver prefix name to qualify options specific to the PETSc KSP object in the current context. Note: A hyphen (-) must NOT be given at the beginning of the prefix name. The first character of all runtime options is AUTOMATICALLY the hyphen.

Definition at line 205 of file petsc_solver.h.


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