22static PetscErrorCode ierr;
45 ierr = PetscOptionsSetValue(NULL,
"-cuda_device",
47 MFEM_VERIFY(!ierr,
"Unable to set initial option value to PETSc");
49 ierr = SlepcInitialize(argc,argv,rc_file,help);
50 MFEM_VERIFY(!ierr,
"Unable to initialize SLEPc");
55 ierr = SlepcFinalize();
56 MFEM_VERIFY(!ierr,
"Unable to finalize SLEPc");
66 ierr = EPSCreate(comm,&eps); CCHKERRQ(comm,ierr);
67 ierr = EPSSetOptionsPrefix(eps, prefix.c_str()); PCHKERRQ(eps, ierr);
76 ierr = PetscObjectGetComm((
PetscObject)eps,&comm); PCHKERRQ(eps,ierr);
77 ierr = EPSDestroy(&eps); CCHKERRQ(comm,ierr);
87 ierr = EPSSetOperators(eps,op,NULL); PCHKERRQ(eps, ierr);
100 ierr = EPSSetOperators(eps,op,opB); PCHKERRQ(eps,ierr);
110 ierr = EPSGetTolerances(eps,NULL,&max_its); PCHKERRQ(eps,ierr);
112 if (max_its==0) { max_its = PETSC_DECIDE; }
113 ierr = EPSSetTolerances(eps,tol,max_its); PCHKERRQ(eps,ierr);
120 ierr = EPSGetTolerances(eps,&tol,NULL); PCHKERRQ(eps,ierr);
121 ierr = EPSSetTolerances(eps,tol,max_its); PCHKERRQ(eps,ierr);
126 ierr = EPSSetDimensions(eps,num_eigs,PETSC_DECIDE,PETSC_DECIDE);
134 ierr = EPSSolve(eps); PCHKERRQ(eps,ierr);
139 if (!customize) {clcustom =
true; }
142 ierr = EPSSetFromOptions(eps); PCHKERRQ(eps,ierr);
149 ierr = EPSGetEigenvalue(eps,i,&lr,NULL); PCHKERRQ(eps,ierr);
155 ierr = EPSGetEigenvalue(eps,i,&lr,&lc); PCHKERRQ(eps,ierr);
160 MFEM_VERIFY(VR,
"Missing real vector");
162 MFEM_ASSERT(vr.
Size() == VR->
Size(),
"invalid vr.Size() = " << vr.
Size()
163 <<
", expected size = " << VR->
Size());
166 ierr = EPSGetEigenvector(eps,i,*VR,NULL); PCHKERRQ(eps,ierr);
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());
183 ierr = EPSGetEigenvector(eps,i,*VR,*VC); PCHKERRQ(eps,ierr);
191 ierr = EPSGetConverged(eps,&num_conv); PCHKERRQ(eps,ierr);
192 return static_cast<int>(num_conv);
200 ierr = EPSSetWhichEigenpairs(eps,EPS_LARGEST_MAGNITUDE); PCHKERRQ(eps,ierr);
203 ierr = EPSSetWhichEigenpairs(eps,EPS_SMALLEST_MAGNITUDE); PCHKERRQ(eps,ierr);
206 ierr = EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL); PCHKERRQ(eps,ierr);
209 ierr = EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL); PCHKERRQ(eps,ierr);
212 ierr = EPSSetWhichEigenpairs(eps,EPS_LARGEST_IMAGINARY); PCHKERRQ(eps,ierr);
215 ierr = EPSSetWhichEigenpairs(eps,EPS_SMALLEST_IMAGINARY); PCHKERRQ(eps,ierr);
218 ierr = EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE); PCHKERRQ(eps,ierr);
221 ierr = EPSSetWhichEigenpairs(eps,EPS_TARGET_REAL); PCHKERRQ(eps,ierr);
224 MFEM_ABORT(
"Which eigenpair not implemented!");
231 ierr = EPSSetTarget(eps,target); PCHKERRQ(eps,ierr);
238 ierr = EPSGetST(eps,&st); PCHKERRQ(eps,ierr);
242 ierr = STSetType(st,STSHIFT); PCHKERRQ(eps,ierr);
245 ierr = STSetType(st,STSINVERT); PCHKERRQ(eps,ierr);
248 MFEM_ABORT(
"Spectral transformation not implemented!");
static bool Allows(unsigned long b_mask)
Return true if any of the backends in the backend mask, b_mask, are allowed.
static int GetId()
Get the device id of the configured device.
Wrapper for PETSc's matrix class.
void ResetMemory()
Completes the operation started with PlaceMemory.
void PlaceMemory(Memory< real_t > &, bool=false)
This requests write access from where the memory is valid and temporarily replaces the corresponding ...
void PlaceArray(PetscScalar *temp_data)
Temporarily replace the data of the PETSc Vec object. To return to the original data array,...
Which
Target spectrum for the eigensolver.
@ LARGEST_REAL
The eigenvalues with the largest real component.
@ SMALLEST_MAGNITUDE
The eigenvalues with the smallest complex magnitude.
@ SMALLEST_REAL
The eigenvalues with the smallest real component.
@ TARGET_MAGNITUDE
The eigenvalues with complex magnitude closest to the target value.
@ LARGEST_MAGNITUDE
The eigenvalues with the largest complex magnitude (default)
@ TARGET_REAL
The eigenvalues with the real component closest to the target value.
@ SMALLEST_IMAGINARY
The eigenvalues with the smallest imaginary component.
@ LARGEST_IMAGINARY
The eigenvalues with the largest imaginary component.
void Customize(bool customize=true) const
Customize object with options set.
void SetOperators(const PetscParMatrix &op, const PetscParMatrix &opB)
Set operators for generalized eigenvalue problem.
void SetWhichEigenpairs(Which which)
Set the which eigenvalues the solver will target and the order they will be indexed in.
SlepcEigenSolver(MPI_Comm comm, const std::string &prefix=std::string())
Constructors.
int GetNumConverged()
Get the number of converged eigenvalues after the call to SlepcEigenSolver::Solve.
void SetTarget(real_t target)
Set the target value for the eigenpairs you want when using SlepcEigenSolver::TARGET_MAGNITUDE or Sle...
void GetEigenvalue(unsigned int i, real_t &lr) const
Get the ith eigenvalue after the system has been solved.
void SetNumModes(int num_eigs)
Set the number of eigenmodes to compute.
virtual ~SlepcEigenSolver()
void SetSpectralTransformation(SpectralTransformation transformation)
Set the spectral transformation strategy for acceletating convergenvce. Both SlepcEigenSolver::SHIFT ...
void SetMaxIter(int max_iter)
Set maximum number of iterations allowed in the call to SlepcEigenSolver::Solve.
SpectralTransformation
Spectral transformations that can be used by the solver in order to accelerate the convergence to the...
@ SHIFT
Utilize the shift of origin strategy.
@ SHIFT_INVERT
Utilize the shift and invert strategy.
void GetEigenvector(unsigned int i, Vector &vr) const
Get the ith eigenvector after the system has been solved.
void SetTol(real_t tol)
Set solver convergence tolerance relative to the magnitude of the eigenvalue.
void SetOperator(const PetscParMatrix &op)
Set operator for standard eigenvalue problem.
void Solve()
Solve the eigenvalue problem for the specified number of eigenvalues.
Memory< real_t > & GetMemory()
Return a reference to the Memory object used by the Vector.
int Size() const
Returns the size of the vector.
void transformation(const Vector &p, Vector &v)
void MFEMInitializeSlepc()
struct _p_PetscObject * PetscObject
@ CUDA_MASK
Biwise-OR of all CUDA backends.