]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Get rid of the confusing ExcInvalidConstructorCall exception. 3015/head
authorWolfgang Bangerth <bangerth@colostate.edu>
Thu, 25 Aug 2016 16:19:03 +0000 (10:19 -0600)
committerWolfgang Bangerth <bangerth@colostate.edu>
Thu, 25 Aug 2016 16:19:03 +0000 (10:19 -0600)
This exception did not print any useful error message, and it was used in a significant
number of places that just happened to check obscure conditions, without really
explaining in the error message what concretely they were doing. Replace all of these
with ExcMessage and provide explanations what exactly it is that's going wrong.

include/deal.II/base/exceptions.h
include/deal.II/base/point.h
include/deal.II/lac/chunk_sparse_matrix.templates.h
include/deal.II/lac/sparse_matrix.templates.h
include/deal.II/lac/sparse_matrix_ez.h
include/deal.II/lac/sparse_matrix_ez.templates.h
source/base/quadrature_lib.cc
source/lac/dynamic_sparsity_pattern.cc

index dd5492f816d8e90c03207bb2048f35b274cb0bb7..507fecd3d653e67c8ec67864267c0ae6fd775424 100644 (file)
@@ -819,13 +819,6 @@ namespace StandardExceptions
                     "been called; if the answer is 'yes', then you need to "
                     "implement the missing override in your class.");
 
-  /**
-   * Used for constructors that are disabled. Examples are copy constructors
-   * and assignment operators of large objects, which are only allowed for
-   * empty objects.
-   */
-  DeclException0 (ExcInvalidConstructorCall);
-
   /**
    * This exception is used if some object is found uninitialized.
    */
index 67a0abba54c4d618526f16ca34979d97de53f66f..8eed55b57e2593cd38037bd8e89874607c611c42 100644 (file)
@@ -263,12 +263,23 @@ template <int dim, typename Number>
 inline
 Point<dim,Number>::Point (const Number x)
 {
+  Assert (dim==1,
+          ExcMessage ("You can only initialize Point<1> objects using the constructor "
+                      "that takes only one argument. Point<dim> objects with dim!=1 "
+                      "require initialization with the constructor that takes 'dim' "
+                      "arguments."));
+
+  // we can only get here if we pass the assertion. use the switch anyway so
+  // as to avoid compiler warnings about uninitialized elements or writing
+  // beyond the end of the 'values' array
   switch (dim)
     {
     case 1:
       this->values[0] = x;
+      break;
+
     default:
-      Assert (dim==1, StandardExceptions::ExcInvalidConstructorCall());
+      ;
     }
 }
 
@@ -276,15 +287,26 @@ Point<dim,Number>::Point (const Number x)
 
 template <int dim, typename Number>
 inline
-Point<dim,Number>::Point (const Number x, const Number y)
+Point<dim,Number>::Point (const Number x,
+                          const Number y)
 {
+  Assert (dim==2,
+          ExcMessage ("You can only initialize Point<2> objects using the constructor "
+                      "that takes two arguments. Point<dim> objects with dim!=2 "
+                      "require initialization with the constructor that takes 'dim' "
+                      "arguments."));
+  // we can only get here if we pass the assertion. use the switch anyway so
+  // as to avoid compiler warnings about uninitialized elements or writing
+  // beyond the end of the 'values' array
   switch (dim)
     {
     case 2:
       this->values[0] = x;
       this->values[1] = y;
+      break;
+
     default:
-      Assert (dim==2, StandardExceptions::ExcInvalidConstructorCall());
+      ;
     }
 }
 
@@ -292,16 +314,29 @@ Point<dim,Number>::Point (const Number x, const Number y)
 
 template <int dim, typename Number>
 inline
-Point<dim,Number>::Point (const Number x, const Number y, const Number z)
+Point<dim,Number>::Point (const Number x,
+                          const Number y,
+                          const Number z)
 {
+  Assert (dim==3,
+          ExcMessage ("You can only initialize Point<3> objects using the constructor "
+                      "that takes three arguments. Point<dim> objects with dim!=3 "
+                      "require initialization with the constructor that takes 'dim' "
+                      "arguments."));
+
+  // we can only get here if we pass the assertion. use the switch anyway so
+  // as to avoid compiler warnings about uninitialized elements or writing
+  // beyond the end of the 'values' array
   switch (dim)
     {
     case 3:
       this->values[0] = x;
       this->values[1] = y;
       this->values[2] = z;
+      break;
+
     default:
-      Assert (dim==3, StandardExceptions::ExcInvalidConstructorCall());
+      ;
     }
 }
 
