]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Adjust error messages in the error estimator class.
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Sun, 15 Mar 2015 06:28:52 +0000 (01:28 -0500)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Tue, 17 Mar 2015 23:27:43 +0000 (18:27 -0500)
include/deal.II/numerics/error_estimator.h
source/numerics/error_estimator.cc
source/numerics/error_estimator_1d.cc

index f92e5bbd475c93ea334c81329e585561e8661228..b3872389838041ddcfe8f08e9194e0650620e85b 100644 (file)
@@ -442,34 +442,51 @@ public:
   /**
    * Exception
    */
-  DeclException0 (ExcInvalidBoundaryIndicator);
+  DeclExceptionMsg (ExcInvalidComponentMask,
+                    "You provided a ComponentMask argument that is invalid. "
+                    "Component masks need to be either default constructed "
+                    "(in which case they indicate that every component is "
+                    "selected) or need to have a length equal to the number "
+                    "of vector components of the finite element in use "
+                    "by the DoFHandler object. In the latter case, at "
+                    "least one component needs to be selected.");
   /**
    * Exception
    */
-  DeclException0 (ExcInvalidComponentMask);
+  DeclExceptionMsg (ExcInvalidCoefficient,
+                    "If you do specify the argument for a (possibly "
+                    "spatially variable) coefficient function for this function, "
+                    "then it needs to refer to a coefficient that is either "
+                    "scalar (has one vector component) or has as many vector "
+                    "components as there are in the finite element used by "
+                    "the DoFHandler argument.");
   /**
    * Exception
    */
-  DeclException0 (ExcInvalidCoefficient);
-  /**
-   * Exception
-   */
-  DeclException0 (ExcInvalidBoundaryFunction);
+  DeclException3 (ExcInvalidBoundaryFunction,
+                  types::boundary_id,
+                  int,
+                  int,
+                  << "You provided a function map that for boundary indicator "
+                  << arg1 << " specifies a function with "
+                  << arg2 << " vector components. However, the finite "
+                  "element in use has "
+                  << arg2 << " components, and these two numbers need to match.");
   /**
    * Exception
    */
   DeclException2 (ExcIncompatibleNumberOfElements,
                   int, int,
-                  << "The number of elements " << arg1 << " and " << arg2
-                  << " of the vectors do not match!");
-  /**
-   * Exception
-   */
-  DeclException0 (ExcInvalidSolutionVector);
+                  << "The number of input vectors, " << arg1
+                  << " needs to be equal to the number of output vectors, "
+                  << arg2
+                  << ". This is not the case in your call of this function.");
   /**
    * Exception
    */
-  DeclException0 (ExcNoSolutions);
+  DeclExceptionMsg (ExcNoSolutions,
+                    "You need to specify at least one solution vector as "
+                    "input.");
 };
 
 
@@ -655,34 +672,51 @@ public:
   /**
    * Exception
    */
-  DeclException0 (ExcInvalidBoundaryIndicator);
+  DeclExceptionMsg (ExcInvalidComponentMask,
+                    "You provided a ComponentMask argument that is invalid. "
+                    "Component masks need to be either default constructed "
+                    "(in which case they indicate that every component is "
+                    "selected) or need to have a length equal to the number "
+                    "of vector components of the finite element in use "
+                    "by the DoFHandler object. In the latter case, at "
+                    "least one component needs to be selected.");
   /**
    * Exception
    */
-  DeclException0 (ExcInvalidComponentMask);
+  DeclExceptionMsg (ExcInvalidCoefficient,
+                    "If you do specify the argument for a (possibly "
+                    "spatially variable) coefficient function for this function, "
+                    "then it needs to refer to a coefficient that is either "
+                    "scalar (has one vector component) or has as many vector "
+                    "components as there are in the finite element used by "
+                    "the DoFHandler argument.");
   /**
    * Exception
    */
