12 #include "../config/config.hpp"
24 static PetscErrorCode ierr;
42 ierr = SlepcInitialize(argc,argv,rc_file,help);
43 MFEM_VERIFY(!ierr,
"Unable to initialize SLEPc");
48 ierr = SlepcFinalize();
49 MFEM_VERIFY(!ierr,
"Unable to finalize SLEPc");
59 ierr = EPSCreate(comm,&eps); CCHKERRQ(comm,ierr);
60 ierr = EPSSetOptionsPrefix(eps, prefix.c_str()); PCHKERRQ(eps, ierr);
69 ierr = PetscObjectGetComm((
PetscObject)eps,&comm); PCHKERRQ(eps,ierr);
70 ierr = EPSDestroy(&eps); CCHKERRQ(comm,ierr);
80 ierr = EPSSetOperators(eps,op,NULL); PCHKERRQ(eps, ierr);
94 ierr = EPSSetOperators(eps,op,opB); PCHKERRQ(eps,ierr);
104 ierr = EPSGetTolerances(eps,NULL,&max_its); PCHKERRQ(eps,ierr);
106 if (max_its==0) { max_its = PETSC_DECIDE; }
107 ierr = EPSSetTolerances(eps,tol,max_its); PCHKERRQ(eps,ierr);
114 ierr = EPSGetTolerances(eps,&tol,NULL); PCHKERRQ(eps,ierr);
115 ierr = EPSSetTolerances(eps,tol,max_its); PCHKERRQ(eps,ierr);
120 ierr = EPSSetDimensions(eps,num_eigs,PETSC_DECIDE,PETSC_DECIDE);
128 ierr = EPSSolve(eps); PCHKERRQ(eps,ierr);
133 if (!customize) {clcustom =
true; }
136 ierr = EPSSetFromOptions(eps); PCHKERRQ(eps,ierr);
143 ierr = EPSGetEigenvalue(eps,i,&lr,NULL); PCHKERRQ(eps,ierr);
149 ierr = EPSGetEigenvalue(eps,i,&lr,&lc); PCHKERRQ(eps,ierr);
154 MFEM_VERIFY(VR,
"Missing real vector");
156 MFEM_ASSERT(vr.
Size() == VR->
Size(),
"invalid vr.Size() = " << vr.
Size()
157 <<
", expected size = " << VR->
Size());
160 ierr = EPSGetEigenvector(eps,i,*VR,NULL); PCHKERRQ(eps,ierr);
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());
177 ierr = EPSGetEigenvector(eps,i,*VR,*VC); PCHKERRQ(eps,ierr);
185 ierr = EPSGetConverged(eps,&num_conv); PCHKERRQ(eps,ierr);
194 ierr = EPSSetWhichEigenpairs(eps,EPS_LARGEST_MAGNITUDE); PCHKERRQ(eps,ierr);
197 ierr = EPSSetWhichEigenpairs(eps,EPS_SMALLEST_MAGNITUDE); PCHKERRQ(eps,ierr);
200 ierr = EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL); PCHKERRQ(eps,ierr);
203 ierr = EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL); PCHKERRQ(eps,ierr);
206 ierr = EPSSetWhichEigenpairs(eps,EPS_LARGEST_IMAGINARY); PCHKERRQ(eps,ierr);
209 ierr = EPSSetWhichEigenpairs(eps,EPS_SMALLEST_IMAGINARY); PCHKERRQ(eps,ierr);
212 ierr = EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE); PCHKERRQ(eps,ierr);
215 ierr = EPSSetWhichEigenpairs(eps,EPS_TARGET_REAL); PCHKERRQ(eps,ierr);
218 MFEM_ABORT(
"Which eigenpair not implemented!");
225 ierr = EPSSetTarget(eps,target); PCHKERRQ(eps,ierr);
232 ierr = EPSGetST(eps,&st); PCHKERRQ(eps,ierr);
233 switch (transformation)
236 ierr = STSetType(st,STSHIFT); PCHKERRQ(eps,ierr);
239 ierr = STSetType(st,STSINVERT); PCHKERRQ(eps,ierr);
242 MFEM_ABORT(
"Spectral transformation not implemented!");
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.
void Customize(bool customize=true) const
Customize object with options set.
void MFEMInitializeSlepc()
Wrapper for PETSc's matrix class.
void GetEigenvalue(unsigned int i, double &lr) const
Get the corresponding eigenvalue.
int Size() const
Returns the size of the vector.
double * GetData() const
Return a pointer to the beginning of the Vector data.
void SetTarget(double target)
Wrapper for PETSc's vector class.
void ResetArray()
Reset the PETSc Vec object to use its default data. Call this method after the use of PlaceArray()...
void SetSpectralTransformation(SpectralTransformation transformation)
virtual ~SlepcEigenSolver()
void transformation(const Vector &p, Vector &v)
void SetMaxIter(int max_iter)
Set maximum number of iterations.
void SetNumModes(int num_eigs)
Set the number of required eigenmodes.
void Solve()
Solve the eigenvalue problem for the specified number of eigenvalues.
void SetTol(double tol)
Set solver tolerance.
SlepcEigenSolver(MPI_Comm comm, const std::string &prefix=std::string())
Constructors.
int GetNumConverged()
Get the number of converged eigenvalues.
void SetWhichEigenpairs(Which which)
struct _p_PetscObject * PetscObject
void PlaceArray(PetscScalar *temp_data)
Temporarily replace the data of the PETSc Vec object. To return to the original data array...
void GetEigenvector(unsigned int i, Vector &vr) const
Get the corresponding eigenvector.
void SetOperators(const PetscParMatrix &op, const PetscParMatrix &opB)
Set operator for generalized eigenvalue problem.