index 76ed57c6282dfa054e86af52911ae13e9c347cbe..9c09c22804b2fccd1367e4e1729d59b38ec706bb 100644 (file)
@@ -300,9 +300,11 @@ ChunkSparseMatrix<number>::ChunkSparseMatrix (const ChunkSparseMatrix &m)
   val(0),
   max_len(0)
 {
-  Assert (m.cols==0, ExcInvalidConstructorCall());
-  Assert (m.val==0, ExcInvalidConstructorCall());
-  Assert (m.max_len==0, ExcInvalidConstructorCall());
+  Assert (m.cols==0 && m.val==0 && m.max_len==0,
+          ExcMessage("This constructor can only be called if the provided argument "
+                     "is an empty matrix. This constructor can not be used to "
+                     "copy-construct a non-empty matrix. Use the "
+                     "ChunkSparseMatrix::copy_from() function for that purpose."));
 }
 
 
@@ -312,9 +314,11 @@ ChunkSparseMatrix<number> &
 ChunkSparseMatrix<number>::operator = (const ChunkSparseMatrix<number> &m)
 {
   (void)m;
-  Assert (m.cols==0, ExcInvalidConstructorCall());
-  Assert (m.val==0, ExcInvalidConstructorCall());
-  Assert (m.max_len==0, ExcInvalidConstructorCall());
+  Assert (m.cols==0 && m.val==0 && m.max_len==0,
+          ExcMessage("This operator can only be called if the provided right "
+                     "hand side is an empty matrix. This operator can not be "
+                     "used to copy a non-empty matrix. Use the "
+                     "ChunkSparseMatrix::copy_from() function for that purpose."));
 
   return *this;
 }
index abc199329fbd272ad1b2979222c889523b055fe5..f07ec8b0af4a7b1b7e48445f9d358e35f96a40d6 100644 (file)
@@ -62,9 +62,11 @@ SparseMatrix<number>::SparseMatrix (const SparseMatrix &m)
   val(0),
   max_len(0)
 {
-  Assert (m.cols==0, ExcInvalidConstructorCall());
-  Assert (m.val==0, ExcInvalidConstructorCall());
-  Assert (m.max_len==0, ExcInvalidConstructorCall());
+  Assert (m.cols==0 && m.val==0 && m.max_len==0,
+          ExcMessage("This constructor can only be called if the provided argument "
+                     "is an empty matrix. This constructor can not be used to "
+                     "copy-construct a non-empty matrix. Use the "
+                     "SparseMatrix::copy_from() function for that purpose."));
 }
 
 
@@ -91,9 +93,11 @@ SparseMatrix<number> &
 SparseMatrix<number>::operator = (const SparseMatrix<number> &m)
 {
   (void)m;
-  Assert (m.cols==0, ExcInvalidConstructorCall());
-  Assert (m.val==0, ExcInvalidConstructorCall());
-  Assert (m.max_len==0, ExcInvalidConstructorCall());
+  Assert (m.cols==0 && m.val==0 && m.max_len==0,
+          ExcMessage("This operator can only be called if the provided right "
+                     "hand side is an empty matrix. This operator can not be "
+                     "used to copy a non-empty matrix. Use the "
+                     "SparseMatrix::copy_from() function for that purpose."));
 
   return *this;
 }
index c78c25ef50b66a2ca54a8edbd9db365855d92310..8d45c0392e556989a5397ae69108cd6fdad5cee7 100644 (file)
@@ -541,7 +541,8 @@ public:
    */
   template <typename MatrixType>
   SparseMatrixEZ<number> &
