]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Build in assertions...
authorToby D. Young <tyoung@ippt.pan.pl>
Thu, 2 May 2013 20:18:40 +0000 (20:18 +0000)
committerToby D. Young <tyoung@ippt.pan.pl>
Thu, 2 May 2013 20:18:40 +0000 (20:18 +0000)
git-svn-id: https://svn.dealii.org/trunk@29432 0785d39b-7218-0410-832d-ea1e28bc413d

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

index 595843a5d08320d6f17bd17fda43d599d5bab77a..402f4ebf69d4cfd42d516306f0867b98ed776658 100644 (file)
@@ -151,9 +151,9 @@ namespace SLEPcWrappers
     template <typename OutputVector>
     void
     solve (const PETScWrappers::MatrixBase &A,
-           std::vector<double>             &kr,
-           std::vector<OutputVector>       &vr,
-           const unsigned int              n_eigenvectors = 1);
+           std::vector<double>             &r_eigenvalues,
+           std::vector<OutputVector>       &r_eigenvectors = std::vector<OutputVector> (),
+           const unsigned int               n_eigenvectors = 1);
 
     /**
      * Same as above, but here a composite method for solving the
@@ -163,8 +163,8 @@ namespace SLEPcWrappers
     void
     solve (const PETScWrappers::MatrixBase &A,
            const PETScWrappers::MatrixBase &B,
-           std::vector<double>             &kr,
-           std::vector<OutputVector>       &vr,
+           std::vector<double>             &r_eigenvalues,
+           std::vector<OutputVector>       &r_eigenvectors = std::vector<OutputVector> (),
            const unsigned int               n_eigenvectors = 1);
 
     /**
@@ -288,12 +288,8 @@ namespace SLEPcWrappers
      * \text{n\_converged}-1$.
      */
     void
-    get_eigenpair (const unsigned int index,
-#ifndef PETSC_USE_COMPLEX
+    get_eigenpair (const unsigned int         index,
                    double                    &kr,
-#else
-                   std::complex<double>      &kr,
-#endif
                    PETScWrappers::VectorBase &vr);
 
     /**
@@ -723,15 +719,15 @@ namespace SLEPcWrappers
                      std::vector<OutputVector>       &vr,
                      const unsigned int               n_eigenvectors)
   {
-    // Panic if no eigenpairs are wanted.
-    AssertThrow (n_eigenvectors != 0, ExcSLEPcWrappersUsageError());
-
-    unsigned int n_converged = 0;
+    // Panic if the number of eigenpairs wanted is out of bounds.
+    AssertThrow ((n_eigenvectors > 0) && (n_eigenvectors <= A.m ()), 
+                ExcSLEPcWrappersUsageError());
 
     // Set the matrices of the problem
     set_matrices (A); 
 
     // and solve
+    unsigned int n_converged = 0;
     solve (n_eigenvectors, &n_converged);
 
     if (n_converged > n_eigenvectors)
@@ -755,15 +751,19 @@ namespace SLEPcWrappers
                       std::vector<OutputVector>       &vr,
                       const unsigned int               n_eigenvectors)
   {
-    // Panic if no eigenpairs are wanted.
-    AssertThrow (n_eigenvectors != 0, ExcSLEPcWrappersUsageError());
+    // Guard against incompatible matrix sizes:
+    AssertThrow (A.m() == B.m (), ExcDimensionMismatch(A.m(), B.m()));
+    AssertThrow (A.n() == B.n (), ExcDimensionMismatch(A.n(), B.n()));
 
-    unsigned int n_converged = 0;
+    // Panic if the number of eigenpairs wanted is out of bounds.
+    AssertThrow ((n_eigenvectors > 0) && (n_eigenvectors <= A.m ()), 
+                ExcSLEPcWrappersUsageError());
 
     // Set the matrices of the problem
     set_matrices (A, B);
 
     // and solve
+    unsigned int n_converged = 0;
     solve (n_eigenvectors, &n_converged);
 
     if (n_converged >= n_eigenvectors)
index ded535c6bcf72acd45bb76dc659356ec50f06a43..9d5877f2aaa7eb8b1eb79f0482c912162952bb37 100644 (file)
@@ -220,11 +220,7 @@ namespace SLEPcWrappers
 
   void
   SolverBase::get_eigenpair (const unsigned int         index,
-#ifndef PETSC_USE_COMPLEX
                              double                    &kr,
-#else
-                             std::complex<double>      &kr,
-#endif
                              PETScWrappers::VectorBase &vr)
   {
     AssertThrow (solver_data.get() != 0, ExcSLEPcWrappersUsageError());

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.