]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
Distinguish between Generalized and Jacobi Davidson solvers. Keep backward comparibil...
authoryoung <young@0785d39b-7218-0410-832d-ea1e28bc413d>
Wed, 9 Jan 2013 20:52:56 +0000 (20:52 +0000)
committeryoung <young@0785d39b-7218-0410-832d-ea1e28bc413d>
Wed, 9 Jan 2013 20:52:56 +0000 (20:52 +0000)
git-svn-id: https://svn.dealii.org/trunk@28000 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/include/deal.II/lac/slepc_solver.h
deal.II/source/lac/slepc_solver.cc

index 6aa1c52cb9af8be8d21f8f7c15bd12d9cf2bd4eb..ed85982912b778342f796d5415c32904a18e9928 100644 (file)
@@ -619,13 +619,74 @@ namespace SLEPcWrappers
 
 #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:
     /**
@@ -645,9 +706,9 @@ namespace SLEPcWrappers
      * 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:
 
@@ -668,6 +729,8 @@ namespace SLEPcWrappers
   };
 #endif
 
+
+
   // --------------------------- inline and template functions -----------
   /**
    * This is declared here to make it possible to take a std::vector
index a0a5c6e865d4ab6b574237e9750e70b450a91bd8..85345bd9a1dcdcda8b47f15233ac63c2347db2d4 100644 (file)
@@ -377,18 +377,18 @@ namespace SLEPcWrappers
   }
 
 #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));
@@ -403,6 +403,33 @@ namespace SLEPcWrappers
                              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
 
 }

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.