]> https://gitweb.dealii.org/ - dealii.git/commitdiff
gracefully fall-back to Lagrange interpolation if conformity is unknown
authorMatthias Maier <tamiko@43-1.org>
Fri, 8 Sep 2017 03:28:53 +0000 (22:28 -0500)
committerMatthias Maier <tamiko@43-1.org>
Sun, 17 Sep 2017 20:26:40 +0000 (15:26 -0500)
include/deal.II/numerics/vector_tools.templates.h

index 27adfdd4fd60fdc9748a02dac778c7bd63fbf455..7ae278f39693b36d5fe69d34fae2fa7222388f43 100644 (file)
@@ -242,16 +242,6 @@ namespace VectorTools
         // element for the correct transformation.
         switch (fe[fe_index].conforming_space)
           {
-          case FiniteElementData<dim>::H1:
-            DEAL_II_FALLTHROUGH;
-          case FiniteElementData<dim>::L2:
-            // See Monk, Finite Element Methods for Maxwell's Equations,
-            // p. 77ff, formula (3.74).
-            // For given mapping F_K: \hat K \to K, we have to transform
-            //  \hat p = p\circ F_K
-            //  i.e., do nothing.
-            break;
-
           case FiniteElementData<dim>::Hcurl:
             // See Monk, Finite Element Methods for Maxwell's Equations,
             // p. 77ff, formula (3.76) and Corollary 3.58.
@@ -295,10 +285,22 @@ namespace VectorTools
               }
             break;
 
+          case FiniteElementData<dim>::H1:
+            DEAL_II_FALLTHROUGH;
+          case FiniteElementData<dim>::L2:
+            // See Monk, Finite Element Methods for Maxwell's Equations,
+            // p. 77ff, formula (3.74).
+            // For given mapping F_K: \hat K \to K, we have to transform
+            //  \hat p = p\circ F_K
+            //  i.e., do nothing.
+            //
+            break;
+
           default:
-            Assert(false,
-                   ExcMessage(
-                     "The supplied finite element has an unknown conformity."));
+            // In case we deal with an unknown conformity, just assume we
+            // deal with a Lagrange element and do nothing.
+            break;
+
           } /*switch*/
 
         FETools::convert_generalized_support_point_values_to_dof_values(

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.