MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
slepc.hpp
Go to the documentation of this file.
1// Copyright (c) 2010-2024, Lawrence Livermore National Security, LLC. Produced
2// at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3// LICENSE and NOTICE for details. LLNL-CODE-806117.
4//
5// This file is part of the MFEM library. For more information and source code
6// availability visit https://mfem.org.
7//
8// MFEM is free software; you can redistribute it and/or modify it under the
9// terms of the BSD-3 license. We welcome feedback and contributions, see file
10// CONTRIBUTING.md for details.
11
12#ifndef MFEM_SLEPC
13#define MFEM_SLEPC
14
15#include "../config/config.hpp"
16
17#ifdef MFEM_USE_SLEPC
18#ifdef MFEM_USE_MPI
19
20#include "petsc.hpp"
21
22// Forward declaration of SLEPc's internal struct _p_EPS:
23struct _p_EPS;
24
25namespace mfem
26{
27
28// Declare an alias of SLEPc's EPS type, mfem::slepc::EPS:
29namespace slepc { typedef struct ::_p_EPS *EPS; }
30
32void MFEMInitializeSlepc(int*,char***);
33void MFEMInitializeSlepc(int*,char***,const char[],const char[]);
35
37{
38private:
39 /// Boolean to handle SetFromOptions calls
40 mutable bool clcustom;
41
42 /// SLEPc linear eigensolver object
43 slepc::EPS eps;
44
45 /// Real and imaginary part of eigenvector
46 mutable PetscParVector *VR, *VC;
47
48public:
49 /// Constructors
50 SlepcEigenSolver(MPI_Comm comm, const std::string &prefix = std::string());
51
52 virtual ~SlepcEigenSolver();
53
54 /** @brief Set solver convergence tolerance relative to the magnitude of the
55 eigenvalue.
56
57 @note Default value is 1e-8
58 */
59 void SetTol(real_t tol);
60
61 /** @brief Set maximum number of iterations allowed in the call to
62 SlepcEigenSolver::Solve */
63 void SetMaxIter(int max_iter);
64
65 /// Set the number of eigenmodes to compute
66 void SetNumModes(int num_eigs);
67
68 /// Set operator for standard eigenvalue problem
69 void SetOperator(const PetscParMatrix &op);
70
71 /// Set operators for generalized eigenvalue problem
72 void SetOperators(const PetscParMatrix &op, const PetscParMatrix &opB);
73
74 /// Customize object with options set
75 void Customize(bool customize = true) const;
76
77 /// Solve the eigenvalue problem for the specified number of eigenvalues
78 void Solve();
79
80 /** @brief Get the number of converged eigenvalues after the call to
81 SlepcEigenSolver::Solve */
82 int GetNumConverged();
83
84 /** @brief Get the ith eigenvalue after the system has been solved
85 @param[in] i The index for the eigenvalue you want ordered by
86 SlepcEigenSolver::SetWhichEigenpairs
87 @param[out] lr The real component of the eigenvalue
88 @note the index @a i must be between 0 and
89 SlepcEigenSolver::GetNumConverged - 1
90 */
91 void GetEigenvalue(unsigned int i, real_t & lr) const;
92
93 /** @brief Get the ith eigenvalue after the system has been solved
94 @param[in] i The index for the eigenvalue you want ordered by
95 SlepcEigenSolver::SetWhichEigenpairs
96 @param[out] lr The real component of the eigenvalue
97 @param[out] lc The imaginary component of the eigenvalue
98 @note the index @a i must be between 0 and
99 SlepcEigenSolver::GetNumConverged - 1
100 */
101 void GetEigenvalue(unsigned int i, real_t & lr, real_t & lc) const;
102
103 /** @brief Get the ith eigenvector after the system has been solved
104 @param[in] i The index for the eigenvector you want ordered by
105 SlepcEigenSolver::SetWhichEigenpairs
106 @param[out] vr The real components of the eigenvector
107 @note the index @a i must be between 0 and
108 SlepcEigenSolver::GetNumConverged - 1
109 */
110 void GetEigenvector(unsigned int i, Vector & vr) const;
111
112 /** @brief Get the ith eigenvector after the system has been solved
113 @param[in] i The index for the eigenvector you want ordered by
114 SlepcEigenSolver::SetWhichEigenpairs
115 @param[out] vr The real components of the eigenvector
116 @param[out] vc The imaginary components of the eigenvector
117 @note the index @a i must be between 0 and
118 SlepcEigenSolver::GetNumConverged - 1
119 */
120 void GetEigenvector(unsigned int i, Vector & vr, Vector & vc) const;
121
122 /** @brief Target spectrum for the eigensolver.
123
124 This will define the order in which the eigenvalues/eigenvectors are
125 indexed after the call to SlepcEigenSolver::Solve.
126 @note Target imaginary is not supported without complex support in SLEPc,
127 and intervals are not implemented.
128 */
129 enum Which
130 {
131 /// The eigenvalues with the largest complex magnitude (default)
133 /// The eigenvalues with the smallest complex magnitude
135 /// The eigenvalues with the largest real component
137 /// The eigenvalues with the smallest real component
139 /// The eigenvalues with the largest imaginary component
141 /// The eigenvalues with the smallest imaginary component
143 /// The eigenvalues with complex magnitude closest to the target value
145 /// The eigenvalues with the real component closest to the target value
147 };
148
149 /** @brief Spectral transformations that can be used by the solver in order
150 to accelerate the convergence to the target eignevalues
151 */
153 {
154 /// Utilize the shift of origin strategy
156 /// Utilize the shift and invert strategy
158 };
159
160 /** @brief Set the which eigenvalues the solver will target and the order
161 they will be indexed in.
162
163 For SlepcEigenSolver::TARGET_MAGNITUDE or SlepcEigenSolver::TARGET_REAL
164 you will also need to set the target value with
165 SlepcEigenSolver::SetTarget.
166 */
167 void SetWhichEigenpairs(Which which);
168
169 /** @brief Set the target value for the eigenpairs you want when using
170 SlepcEigenSolver::TARGET_MAGNITUDE or SlepcEigenSolver::TARGET_REAL in
171 the SlepcEigenSolver::SetWhichEigenpairs method.
172 */
173 void SetTarget(real_t target);
174
175 /** @brief Set the spectral transformation strategy for acceletating
176 convergenvce. Both SlepcEigenSolver::SHIFT and
177 SlepcEigenSolver::SHIFT_INVERT are available.
178 */
180
181 /// Conversion function to SLEPc's EPS type.
182 operator slepc::EPS() const { return eps; }
183
184 /// Conversion function to PetscObject
185 operator PetscObject() const {return (PetscObject)eps; }
186};
187
188}
189
190#endif // MFEM_USE_MPI
191#endif // MFEM_USE_SLEPC
192
193#endif // MFEM_SLEPC
Wrapper for PETSc's matrix class.
Definition petsc.hpp:319
Which
Target spectrum for the eigensolver.
Definition slepc.hpp:130
@ LARGEST_REAL
The eigenvalues with the largest real component.
Definition slepc.hpp:136
@ SMALLEST_MAGNITUDE
The eigenvalues with the smallest complex magnitude.
Definition slepc.hpp:134
@ SMALLEST_REAL
The eigenvalues with the smallest real component.
Definition slepc.hpp:138
@ TARGET_MAGNITUDE
The eigenvalues with complex magnitude closest to the target value.
Definition slepc.hpp:144
@ LARGEST_MAGNITUDE
The eigenvalues with the largest complex magnitude (default)
Definition slepc.hpp:132
@ TARGET_REAL
The eigenvalues with the real component closest to the target value.
Definition slepc.hpp:146
@ SMALLEST_IMAGINARY
The eigenvalues with the smallest imaginary component.
Definition slepc.hpp:142
@ LARGEST_IMAGINARY
The eigenvalues with the largest imaginary component.
Definition slepc.hpp:140
void Customize(bool customize=true) const
Customize object with options set.
Definition slepc.cpp:137
void SetOperators(const PetscParMatrix &op, const PetscParMatrix &opB)
Set operators for generalized eigenvalue problem.
Definition slepc.cpp:93
void SetWhichEigenpairs(Which which)
Set the which eigenvalues the solver will target and the order they will be indexed in.
Definition slepc.cpp:195
SlepcEigenSolver(MPI_Comm comm, const std::string &prefix=std::string())
Constructors.
Definition slepc.cpp:60
int GetNumConverged()
Get the number of converged eigenvalues after the call to SlepcEigenSolver::Solve.
Definition slepc.cpp:188
void SetTarget(real_t target)
Set the target value for the eigenpairs you want when using SlepcEigenSolver::TARGET_MAGNITUDE or Sle...
Definition slepc.cpp:229
void GetEigenvalue(unsigned int i, real_t &lr) const
Get the ith eigenvalue after the system has been solved.
Definition slepc.cpp:147
void SetNumModes(int num_eigs)
Set the number of eigenmodes to compute.
Definition slepc.cpp:124
virtual ~SlepcEigenSolver()
Definition slepc.cpp:70
void SetSpectralTransformation(SpectralTransformation transformation)
Set the spectral transformation strategy for acceletating convergenvce. Both SlepcEigenSolver::SHIFT ...
Definition slepc.cpp:234
void SetMaxIter(int max_iter)
Set maximum number of iterations allowed in the call to SlepcEigenSolver::Solve.
Definition slepc.cpp:116
SpectralTransformation
Spectral transformations that can be used by the solver in order to accelerate the convergence to the...
Definition slepc.hpp:153
@ SHIFT
Utilize the shift of origin strategy.
Definition slepc.hpp:155
@ SHIFT_INVERT
Utilize the shift and invert strategy.
Definition slepc.hpp:157
void GetEigenvector(unsigned int i, Vector &vr) const
Get the ith eigenvector after the system has been solved.
Definition slepc.cpp:158
void SetTol(real_t tol)
Set solver convergence tolerance relative to the magnitude of the eigenvalue.
Definition slepc.cpp:106
void SetOperator(const PetscParMatrix &op)
Set operator for standard eigenvalue problem.
Definition slepc.cpp:81
void Solve()
Solve the eigenvalue problem for the specified number of eigenvalues.
Definition slepc.cpp:130
Vector data type.
Definition vector.hpp:80
void transformation(const Vector &p, Vector &v)
struct ::_p_EPS * EPS
Definition slepc.hpp:29
void MFEMFinalizeSlepc()
Definition slepc.cpp:53
void MFEMInitializeSlepc()
Definition slepc.cpp:29
float real_t
Definition config.hpp:43
struct _p_PetscObject * PetscObject
Definition petsc.hpp:54