#if DEAL_II_PETSC_VERSION_GTE(3,1,0)
/**
- * An implementation of the solver interface using the SLEPc Davidson
- * solver. Usage (incomplete/untested): All problem types.
+ * An implementation of the solver interface using the SLEPc
+ * Generalized Davidson solver. Usage: All problem types.
+ *
+ * @deprecated The use of this class has been superceded by the
+ * class SolverGeneralizedDavidson.
*
* @ingroup SLEPcWrappers
- * @author Toby D. Young 2010
+ * @author Toby D. Young 2011
+ */
+#define SolverDavidson SolverGeneralizedDavidson;
+
+ /**
+ * An implementation of the solver interface using the SLEPc
+ * Generalized Davidson solver. Usage: All problem types.
+ *
+ * @ingroup SLEPcWrappers
+ * @author Toby D. Young 2013
+ */
+ class SolverGeneralizedDavidson : 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.
+ */
+ SolverGeneralizedDavidson (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;
+
+ };
+
+ /**
+ * An implementation of the solver interface using the SLEPc
+ * Jacobi-Davidson solver. Usage: All problem types.
+ *
+ * @ingroup SLEPcWrappers
+ * @author Toby D. Young 2013
*/
- class SolverDavidson : public SolverBase
+ class SolverJacobiDavidson : public SolverBase
{
public:
/**
* the PETScWrappers, but you can
* change that.
*/
- SolverDavidson (SolverControl &cn,
- const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
- const AdditionalData &data = AdditionalData());
+ SolverJacobiDavidson (SolverControl &cn,
+ const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
+ const AdditionalData &data = AdditionalData());
protected:
};
#endif
+
+
// --------------------------- inline and template functions -----------
/**
* This is declared here to make it possible to take a std::vector
}
#if DEAL_II_PETSC_VERSION_GTE(3,1,0)
- /* ---------------------- Davidson ----------------------- */
- SolverDavidson::SolverDavidson (SolverControl &cn,
- const MPI_Comm &mpi_communicator,
- const AdditionalData &data)
+ /* ---------------- Generalized Davidson ----------------- */
+ SolverGeneralizedDavidson::SolverGeneralizedDavidson (SolverControl &cn,
+ const MPI_Comm &mpi_communicator,
+ const AdditionalData &data)
:
SolverBase (cn, mpi_communicator),
additional_data (data)
{}
void
- SolverDavidson::set_solver_type (EPS &eps) const
+ SolverGeneralizedDavidson::set_solver_type (EPS &eps) const
{
int ierr;
ierr = EPSSetType (eps, const_cast<char *>(EPSGD));
this->solver_control.max_steps());
AssertThrow (ierr == 0, ExcSLEPcError(ierr));
}
+
+ /* ------------------ Jacobi Davidson -------------------- */
+
+ SolverJacobiDavidson::SolverJacobiDavidson (SolverControl &cn,
+ const MPI_Comm &mpi_communicator,
+ const AdditionalData &data)
+ :
+ SolverBase (cn, mpi_communicator),
+ additional_data (data)
+ {}
+
+ void
+ SolverJacobiDavidson::set_solver_type (EPS &eps) const
+ {
+ int ierr;
+ ierr = EPSSetType (eps, const_cast<char *>(EPSJD));
+ 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));
+ }
#endif
}