]> https://gitweb.dealii.org/ - dealii.git/commitdiff
add set_initial_space() to SLEPc solvers
authorDenis Davydov <davydden@gmail.com>
Mon, 9 Nov 2015 14:11:15 +0000 (15:11 +0100)
committerDenis Davydov <davydden@gmail.com>
Mon, 9 Nov 2015 16:38:18 +0000 (17:38 +0100)
include/deal.II/lac/slepc_solver.h
source/lac/slepc_solver.cc

index 6fac9094b48ff441737dd115a269ce78f1b7adf7..0e6770a79f72c29c5034ff68cfce6993b78093fb 100644 (file)
@@ -186,7 +186,15 @@ namespace SLEPcWrappers
      */
     void
     set_initial_vector
-    (const PETScWrappers::VectorBase &this_initial_vector);
+    (const PETScWrappers::VectorBase &this_initial_vector) DEAL_II_DEPRECATED;
+
+    /**
+     * Set the initial vector space for the solver.
+     */
+    template <typename Vector>
+    void
+    set_initial_space
+    (const std::vector<Vector> &initial_vectors);
 
     /**
      * Set the spectral transformation to be used.
@@ -361,9 +369,9 @@ namespace SLEPcWrappers
     const PETScWrappers::MatrixBase *opB;
 
     /**
-     * An initial vector used to "feed" some SLEPc solvers.
+     * An initial vector space used in SLEPc solvers.
      */
-    const PETScWrappers::VectorBase *initial_vector;
+    std::vector<Vec> initial_vectors;
 
     /**
      * Pointer to an an object that describes transformations that can be
@@ -844,6 +852,17 @@ namespace SLEPcWrappers
                      real_eigenvectors[index], imag_eigenvectors[index]);
   }
 
+  template <typename Vector>
+  void
+  SolverBase::set_initial_space(const std::vector<Vector> &vectors)
+  {
+    initial_vectors.resize(vectors.size());
+    for (unsigned int i = 0; i < initial_vectors.size(); i++)
+      {
+        initial_vectors[i] = vectors[i];
+      }
+  }
+
 }
 
 DEAL_II_NAMESPACE_CLOSE
index 3953ab2bad32d909892aa2244954eb5dc04ce2e0..c442b180559db0911f1a96e0a0e304bc3e4983d1 100644 (file)
@@ -53,7 +53,6 @@ namespace SLEPcWrappers
     set_problem (EPS_GNHEP),
     opA (NULL),
     opB (NULL),
-    initial_vector (NULL),
     transformation (NULL)
   {}
 
@@ -80,7 +79,8 @@ namespace SLEPcWrappers
   void
   SolverBase::set_initial_vector (const PETScWrappers::VectorBase &this_initial_vector)
   {
-    initial_vector = (&this_initial_vector);
+    initial_vectors.resize(1);
+    initial_vectors[0] = this_initial_vector;
   }
 
   void
@@ -140,14 +140,13 @@ namespace SLEPcWrappers
     AssertThrow (ierr == 0, ExcSLEPcError(ierr));
 
     // set the initial vector(s) if any
-    if (initial_vector && initial_vector->size() != 0)
+    if (initial_vectors.size() != 0)
       {
 
 #if DEAL_II_PETSC_VERSION_LT(3,1,0)
-        ierr = EPSSetInitialVector (solver_data->eps, *initial_vector);
+        ierr = EPSSetInitialVector (solver_data->eps, &initial_vectors[0]);
 #else
-        Vec this_vector = *initial_vector;
-        ierr = EPSSetInitialSpace (solver_data->eps, 1, &this_vector);
+        ierr = EPSSetInitialSpace (solver_data->eps, initial_vectors.size(), &initial_vectors[0]);
 #endif
 
         AssertThrow (ierr == 0, ExcSLEPcError(ierr));

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.