]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Fixes in relation to #1673. 1784/head
authorJean-Paul Pelteret <jppelteret@gmail.com>
Mon, 26 Oct 2015 21:11:29 +0000 (22:11 +0100)
committerJean-Paul Pelteret <jppelteret@gmail.com>
Mon, 26 Oct 2015 21:11:29 +0000 (22:11 +0100)
Changed size_type in PreconditionChebyshev. Added static_asserts that
compare matrix and vector types in precondtioners. Updated
linear_operator_08 test to account for inability to wrap
IterativeInverse in a linear_operator.

include/deal.II/lac/precondition.h
tests/lac/linear_operator_08.cc
tests/lac/linear_operator_08.with_cxx11=on.output

index 536d4d097faf40cce558017dda0928e2e92f9ad2..6c0285b23aad96d242784661ea2e1e714ca882d1 100644 (file)
@@ -810,7 +810,7 @@ public:
   /**
    * Declare type for container size.
    */
-  typedef typename MATRIX::size_type size_type;
+  typedef types::global_dof_index size_type;
 
   /**
    * Standardized data struct to pipe additional parameters to the
@@ -1092,6 +1092,12 @@ template<class VECTOR>
 inline void
 PreconditionRichardson::vmult (VECTOR &dst, const VECTOR &src) const
 {
+#ifdef DEAL_II_WITH_CXX11
+  static_assert(
+    std::is_same<size_type, typename VECTOR::size_type>::value,
+    "PreconditionRichardson and VECTOR must have the same size_type.");
+#endif // DEAL_II_WITH_CXX11
+
   dst.equ(relaxation,src);
 }
 
@@ -1101,6 +1107,12 @@ template<class VECTOR>
 inline void
 PreconditionRichardson::Tvmult (VECTOR &dst, const VECTOR &src) const
 {
+#ifdef DEAL_II_WITH_CXX11
+  static_assert(
+    std::is_same<size_type, typename VECTOR::size_type>::value,
+    "PreconditionRichardson and VECTOR must have the same size_type.");
+#endif // DEAL_II_WITH_CXX11
+
   dst.equ(relaxation,src);
 }
 
@@ -1108,6 +1120,12 @@ template<class VECTOR>
 inline void
 PreconditionRichardson::vmult_add (VECTOR &dst, const VECTOR &src) const
 {
+#ifdef DEAL_II_WITH_CXX11
+  static_assert(
+    std::is_same<size_type, typename VECTOR::size_type>::value,
+    "PreconditionRichardson and VECTOR must have the same size_type.");
+#endif // DEAL_II_WITH_CXX11
+
   dst.add(relaxation,src);
 }
 
@@ -1117,6 +1135,12 @@ template<class VECTOR>
 inline void
 PreconditionRichardson::Tvmult_add (VECTOR &dst, const VECTOR &src) const
 {
+#ifdef DEAL_II_WITH_CXX11
+  static_assert(
+    std::is_same<size_type, typename VECTOR::size_type>::value,
+    "PreconditionRichardson and VECTOR must have the same size_type.");
+#endif // DEAL_II_WITH_CXX11
+
   dst.add(relaxation,src);
 }
 
@@ -1176,6 +1200,12 @@ template<class VECTOR>
 inline void
 PreconditionJacobi<MATRIX>::vmult (VECTOR &dst, const VECTOR &src) const
 {
+#ifdef DEAL_II_WITH_CXX11
+  static_assert(
+    std::is_same<typename PreconditionJacobi<MATRIX>::size_type, typename VECTOR::size_type>::value,
+    "PreconditionJacobi and VECTOR must have the same size_type.");
+#endif // DEAL_II_WITH_CXX11
+
   Assert (this->A!=0, ExcNotInitialized());
   this->A->precondition_Jacobi (dst, src, this->relaxation);
 }
@@ -1187,6 +1217,12 @@ template<class VECTOR>
 inline void
 PreconditionJacobi<MATRIX>::Tvmult (VECTOR &dst, const VECTOR &src) const
 {
+#ifdef DEAL_II_WITH_CXX11
+  static_assert(
+    std::is_same<typename PreconditionJacobi<MATRIX>::size_type, typename VECTOR::size_type>::value,
+    "PreconditionJacobi and VECTOR must have the same size_type.");
+#endif // DEAL_II_WITH_CXX11
+
   Assert (this->A!=0, ExcNotInitialized());
   this->A->precondition_Jacobi (dst, src, this->relaxation);
 }
@@ -1198,6 +1234,12 @@ template<class VECTOR>
 inline void
 PreconditionJacobi<MATRIX>::step (VECTOR &dst, const VECTOR &src) const
 {
+#ifdef DEAL_II_WITH_CXX11
+  static_assert(
+    std::is_same<typename PreconditionJacobi<MATRIX>::size_type, typename VECTOR::size_type>::value,
+    "PreconditionJacobi and VECTOR must have the same size_type.");
+#endif // DEAL_II_WITH_CXX11
+
   Assert (this->A!=0, ExcNotInitialized());
   this->A->Jacobi_step (dst, src, this->relaxation);
 }
@@ -1209,6 +1251,12 @@ template<class VECTOR>
 inline void
 PreconditionJacobi<MATRIX>::Tstep (VECTOR &dst, const VECTOR &src) const
 {
+#ifdef DEAL_II_WITH_CXX11
+  static_assert(
+    std::is_same<typename PreconditionJacobi<MATRIX>::size_type, typename VECTOR::size_type>::value,
+    "PreconditionJacobi and VECTOR must have the same size_type.");
+#endif // DEAL_II_WITH_CXX11
+
   step (dst, src);
 }
 
@@ -1221,6 +1269,12 @@ template<class VECTOR>
 inline void
 PreconditionSOR<MATRIX>::vmult (VECTOR &dst, const VECTOR &src) const
 {
+#ifdef DEAL_II_WITH_CXX11
+  static_assert(
+    std::is_same<typename PreconditionSOR<MATRIX>::size_type, typename VECTOR::size_type>::value,
+    "PreconditionSOR and VECTOR must have the same size_type.");
+#endif // DEAL_II_WITH_CXX11
+
   Assert (this->A!=0, ExcNotInitialized());
   this->A->precondition_SOR (dst, src, this->relaxation);
 }
@@ -1232,6 +1286,12 @@ template<class VECTOR>
 inline void
 PreconditionSOR<MATRIX>::Tvmult (VECTOR &dst, const VECTOR &src) const
 {
+#ifdef DEAL_II_WITH_CXX11
+  static_assert(
+    std::is_same<typename PreconditionSOR<MATRIX>::size_type, typename VECTOR::size_type>::value,
+    "PreconditionSOR and VECTOR must have the same size_type.");
+#endif // DEAL_II_WITH_CXX11
+
   Assert (this->A!=0, ExcNotInitialized());
   this->A->precondition_TSOR (dst, src, this->relaxation);
 }
@@ -1243,6 +1303,12 @@ template<class VECTOR>
 inline void
 PreconditionSOR<MATRIX>::step (VECTOR &dst, const VECTOR &src) const
 {
+#ifdef DEAL_II_WITH_CXX11
+  static_assert(
+    std::is_same<typename PreconditionSOR<MATRIX>::size_type, typename VECTOR::size_type>::value,
+    "PreconditionSOR and VECTOR must have the same size_type.");
+#endif // DEAL_II_WITH_CXX11
+
   Assert (this->A!=0, ExcNotInitialized());
   this->A->SOR_step (dst, src, this->relaxation);
 }
@@ -1254,6 +1320,12 @@ template<class VECTOR>
 inline void
 PreconditionSOR<MATRIX>::Tstep (VECTOR &dst, const VECTOR &src) const
 {
+#ifdef DEAL_II_WITH_CXX11
+  static_assert(
+    std::is_same<typename PreconditionSOR<MATRIX>::size_type, typename VECTOR::size_type>::value,
+    "PreconditionSOR and VECTOR must have the same size_type.");
+#endif // DEAL_II_WITH_CXX11
+
   Assert (this->A!=0, ExcNotInitialized());
   this->A->TSOR_step (dst, src, this->relaxation);
 }
@@ -1301,6 +1373,12 @@ template<class VECTOR>
 inline void
 PreconditionSSOR<MATRIX>::vmult (VECTOR &dst, const VECTOR &src) const
 {
+#ifdef DEAL_II_WITH_CXX11
+  static_assert(
+    std::is_same<typename PreconditionSSOR<MATRIX>::size_type, typename VECTOR::size_type>::value,
+    "PreconditionSSOR and VECTOR must have the same size_type.");
+#endif // DEAL_II_WITH_CXX11
+
   Assert (this->A!=0, ExcNotInitialized());
   this->A->precondition_SSOR (dst, src, this->relaxation, pos_right_of_diagonal);
 }
@@ -1312,6 +1390,12 @@ template<class VECTOR>
 inline void
 PreconditionSSOR<MATRIX>::Tvmult (VECTOR &dst, const VECTOR &src) const
 {
+#ifdef DEAL_II_WITH_CXX11
+  static_assert(
+    std::is_same<typename PreconditionSSOR<MATRIX>::size_type, typename VECTOR::size_type>::value,
+    "PreconditionSSOR and VECTOR must have the same size_type.");
+#endif // DEAL_II_WITH_CXX11
+
   Assert (this->A!=0, ExcNotInitialized());
   this->A->precondition_SSOR (dst, src, this->relaxation, pos_right_of_diagonal);
 }
@@ -1323,6 +1407,12 @@ template<class VECTOR>
 inline void
 PreconditionSSOR<MATRIX>::step (VECTOR &dst, const VECTOR &src) const
 {
+#ifdef DEAL_II_WITH_CXX11
+  static_assert(
+    std::is_same<typename PreconditionSSOR<MATRIX>::size_type, typename VECTOR::size_type>::value,
+    "PreconditionSSOR and VECTOR must have the same size_type.");
+#endif // DEAL_II_WITH_CXX11
+
   Assert (this->A!=0, ExcNotInitialized());
   this->A->SSOR_step (dst, src, this->relaxation);
 }
@@ -1334,6 +1424,12 @@ template<class VECTOR>
 inline void
 PreconditionSSOR<MATRIX>::Tstep (VECTOR &dst, const VECTOR &src) const
 {
+#ifdef DEAL_II_WITH_CXX11
+  static_assert(
+    std::is_same<typename PreconditionSSOR<MATRIX>::size_type, typename VECTOR::size_type>::value,
+    "PreconditionSSOR and VECTOR must have the same size_type.");
+#endif // DEAL_II_WITH_CXX11
+
   step (dst, src);
 }
 
@@ -1373,6 +1469,12 @@ template<class VECTOR>
 inline void
 PreconditionPSOR<MATRIX>::vmult (VECTOR &dst, const VECTOR &src) const
 {
+#ifdef DEAL_II_WITH_CXX11
+  static_assert(
+    std::is_same<typename PreconditionPSOR<MATRIX>::size_type, typename VECTOR::size_type>::value,
+    "PreconditionPSOR and VECTOR must have the same size_type.");
+#endif // DEAL_II_WITH_CXX11
+
   Assert (this->A!=0, ExcNotInitialized());
   dst = src;
   this->A->PSOR (dst, *permutation, *inverse_permutation, this->relaxation);
@@ -1385,6 +1487,12 @@ template<class VECTOR>
 inline void
 PreconditionPSOR<MATRIX>::Tvmult (VECTOR &dst, const VECTOR &src) const
 {
+#ifdef DEAL_II_WITH_CXX11
+  static_assert(
+    std::is_same<typename PreconditionPSOR<MATRIX>::size_type, typename VECTOR::size_type>::value,
+    "PreconditionPSOR and VECTOR must have the same size_type.");
+#endif // DEAL_II_WITH_CXX11
+
   Assert (this->A!=0, ExcNotInitialized());
   dst = src;
   this->A->TPSOR (dst, *permutation, *inverse_permutation, this->relaxation);
@@ -1678,7 +1786,13 @@ inline
 PreconditionChebyshev<MATRIX,VECTOR>::PreconditionChebyshev ()
   :
   is_initialized (false)
-{}
+{
+#ifdef DEAL_II_WITH_CXX11
+  static_assert(
+    std::is_same<size_type, typename VECTOR::size_type>::value,
+    "PreconditionChebyshev and VECTOR must have the same size_type.");
+#endif // DEAL_II_WITH_CXX11
+}
 
 
 
index 965fda2f21206eb24c30eb650acd40396401817a..6ba3f68e7de732d4e225fb06fa37e8ae60bef806 100644 (file)
@@ -396,21 +396,21 @@ int main()
       const auto lo_A_inv = linear_operator(solver);
       const Vector<double> x = lo_A_inv*b;
     }
-    {
-      deallog << "IterativeInverse" << std::endl;
-
-      PreconditionJacobi< SparseMatrix<double> > preconditioner;
-      preconditioner.initialize(A);
-
-      ReductionControl solver_control (10, 1.e-30, 1.e-2);
-      IterativeInverse< Vector<double> > A_inv;
-      A_inv.initialize(A,preconditioner);
-      A_inv.solver.select("cg");
-      A_inv.solver.set_control(solver_control);
-
-      const auto lo_A_inv = linear_operator(A_inv);
-      const Vector<double> x = lo_A_inv*b;
-    }
+//    { // See #1673 and #1784
+//      deallog << "IterativeInverse" << std::endl;
+//
+//      PreconditionJacobi< SparseMatrix<double> > preconditioner;
+//      preconditioner.initialize(A);
+//
+//      ReductionControl solver_control (10, 1.e-30, 1.e-2);
+//      IterativeInverse< Vector<double> > A_inv;
+//      A_inv.initialize(A,preconditioner);
+//      A_inv.solver.select("cg");
+//      A_inv.solver.set_control(solver_control);
+//
+//      const auto lo_A_inv = linear_operator(A_inv);
+//      const Vector<double> x = lo_A_inv*b;
+//    }
     deallog.pop();
 
 
index 0aeb146b8214355c6251d9eb48c82a36f92c62f5..459bc5bca8d87e695c2138450ef91c94c536969a 100644 (file)
@@ -80,9 +80,6 @@ DEAL:Solvers:Linear operator:Richardson::Convergence step 1 value 0.000000000
 DEAL:Solvers::SolverSelector
 DEAL:Solvers:cg::Convergence step 0 value 0.000000000
 DEAL:Solvers::SparseDirectUMFPACK
-DEAL:Solvers::IterativeInverse
-DEAL:Solvers:cg::Starting value 16.88194302
-DEAL:Solvers:cg::Convergence step 1 value 0.000000000
 DEAL::SparseMatrix OK
 DEAL::PreconditionBlockIdentity
 DEAL:cg::Starting value 49.69909456

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.