MFEM  v4.2.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
slepc.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2020, 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 
20 #include "slepc.h"
21 
22 #include "petscinternals.hpp"
23 
24 static PetscErrorCode ierr;
25 
26 namespace mfem
27 {
28 
30 {
31  MFEMInitializeSlepc(NULL,NULL,NULL,NULL);
32 }
33 
34 void MFEMInitializeSlepc(int *argc,char*** argv)
35 {
36  MFEMInitializeSlepc(argc,argv,NULL,NULL);
37 }
38 
39 void MFEMInitializeSlepc(int *argc,char ***argv,const char rc_file[],
40  const char help[])
41 {
42  ierr = SlepcInitialize(argc,argv,rc_file,help);
43  MFEM_VERIFY(!ierr,"Unable to initialize SLEPc");
44 }
45 
47 {
48  ierr = SlepcFinalize();
49  MFEM_VERIFY(!ierr,"Unable to finalize SLEPc");
50 }
51 
52 
53 SlepcEigenSolver::SlepcEigenSolver(MPI_Comm comm, const std::string &prefix)
54 {
55  clcustom = false;
56  VR = NULL;
57  VC = NULL;
58 
59  ierr = EPSCreate(comm,&eps); CCHKERRQ(comm,ierr);
60  ierr = EPSSetOptionsPrefix(eps, prefix.c_str()); PCHKERRQ(eps, ierr);
61 }
62 
64 {
65  delete VR;
66  delete VC;
67 
68  MPI_Comm comm;
69  ierr = PetscObjectGetComm((PetscObject)eps,&comm); PCHKERRQ(eps,ierr);
70  ierr = EPSDestroy(&eps); CCHKERRQ(comm,ierr);
71 }
72 
73 
75 {
76  delete VR;
77  delete VC;
78  VR = VC = NULL;
79 
80  ierr = EPSSetOperators(eps,op,NULL); PCHKERRQ(eps, ierr);
81 
82  VR = new PetscParVector(op, true, false);
83  VC = new PetscParVector(op, true, false);
84 
85 }
86 
88  const PetscParMatrix&opB)
89 {
90  delete VR;
91  delete VC;
92  VR = VC = NULL;
93 
94  ierr = EPSSetOperators(eps,op,opB); PCHKERRQ(eps,ierr);
95 
96  VR = new PetscParVector(op, true, false);
97  VC = new PetscParVector(op, true, false);
98 }
99 
100 void SlepcEigenSolver::SetTol(double tol)
101 {
102  int max_its;
103 
104  ierr = EPSGetTolerances(eps,NULL,&max_its); PCHKERRQ(eps,ierr);
105  // Work around uninitialized maximum iterations
106  if (max_its==0) { max_its = PETSC_DECIDE; }
107  ierr = EPSSetTolerances(eps,tol,max_its); PCHKERRQ(eps,ierr);
108 }
109 
111 {
112  double tol;
113 
114  ierr = EPSGetTolerances(eps,&tol,NULL); PCHKERRQ(eps,ierr);
115  ierr = EPSSetTolerances(eps,tol,max_its); PCHKERRQ(eps,ierr);
116 }
117 
119 {
120  ierr = EPSSetDimensions(eps,num_eigs,PETSC_DECIDE,PETSC_DECIDE);
121  PCHKERRQ(eps,ierr);
122 }
123 
125 {
126  Customize();
127 
128  ierr = EPSSolve(eps); PCHKERRQ(eps,ierr);
129 }
130 
131 void SlepcEigenSolver::Customize(bool customize) const
132 {
133  if (!customize) {clcustom = true; }
134  if (!clcustom)
135  {
136  ierr = EPSSetFromOptions(eps); PCHKERRQ(eps,ierr);
137  }
138  clcustom = true;
139 }
140 
141 void SlepcEigenSolver::GetEigenvalue(unsigned int i, double & lr) const
142 {
143  ierr = EPSGetEigenvalue(eps,i,&lr,NULL); PCHKERRQ(eps,ierr);
144 }
145 
146 void SlepcEigenSolver::GetEigenvalue(unsigned int i, double & lr,
147  double & lc) const
148 {
149  ierr = EPSGetEigenvalue(eps,i,&lr,&lc); PCHKERRQ(eps,ierr);
150 }
151 
152 void SlepcEigenSolver::GetEigenvector(unsigned int i, Vector & vr) const
153 {
154  MFEM_VERIFY(VR,"Missing real vector");
155 
156  MFEM_ASSERT(vr.Size() == VR->Size(), "invalid vr.Size() = " << vr.Size()
157  << ", expected size = " << VR->Size());
158 
159  VR->PlaceArray(vr.GetData());
160  ierr = EPSGetEigenvector(eps,i,*VR,NULL); PCHKERRQ(eps,ierr);
161  VR->ResetArray();
162 
163 }
164 
165 void SlepcEigenSolver::GetEigenvector(unsigned int i, Vector & vr,
166  Vector & vc) const
167 {
168  MFEM_VERIFY(VR,"Missing real vector");
169  MFEM_VERIFY(VC,"Missing imaginary vector");
170  MFEM_ASSERT(vr.Size() == VR->Size(), "invalid vr.Size() = " << vr.Size()
171  << ", expected size = " << VR->Size());
172  MFEM_ASSERT(vc.Size() == VC->Size(), "invalid vc.Size() = " << vc.Size()
173  << ", expected size = " << VC->Size());
174 
175  VR->PlaceArray(vr.GetData());
176  VC->PlaceArray(vc.GetData());
177  ierr = EPSGetEigenvector(eps,i,*VR,*VC); PCHKERRQ(eps,ierr);
178  VR->ResetArray();
179  VC->ResetArray();
180 }
181 
183 {
184  int num_conv;
185  ierr = EPSGetConverged(eps,&num_conv); PCHKERRQ(eps,ierr);
186  return num_conv;
187 }
188 
190 {
191  switch (which)
192  {
194  ierr = EPSSetWhichEigenpairs(eps,EPS_LARGEST_MAGNITUDE); PCHKERRQ(eps,ierr);
195  break;
197  ierr = EPSSetWhichEigenpairs(eps,EPS_SMALLEST_MAGNITUDE); PCHKERRQ(eps,ierr);
198  break;
200  ierr = EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL); PCHKERRQ(eps,ierr);
201  break;
203  ierr = EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL); PCHKERRQ(eps,ierr);
204  break;
206  ierr = EPSSetWhichEigenpairs(eps,EPS_LARGEST_IMAGINARY); PCHKERRQ(eps,ierr);
207  break;
209  ierr = EPSSetWhichEigenpairs(eps,EPS_SMALLEST_IMAGINARY); PCHKERRQ(eps,ierr);
210  break;
212  ierr = EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE); PCHKERRQ(eps,ierr);
213  break;
215  ierr = EPSSetWhichEigenpairs(eps,EPS_TARGET_REAL); PCHKERRQ(eps,ierr);
216  break;
217  default:
218  MFEM_ABORT("Which eigenpair not implemented!");
219  break;
220  }
221 }
222 
223 void SlepcEigenSolver::SetTarget(double target)
224 {
225  ierr = EPSSetTarget(eps,target); PCHKERRQ(eps,ierr);
226 }
227 
230 {
231  ST st;
232  ierr = EPSGetST(eps,&st); PCHKERRQ(eps,ierr);
233  switch (transformation)
234  {
236  ierr = STSetType(st,STSHIFT); PCHKERRQ(eps,ierr);
237  break;
239  ierr = STSetType(st,STSINVERT); PCHKERRQ(eps,ierr);
240  break;
241  default:
242  MFEM_ABORT("Spectral transformation not implemented!");
243  break;
244  }
245 }
246 
247 }
248 
249 #endif // MFEM_USE_SLEPC
250 #endif // MFEM_USE_PETSC
251 #endif // MFEM_USE_MPI
void SetOperator(const PetscParMatrix &op)
Set operator for standard eigenvalue problem.
Definition: slepc.cpp:74
void Customize(bool customize=true) const
Customize object with options set.
Definition: slepc.cpp:131
void MFEMInitializeSlepc()
Definition: slepc.cpp:29
Wrapper for PETSc&#39;s matrix class.
Definition: petsc.hpp:197
void GetEigenvalue(unsigned int i, double &lr) const
Get the corresponding eigenvalue.
Definition: slepc.cpp:141
int Size() const
Returns the size of the vector.
Definition: vector.hpp:160
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:169
void SetTarget(double target)
Definition: slepc.cpp:223
Wrapper for PETSc&#39;s vector class.
Definition: petsc.hpp:75
void ResetArray()
Reset the PETSc Vec object to use its default data. Call this method after the use of PlaceArray()...
Definition: petsc.cpp:415
void SetSpectralTransformation(SpectralTransformation transformation)
Definition: slepc.cpp:228
virtual ~SlepcEigenSolver()
Definition: slepc.cpp:63
void transformation(const Vector &p, Vector &v)
void SetMaxIter(int max_iter)
Set maximum number of iterations.
Definition: slepc.cpp:110
void SetNumModes(int num_eigs)
Set the number of required eigenmodes.
Definition: slepc.cpp:118
void Solve()
Solve the eigenvalue problem for the specified number of eigenvalues.
Definition: slepc.cpp:124
void MFEMFinalizeSlepc()
Definition: slepc.cpp:46
void SetTol(double tol)
Set solver tolerance.
Definition: slepc.cpp:100
SlepcEigenSolver(MPI_Comm comm, const std::string &prefix=std::string())
Constructors.
Definition: slepc.cpp:53
int GetNumConverged()
Get the number of converged eigenvalues.
Definition: slepc.cpp:182
void SetWhichEigenpairs(Which which)
Definition: slepc.cpp:189
struct _p_PetscObject * PetscObject
Definition: petsc.hpp:50
Vector data type.
Definition: vector.hpp:51
void PlaceArray(PetscScalar *temp_data)
Temporarily replace the data of the PETSc Vec object. To return to the original data array...
Definition: petsc.cpp:410
void GetEigenvector(unsigned int i, Vector &vr) const
Get the corresponding eigenvector.
Definition: slepc.cpp:152
void SetOperators(const PetscParMatrix &op, const PetscParMatrix &opB)
Set operator for generalized eigenvalue problem.
Definition: slepc.cpp:87