]> https://gitweb.dealii.org/ - dealii.git/commitdiff
add a test
authorMatthias Maier <tamiko@43-1.org>
Tue, 26 Sep 2017 16:40:25 +0000 (11:40 -0500)
committerMatthias Maier <tamiko@43-1.org>
Wed, 27 Sep 2017 17:13:36 +0000 (12:13 -0500)
include/deal.II/numerics/vector_tools.templates.h
tests/numerics/generalized_interpolation_04.output
tests/numerics/generalized_interpolation_05.cc [new file with mode: 0644]
tests/numerics/generalized_interpolation_05.output [new file with mode: 0644]

index 04966b4bffe5c914ea341564e55b0c55c70ae92b..c1c98a6208e661fb01e5a93a91b77b57d22349ea 100644 (file)
@@ -238,15 +238,21 @@ namespace VectorTools
 
         // A small helper function to transform a component range starting
         // at offset from the real to the unit cell according to the
-        // supplied conformity.
+        // supplied conformity. The function_values vector is transformed
+        // in place.
         //
         // FIXME: This should be refactored into the mapping (i.e.
-        // implement the inverse function of Mapping<dim,
-        // spacedim>::transform). Further, the finite element should make
-        // the information about the correct mapping directly accessible -
-        // fe.conforming_space is not the right call (thing about BDM).
-        const auto transform = [&](const typename FiniteElementData<dim>::Conformity conformity,
-                                   unsigned int offset)
+        // implement the inverse function of Mapping<dim, spacedim>::transform).
+        // Further, the finite element should make the information about
+        // the correct mapping directly accessible (i.e. which MappingType
+        // should be used). Using fe.conforming_space might be a bit of a
+        // problem because we only support doing nothing, Hcurl, and Hdiv
+        // conforming mappings.
+        //
+
+        const auto transform = [&function_values, &fe_values_jacobians, &cell](
+            const typename FiniteElementData<dim>::Conformity conformity,
+            const unsigned int offset)
         {
           switch (conformity)
             {
index 59027680b3a36e465e62e96ebc487a321d02d257..ce616a49ad48b7dbb7fecc20596d5b4f4a6ac8d6 100644 (file)
@@ -1,3 +1,3 @@
 
-DEAL::dim 2 FESystem<2>[FE_RaviartThomas<2>(1)^2-FE_Q<2>(1)-FESystem<2>[FE_Nedelec<2>(2)^2]-FESystem<2>[FE_Q<2>(1)^2]]
-DEAL::Check projection property: 20.8500
+DEAL::dim 2 FESystem<2>[FE_RaviartThomas<2>(1)-FE_Q<2>(1)-FE_Nedelec<2>(2)^2-FESystem<2>[FE_Q<2>(1)^2]]
+DEAL::Check projection property: 0.00000
diff --git a/tests/numerics/generalized_interpolation_05.cc b/tests/numerics/generalized_interpolation_05.cc
new file mode 100644 (file)
index 0000000..7d82fe3
--- /dev/null
@@ -0,0 +1,45 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2017 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+// Check projection property of VectorTools::interpolate for a complex,
+// staggered system of Hdiv / Hcurl / L2 conforming spaces.
+
+#include <deal.II/fe/fe_nedelec.h>
+#include <deal.II/fe/fe_q.h>
+#include <deal.II/fe/fe_raviart_thomas.h>
+#include <deal.II/fe/fe_system.h>
+
+#include "../tests.h"
+#include "generalized_interpolation.h"
+
+int main ()
+{
+  initlog();
+
+  static const int dim = 2;
+
+  FESystem<dim> fe(FE_RaviartThomas<dim>(1),
+                   2,
+                   FE_Q<dim>(1),
+                   1,
+                   FESystem<dim>(FE_Nedelec<dim>(2), 2),
+                   1,
+                   FESystem<dim>(FE_Q<dim>(1), dim),
+                   1);
+
+  const unsigned int n_comp = fe.n_components();
+
+  test<2>(fe, F<2>(n_comp, 1), 1, false);
+}
diff --git a/tests/numerics/generalized_interpolation_05.output b/tests/numerics/generalized_interpolation_05.output
new file mode 100644 (file)
index 0000000..e6dc071
--- /dev/null
@@ -0,0 +1,3 @@
+
+DEAL::dim 2 FESystem<2>[FE_RaviartThomas<2>(1)^2-FE_Q<2>(1)-FESystem<2>[FE_Nedelec<2>(2)^2]-FESystem<2>[FE_Q<2>(1)^2]]
+DEAL::Check projection property: 0.00000

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.