Reference documentation for deal.II version GIT relicensing-136-gb80d0be4af 2024-03-18 08:20: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

Classes

class  SolverArnoldi
 
class  SolverBase
 
class  SolverGeneralizedDavidson
 
class  SolverJacobiDavidson
 
class  SolverKrylovSchur
 
class  SolverLanczos
 
class  SolverLAPACK
 
class  SolverPower
 
class  TransformationBase
 
class  TransformationCayley
 
class  TransformationShift
 
class  TransformationShiftInvert
 
class  TransformationSpectrumFolding
 

Detailed Description

Base namespace for solver classes using the SLEPc solvers which are selected based on flags passed to the eigenvalue problem solver context. Derived classes set the right flags to set the right solver.

The SLEPc solvers are intended to be used for solving the generalized eigenspectrum problem \((A-\lambda B)x=0\), for \(x\neq0\); where \(A\) is a system matrix, \(B\) is a mass matrix, and \(\lambda, x\) are a set of eigenvalues and eigenvectors respectively. The emphasis is on methods and techniques appropriate for problems in which the associated matrices are sparse. Most of the methods offered by the SLEPc library are projection methods or other methods with similar properties; and wrappers are provided to interface to SLEPc solvers that handle both of these problem sets.

SLEPcWrappers can be implemented in application codes in the following way:

SolverControl solver_control (1000, 1e-9);
SolverArnoldi system (solver_control, mpi_communicator);
system.solve (A, B, lambda, x, size_of_spectrum);

for the generalized eigenvalue problem \(Ax=B\lambda x\), where the variable const unsigned int size_of_spectrum tells SLEPc the number of eigenvector/eigenvalue pairs to solve for. Additional options and solver parameters can be passed to the SLEPc solvers before calling solve(). For example, if the matrices of the general eigenspectrum problem are not hermitian and the lower eigenvalues are wanted only, the following code can be implemented before calling solve():

system.set_problem_type (EPS_NHEP);
system.set_which_eigenpairs (EPS_SMALLEST_REAL);

These options can also be set at the command line.

See also step-36 for a hands-on example.

For cases when spectral transformations are used in conjunction with Krylov-type solvers or Davidson-type eigensolvers are employed one can additionally specify which linear solver and preconditioner to use. This can be achieved as follows

data.symmetric_operator = true;
PETScWrappers::PreconditionBoomerAMG preconditioner(mpi_communicator, data);
SolverControl linear_solver_control (dof_handler.n_dofs(),
1e-12, false, false);
PETScWrappers::SolverCG linear_solver(linear_solver_control,
mpi_communicator);
linear_solver.initialize(preconditioner);
SolverControl solver_control (100, 1e-9,false,false);
SLEPcWrappers::SolverKrylovSchur eigensolver(solver_control,
mpi_communicator);
SLEPcWrappers::TransformationShift spectral_transformation(mpi_communicator);
spectral_transformation.set_solver(linear_solver);
eigensolver.set_transformation(spectral_transformation);
eigensolver.solve(stiffness_matrix, mass_matrix,
eigenvalues, eigenfunctions, eigenfunctions.size());
std::array< Number, 1 > eigenvalues(const SymmetricTensor< 2, 1, Number > &T)

In order to support this usage case, different from PETSc wrappers, the classes in this namespace are written in such a way that the underlying SLEPc objects are initialized in constructors. By doing so one also avoid caching of different settings (such as target eigenvalue or type of the problem); instead those are applied straight away when the corresponding functions of the wrapper classes are called.

An alternative implementation to the one above is to use the API internals directly within the application code. In this way the calling sequence requires calling several of SolverBase functions rather than just one. This freedom is intended for use of the SLEPcWrappers that require a greater handle on the eigenvalue problem solver context. See also the API of, for example:

template <typename OutputVector>
void
std::vector<PetscScalar> &eigenvalues,
std::vector<OutputVector> &eigenvectors,
const unsigned int n_eigenpairs)
{
...
}
void solve(const PETScWrappers::MatrixBase &A, std::vector< PetscScalar > &eigenvalues, std::vector< OutputVector > &eigenvectors, const unsigned int n_eigenpairs=1)
std::array< std::pair< Number, Tensor< 1, dim, Number > >, std::integral_constant< int, dim >::value > eigenvectors(const SymmetricTensor< 2, dim, Number > &T, const SymmetricTensorEigenvectorMethod method=SymmetricTensorEigenvectorMethod::ql_implicit_shifts)

as an example on how to do this.

For further information and explanations on handling the SLEPcWrappers, see also the PETScWrappers, on which they depend.