-  copy_from (const MatrixType &source, const bool elide_zero_values = true);
+  copy_from (const MatrixType &source,
+             const bool elide_zero_values = true);
 
   /**
    * Add @p matrix scaled by @p factor to this matrix.
index db928932d305e44404ff37521f10244373ca126f..4b522b166edc7fd539b5cf1e15dea6fb1fdf36e4 100644 (file)
@@ -43,8 +43,11 @@ SparseMatrixEZ<number>::SparseMatrixEZ(const SparseMatrixEZ<number> &m)
   Subscriptor (m),
   n_columns (0)
 {
-  Assert(m.n() == 0, ExcNotImplemented());
-  Assert(m.m() == 0, ExcNotImplemented());
+  Assert (m.m()==0 && m.n()==0,
+          ExcMessage("This constructor can only be called if the provided argument "
+                     "is an empty matrix. This constructor can not be used to "
+                     "copy-construct a non-empty matrix. Use the "
+                     "SparseMatrixEZ::copy_from() function for that purpose."));
 }
 
 
@@ -68,7 +71,11 @@ SparseMatrixEZ<number> &
 SparseMatrixEZ<number>::operator= (const SparseMatrixEZ<number> &m)
 {
   (void)m;
-  Assert (m.empty(), ExcInvalidConstructorCall());
+  Assert (m.empty(),
+          ExcMessage("This operator can only be called if the provided right "
+                     "hand side is an empty matrix. This operator can not be "
+                     "used to copy a non-empty matrix. Use the "
+                     "SparseMatrixEZ::copy_from() function for that purpose."));
   return *this;
 }
 
index 107f851fbb6fc82e14163d3621dacbb897668170..1007fb52deb59a9988d6692c023d1cb43c41c9e1 100644 (file)
@@ -1208,12 +1208,24 @@ QGaussRadauChebyshev<1>::get_quadrature_points(const unsigned int n,
     // would be -cos(2i Pi/(2N+1))
     // put + Pi so we start from the smallest point
     // then map from [-1,1] to [0,1]
-    if (ep == QGaussRadauChebyshev::left)
-      points[i] = 1./2.*(1.-std::cos(numbers::PI*(1+2*double(i)/(2*double(n-1)+1.))));
-    else
+    switch (ep)
+      {
+      case QGaussRadauChebyshev::left:
+      {
+        points[i] = 1./2.*(1.-std::cos(numbers::PI*(1+2*double(i)/(2*double(n-1)+1.))));
+        break;
+      }
+
+      case QGaussRadauChebyshev::right:
       {
-        Assert(ep==QGaussRadauChebyshev::right,ExcInvalidConstructorCall());
         points[i] = 1./2.*(1.-std::cos(numbers::PI*(2*double(n-1-i)/(2*double(n-1)+1.))));
+        break;
+      }
+
+      default:
+        Assert (false, ExcMessage ("This constructor can only be called with either "
+                                   "QGaussRadauChebyshev::left or QGaussRadauChebyshev::right as "
+                                   "second argument."));
       }
 
   return points;
index d097f52f07c5fae05fbffe79ff9b95001a57f2a0..e1b84a58fa0dca893ae0b364e1ca9057a2fd67d0 100644 (file)
@@ -234,8 +234,10 @@ DynamicSparsityPattern (const DynamicSparsityPattern &s)
   rowset (0)
 {
   (void)s;
-  Assert (s.rows == 0, ExcInvalidConstructorCall());
-  Assert (s.cols == 0, ExcInvalidConstructorCall());
+  Assert (s.rows==0 && s.cols==0,
+          ExcMessage("This constructor can only be called if the provided argument "
+                     "is the sparsity pattern for an empty matrix. This constructor can "
+                     "not be used to copy-construct a non-empty sparsity pattern."));
 }
 
 
@@ -281,11 +283,14 @@ DynamicSparsityPattern &
 DynamicSparsityPattern::operator = (const DynamicSparsityPattern &s)
 {
   (void)s;
-  Assert (s.rows == 0, ExcInvalidConstructorCall());
-  Assert (s.cols == 0, ExcInvalidConstructorCall());
+  Assert (s.rows==0 && s.cols==0,
+          ExcMessage("This operator can only be called if the provided argument "
+                     "is the sparsity pattern for an empty matrix. This operator can "
+                     "not be used to copy a non-empty sparsity pattern."));
 
-  Assert (rows == 0, ExcInvalidConstructorCall());
-  Assert (cols == 0, ExcInvalidConstructorCall());
+  Assert (rows==0 && cols==0,
+          ExcMessage("This operator can only be called if the current object is"
+                     "empty."));
 
   return *this;
 }
@@ -302,7 +307,13 @@ DynamicSparsityPattern::reinit (const size_type m,
   cols = n;
   rowset=rowset_;
 
-  Assert(rowset.size()==0 || rowset.size() == m, ExcInvalidConstructorCall());
+  Assert(rowset.size()==0 || rowset.size() == m,
+         ExcMessage("The IndexSet argument to this function needs to either "
+                    "be empty (indicating the complete set of rows), or have size "
+                    "equal to the desired number of rows as specified by the "
+                    "first argument to this function. (Of course, the number "
+                    "of indices in this IndexSet may be less than the number "
+                    "of rows, but the *size* of the IndexSet must be equal.)"));
 
   std::vector<Line> new_lines (rowset.size()==0 ? rows : rowset.n_elements());
   lines.swap (new_lines);

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.