MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
slepc.cpp
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#include "../config/config.hpp"
13
14#ifdef MFEM_USE_MPI
15#ifdef MFEM_USE_PETSC
16#ifdef MFEM_USE_SLEPC
17
18#include "linalg.hpp"
19#include "petscinternals.hpp"
20#include "slepc.h"
21
22static PetscErrorCode ierr;
23
24using namespace std;
25
26namespace mfem
27{
28
30{
31 MFEMInitializeSlepc(NULL,NULL,NULL,NULL);
32}
33
34void MFEMInitializeSlepc(int *argc,char*** argv)
35{
36 MFEMInitializeSlepc(argc,argv,NULL,NULL);
37}
38
39void MFEMInitializeSlepc(int *argc,char ***argv,const char rc_file[],
40 const char help[])
41{
43 {
44 // Tell PETSc to use the same CUDA device as MFEM:
45 ierr = PetscOptionsSetValue(NULL,"-cuda_device",
46 to_string(mfem::Device::GetId()).c_str());
47 MFEM_VERIFY(!ierr,"Unable to set initial option value to PETSc");
48 }
49 ierr = SlepcInitialize(argc,argv,rc_file,help);
50 MFEM_VERIFY(!ierr,"Unable to initialize SLEPc");
51}
52
54{
55 ierr = SlepcFinalize();
56 MFEM_VERIFY(!ierr,"Unable to finalize SLEPc");
57}
58
59
60SlepcEigenSolver::SlepcEigenSolver(MPI_Comm comm, const std::string &prefix)
61{
62 clcustom = false;
63 VR = NULL;
64 VC = NULL;
65
66 ierr = EPSCreate(comm,&eps); CCHKERRQ(comm,ierr);
67 ierr = EPSSetOptionsPrefix(eps, prefix.c_str()); PCHKERRQ(eps, ierr);
68}
69
71{
72 delete VR;
73 delete VC;
74
75 MPI_Comm comm;
76 ierr = PetscObjectGetComm((PetscObject)eps,&comm); PCHKERRQ(eps,ierr);
77 ierr = EPSDestroy(&eps); CCHKERRQ(comm,ierr);
78}
79
80
82{
83 delete VR;
84 delete VC;
85 VR = VC = NULL;
86
87 ierr = EPSSetOperators(eps,op,NULL); PCHKERRQ(eps, ierr);
88
89 VR = new PetscParVector(op, true, false);
90 VC = new PetscParVector(op, true, false);
91}
92
94 const PetscParMatrix&opB)
95{
96 delete VR;
97 delete VC;
98 VR = VC = NULL;
99
100 ierr = EPSSetOperators(eps,op,opB); PCHKERRQ(eps,ierr);
101
102 VR = new PetscParVector(op, true, false);
103 VC = new PetscParVector(op, true, false);
104}
105
107{
108 PetscInt max_its;
109
110 ierr = EPSGetTolerances(eps,NULL,&max_its); PCHKERRQ(eps,ierr);
111 // Work around uninitialized maximum iterations
112 if (max_its==0) { max_its = PETSC_DECIDE; }
113 ierr = EPSSetTolerances(eps,tol,max_its); PCHKERRQ(eps,ierr);
114}
115
117{
118 real_t tol;
119
120 ierr = EPSGetTolerances(eps,&tol,NULL); PCHKERRQ(eps,ierr);
121 ierr = EPSSetTolerances(eps,tol,max_its); PCHKERRQ(eps,ierr);
122}
123
125{
126 ierr = EPSSetDimensions(eps,num_eigs,PETSC_DECIDE,PETSC_DECIDE);
127 PCHKERRQ(eps,ierr);
128}
129
131{
132 Customize();
133
134 ierr = EPSSolve(eps); PCHKERRQ(eps,ierr);
135}
136
137void SlepcEigenSolver::Customize(bool customize) const
138{
139 if (!customize) {clcustom = true; }
140 if (!clcustom)
141 {
142 ierr = EPSSetFromOptions(eps); PCHKERRQ(eps,ierr);
143 }
144 clcustom = true;
145}
146
147void SlepcEigenSolver::GetEigenvalue(unsigned int i, real_t & lr) const
148{
149 ierr = EPSGetEigenvalue(eps,i,&lr,NULL); PCHKERRQ(eps,ierr);
150}
151
152void SlepcEigenSolver::GetEigenvalue(unsigned int i, real_t & lr,
153 real_t & lc) const
154{
155 ierr = EPSGetEigenvalue(eps,i,&lr,&lc); PCHKERRQ(eps,ierr);
156}
157
158void SlepcEigenSolver::GetEigenvector(unsigned int i, Vector & vr) const
159{
160 MFEM_VERIFY(VR,"Missing real vector");
161
162 MFEM_ASSERT(vr.Size() == VR->Size(), "invalid vr.Size() = " << vr.Size()
163 << ", expected size = " << VR->Size());
164
165 VR->PlaceMemory(vr.GetMemory());
166 ierr = EPSGetEigenvector(eps,i,*VR,NULL); PCHKERRQ(eps,ierr);
167 VR->ResetMemory();
168
169}
170
171void SlepcEigenSolver::GetEigenvector(unsigned int i, Vector & vr,
172 Vector & vc) const
173{
174 MFEM_VERIFY(VR,"Missing real vector");
175 MFEM_VERIFY(VC,"Missing imaginary vector");
176 MFEM_ASSERT(vr.Size() == VR->Size(), "invalid vr.Size() = " << vr.Size()
177 << ", expected size = " << VR->Size());
178 MFEM_ASSERT(vc.Size() == VC->Size(), "invalid vc.Size() = " << vc.Size()
179 << ", expected size = " << VC->Size());
180
181 VR->PlaceArray(vr.GetMemory());
182 VC->PlaceArray(vc.GetMemory());
183 ierr = EPSGetEigenvector(eps,i,*VR,*VC); PCHKERRQ(eps,ierr);
184 VR->ResetMemory();
185 VC->ResetMemory();
186}
187
189{
190 PetscInt num_conv;
191 ierr = EPSGetConverged(eps,&num_conv); PCHKERRQ(eps,ierr);
192 return static_cast<int>(num_conv);
193}
194
196{
197 switch (which)
198 {
200 ierr = EPSSetWhichEigenpairs(eps,EPS_LARGEST_MAGNITUDE); PCHKERRQ(eps,ierr);
201 break;
203 ierr = EPSSetWhichEigenpairs(eps,EPS_SMALLEST_MAGNITUDE); PCHKERRQ(eps,ierr);
204 break;
206 ierr = EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL); PCHKERRQ(eps,ierr);
207 break;
209 ierr = EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL); PCHKERRQ(eps,ierr);
210 break;
212 ierr = EPSSetWhichEigenpairs(eps,EPS_LARGEST_IMAGINARY); PCHKERRQ(eps,ierr);
213 break;
215 ierr = EPSSetWhichEigenpairs(eps,EPS_SMALLEST_IMAGINARY); PCHKERRQ(eps,ierr);
216 break;
218 ierr = EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE); PCHKERRQ(eps,ierr);
219 break;
221 ierr = EPSSetWhichEigenpairs(eps,EPS_TARGET_REAL); PCHKERRQ(eps,ierr);
222 break;
223 default:
224 MFEM_ABORT("Which eigenpair not implemented!");
225 break;
226 }
227}
228
230{
231 ierr = EPSSetTarget(eps,target); PCHKERRQ(eps,ierr);
232}
233
236{
237 ST st;
238 ierr = EPSGetST(eps,&st); PCHKERRQ(eps,ierr);
239 switch (transformation)
240 {
242 ierr = STSetType(st,STSHIFT); PCHKERRQ(eps,ierr);
243 break;
245 ierr = STSetType(st,STSINVERT); PCHKERRQ(eps,ierr);
246 break;
247 default:
248 MFEM_ABORT("Spectral transformation not implemented!");
249 break;
250 }
251}
252
253}
254
255#endif // MFEM_USE_SLEPC
256#endif // MFEM_USE_PETSC
257#endif // MFEM_USE_MPI
static bool Allows(unsigned long b_mask)
Return true if any of the backends in the backend mask, b_mask, are allowed.
Definition device.hpp:259
static int GetId()
Get the device id of the configured device.
Definition device.hpp:253
Wrapper for PETSc's matrix class.
Definition petsc.hpp:319
void ResetMemory()
Completes the operation started with PlaceMemory.
Definition petsc.cpp:891
void PlaceMemory(Memory< real_t > &, bool=false)
This requests write access from where the memory is valid and temporarily replaces the corresponding ...
Definition petsc.cpp:804
void PlaceArray(PetscScalar *temp_data)
Temporarily replace the data of the PETSc Vec object. To return to the original data array,...
Definition petsc.cpp:794
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
Memory< real_t > & GetMemory()
Return a reference to the Memory object used by the Vector.
Definition vector.hpp:249
int Size() const
Returns the size of the vector.
Definition vector.hpp:218
void transformation(const Vector &p, Vector &v)
void MFEMFinalizeSlepc()
Definition slepc.cpp:53
void MFEMInitializeSlepc()
Definition slepc.cpp:29
float real_t
Definition config.hpp:43
HYPRE_Int PetscInt
Definition petsc.hpp:50
struct _p_PetscObject * PetscObject
Definition petsc.hpp:54
@ CUDA_MASK
Biwise-OR of all CUDA backends.
Definition device.hpp:89