]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
Propagate exceptions from the member functions of InverseMatrixRichardson.
authorbangerth <bangerth@0785d39b-7218-0410-832d-ea1e28bc413d>
Sun, 16 Feb 2014 17:55:36 +0000 (17:55 +0000)
committerbangerth <bangerth@0785d39b-7218-0410-832d-ea1e28bc413d>
Sun, 16 Feb 2014 17:55:36 +0000 (17:55 +0000)
git-svn-id: https://svn.dealii.org/trunk@32494 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/doc/news/changes.h
deal.II/include/deal.II/lac/matrix_lib.h
deal.II/include/deal.II/lac/matrix_lib.templates.h

index 83cd398594ec0f3e20bb45fc7ea67dcf04a4a6cc..87185c59505da6a1376986275e1a4a3ee45f4cf6 100644 (file)
@@ -135,6 +135,14 @@ inconvenience this causes.
 <h3>Specific improvements</h3>
 
 <ol>
+  <li>Changed: The InverseMatrixRichardson used to eat all exceptions
+  that may have been produced by the underlying Richardson solver, leaving
+  no trace that the underlying solver may have failed when you call functions
+  such as InverseMatrixRichardson::vmult(). These exceptions are now propagated
+  out to the caller.
+  <br>
+  (Wolfgang Bangerth, 2014/02/16)
+
 
   <li>New: FE_TraceQ implements finite elements on faces, which
   correspond to the traces of H<sup>1</sup>-conforming elements.
index 03eb9b7a3d35a15f5d9d16f8938f6511385a3c63..3bd6318ea5a7b44192d45024278dece5d671065c 100644 (file)
@@ -415,15 +415,29 @@ private:
 
 
 /**
- * Inverse matrix computed approximately by using the SolverRichardson
- * iterative solver. In particular, the function
- * SolverRichardson::Tsolve() allows for the implementation of
- * transpose matrix vector products.
+ * Objects of this type represent the inverse of a matrix as
+ * computed approximately by using the SolverRichardson
+ * iterative solver. In other words, if you set up an object
+ * of the current type for a matrix $A$, then calling the
+ * vmult() function with arguments $v,w$ amounts to setting
+ * $w=A^{-1}v$ by solving the linear system $Aw=v$ using the
+ * Richardson solver with a preconditioner that can be chosen. Similarly,
+ * this class allows to also multiple with the transpose of the
+ * inverse (i.e., the inverse of the transpose) using the function
+ * SolverRichardson::Tsolve().
  *
  * The functions vmult() and Tvmult() approximate the inverse
  * iteratively starting with the vector <tt>dst</tt>. Functions
  * vmult_add() and Tvmult_add() start the iteration with a zero
- * vector.
+ * vector. All of the matrix-vector multiplication functions
+ * expect that the Richardson solver with the given preconditioner
+ * actually converge. If the Richardson solver does not converge
+ * within the specified number of iterations, the exception that will
+ * result in the solver will simply be propagated to the caller of
+ * the member function of the current class.
+ *
+ * @note A more powerful version of this class is provided by the
+ * IterativeInverse class.
  *
  * @note Instantiations for this template are provided for <tt>@<float@> and
  * @<double@></tt>; others can be generated in application programs (see the
index 81e0a09157467adb68151553d290b2ca0ab2465a..7d5b067d2da7405219615b85b6b0f7243a87d230 100644 (file)
@@ -152,12 +152,7 @@ InverseMatrixRichardson<VECTOR>::vmult(VECTOR &dst, const VECTOR &src) const
   Assert (matrix != 0, ExcNotInitialized());
   Assert (precondition != 0, ExcNotInitialized());
   dst = 0.;
-  try
-    {
-      solver.solve(*matrix, dst, src, *precondition);
-    }
-  catch (...)
-    {}
+  solver.solve(*matrix, dst, src, *precondition);
 }
 
 
@@ -170,12 +165,9 @@ InverseMatrixRichardson<VECTOR>::vmult_add(VECTOR &dst, const VECTOR &src) const
   Assert (precondition != 0, ExcNotInitialized());
   VECTOR *aux = mem.alloc();
   aux->reinit(dst);
-  try
-    {
-      solver.solve(*matrix, *aux, src, *precondition);
-    }
-  catch (...)
-    {}
+
+  solver.solve(*matrix, *aux, src, *precondition);
+
   dst += *aux;
   mem.free(aux);
 }
@@ -189,12 +181,7 @@ InverseMatrixRichardson<VECTOR>::Tvmult(VECTOR &dst, const VECTOR &src) const
   Assert (matrix != 0, ExcNotInitialized());
   Assert (precondition != 0, ExcNotInitialized());
   dst = 0.;
-  try
-    {
-      solver.Tsolve(*matrix, dst, src, *precondition);
-    }
-  catch (...)
-    {}
+  solver.Tsolve(*matrix, dst, src, *precondition);
 }
 
 
@@ -207,12 +194,9 @@ InverseMatrixRichardson<VECTOR>::Tvmult_add(VECTOR &dst, const VECTOR &src) cons
   Assert (precondition != 0, ExcNotInitialized());
   VECTOR *aux = mem.alloc();
   aux->reinit(dst);
-  try
-    {
-      solver.Tsolve(*matrix, *aux, src, *precondition);
-    }
-  catch (...)
-    {}
+
+  solver.Tsolve(*matrix, *aux, src, *precondition);
+
   dst += *aux;
   mem.free(aux);
 }

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.