virtual void set_solver_type (EPS &eps) const;
};
+
+ /**
+ * An implementation of the solver interface using the SLEPc LAPACK
+ * direct solver.
+ *
+ * @ingroup SLEPcWrappers
+ *
+ * @author Toby D. Young 2013
+ */
+ class SolverLAPACK : public SolverBase
+ {
+ public:
+
+ /**
+ * Standardized data struct to pipe additional data to the solver,
+ * should it be needed.
+ */
+ struct AdditionalData
+ {};
+
+ /**
+ * SLEPc solvers will want to have an MPI communicator context
+ * over which computations are parallelized. By default, this
+ * carries the same behaviour has the PETScWrappers, but you can
+ * change that.
+ */
+ SolverLAPACK (SolverControl &cn,
+ const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
+ const AdditionalData &data = AdditionalData());
+
+ protected:
+
+ /**
+ * Store a copy of the flags for this particular solver.
+ */
+ const AdditionalData additional_data;
+
+ /**
+ * Function that takes a Eigenvalue Problem Solver context object,
+ * and sets the type of solver that is appropriate for this class.
+ */
+ virtual void set_solver_type (EPS &eps) const;
+ };
+
// --------------------------- inline and template functions -----------
/**
* This is declared here to make it possible to take a std::vector
#endif
}
+
+ /* ---------------------- LAPACK ------------------------- */
+ SolverLAPACK::SolverLAPACK (SolverControl &cn,
+ const MPI_Comm &mpi_communicator,
+ const AdditionalData &data)
+ :
+ SolverBase (cn, mpi_communicator),
+ additional_data (data)
+ {}
+
+ void
+ SolverLAPACK::set_solver_type (EPS &eps) const
+ {
+ // 'Tis overwhelmingly likely that PETSc/SLEPc *always* has
+ // BLAS/LAPACK, but let's be defensive.
+#if PETSC_HAVE_BLASLAPACK
+ int ierr;
+ ierr = EPSSetType (eps, const_cast<char *>(EPSLAPACK));
+ AssertThrow (ierr == 0, ExcSLEPcError(ierr));
+
+ // hand over the absolute tolerance and the maximum number of
+ // iteration steps to the SLEPc convergence criterion.
+ ierr = EPSSetTolerances (eps, this->solver_control.tolerance(),
+ this->solver_control.max_steps());
+ AssertThrow (ierr == 0, ExcSLEPcError(ierr));
+#else
+ // Supress compiler warnings about unused paameters.
+ (void) eps;
+
+ Assert ((false),
+ ExcMessage ("Your PETSc/SLEPc installation was not configured with BLAS/LAPACK "
+ "but this is needed to use the LAPACK solver."));
+#endif
+ }
+
}
DEAL_II_NAMESPACE_CLOSE