-  DeclException0 (ExcInvalidCoefficient);
-  /**
-   * Exception
-   */
-  DeclException0 (ExcInvalidBoundaryFunction);
+  DeclException3 (ExcInvalidBoundaryFunction,
+                  types::boundary_id,
+                  int,
+                  int,
+                  << "You provided a function map that for boundary indicator "
+                  << arg1 << " specifies a function with "
+                  << arg2 << " vector components. However, the finite "
+                  "element in use has "
+                  << arg3 << " components, and these two numbers need to match.");
   /**
    * Exception
    */
   DeclException2 (ExcIncompatibleNumberOfElements,
                   int, int,
-                  << "The number of elements " << arg1 << " and " << arg2
-                  << " of the vectors do not match!");
-  /**
-   * Exception
-   */
-  DeclException0 (ExcInvalidSolutionVector);
+                  << "The number of input vectors, " << arg1
+                  << " needs to be equal to the number of output vectors, "
+                  << arg2
+                  << ". This is not the case in your call of this function.");
   /**
    * Exception
    */
-  DeclException0 (ExcNoSolutions);
+  DeclExceptionMsg (ExcNoSolutions,
+                    "You need to specify at least one solution vector as "
+                    "input.");
 };
 
 
index 19ba50108a71022ea27546730fbb6d2cfabfb98a..2d6e2fe6fd32979485d09db5ed1514662f859dbf 100644 (file)
@@ -934,7 +934,9 @@ estimate (const Mapping<dim, spacedim>                  &mapping,
   for (typename FunctionMap<spacedim>::type::const_iterator i=neumann_bc.begin();
        i!=neumann_bc.end(); ++i)
     Assert (i->second->n_components == n_components,
-            ExcInvalidBoundaryFunction());
+            ExcInvalidBoundaryFunction(i->first,
+                                       i->second->n_components,
+                                       n_components));
 
   Assert (component_mask.represents_n_components(n_components),
           ExcInvalidComponentMask());
@@ -948,7 +950,8 @@ estimate (const Mapping<dim, spacedim>                  &mapping,
 
   for (unsigned int n=0; n<solutions.size(); ++n)
     Assert (solutions[n]->size() == dof_handler.n_dofs(),
-            ExcInvalidSolutionVector());
+            ExcDimensionMismatch(solutions[n]->size(),
+                                 dof_handler.n_dofs()));
 
   const unsigned int n_solution_vectors = solutions.size();
 
index e91bfe0ab18476efd3b6aac3b758d34b09c9aba6..6e95f5dad454b3030e88c4019960f5b5830661df 100644 (file)
@@ -250,11 +250,16 @@ estimate (const Mapping<1,spacedim>                    &mapping,
 
   // sanity checks
   Assert (neumann_bc.find(numbers::internal_face_boundary_id) == neumann_bc.end(),
-          ExcInvalidBoundaryIndicator());
+          ExcMessage("You are not allowed to list the special boundary "
+                     "indicator for internal boundaries in your boundary "
+                     "value map."));
 
   for (typename FunctionMap<spacedim>::type::const_iterator i=neumann_bc.begin();
        i!=neumann_bc.end(); ++i)
-    Assert (i->second->n_components == n_components, ExcInvalidBoundaryFunction());
+    Assert (i->second->n_components == n_components,
+            ExcInvalidBoundaryFunction(i->first,
+                                       i->second->n_components,
+                                       n_components));
 
   Assert (component_mask.represents_n_components(n_components),
           ExcInvalidComponentMask());
@@ -272,7 +277,8 @@ estimate (const Mapping<1,spacedim>                    &mapping,
           ExcIncompatibleNumberOfElements(solutions.size(), errors.size()));
   for (unsigned int n=0; n<solutions.size(); ++n)
     Assert (solutions[n]->size() == dof_handler.n_dofs(),
-            ExcInvalidSolutionVector());
+            ExcDimensionMismatch(solutions[n]->size(),
+                                 dof_handler.n_dofs()));
 
   Assert ((coefficient == 0) ||
           (coefficient->n_components == n_components) ||
@@ -282,7 +288,9 @@ estimate (const Mapping<1,spacedim>                    &mapping,
   for (typename FunctionMap<spacedim>::type::const_iterator i=neumann_bc.begin();
        i!=neumann_bc.end(); ++i)
     Assert (i->second->n_components == n_components,
-            ExcInvalidBoundaryFunction());
+            ExcInvalidBoundaryFunction(i->first,
+                                       i->second->n_components,
+                                       n_components));
 
   // reserve one slot for each cell and set it to zero
   for (unsigned int n=0; n<n_solution_vectors; ++n)

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.