]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Fix FE_ABF::convert_generalized_support_point_values_to_nodal_values 4114/head
authorDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Sun, 26 Mar 2017 15:44:07 +0000 (17:44 +0200)
committerDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Sun, 26 Mar 2017 15:51:41 +0000 (17:51 +0200)
source/fe/fe_abf.cc
tests/fe/abf_03.cc [new file with mode: 0644]
tests/fe/abf_03.output [new file with mode: 0644]

index 86d0b0b8ac1fb8d04e6352ce5ed9ed50fb0fd6e0..4ad9967d18e0acc6563efcd571782a064fc32580 100644 (file)
@@ -558,7 +558,7 @@ convert_generalized_support_point_values_to_nodal_values (const std::vector<Vect
           unsigned int k = QProjector<dim>::DataSetDescriptor::face (face, false, false, false, n_face_points);
           for (unsigned int i=0; i<boundary_weights_abf.size(1); ++i)
             nodal_values[start_abf_dofs+i] += n_orient * boundary_weights_abf(k + fp, i)
-                                              * support_point_values[k + fp][GeometryInfo<dim>::unit_normal_direction[face]];
+                                              * support_point_values[face*n_face_points+fp][GeometryInfo<dim>::unit_normal_direction[face]];
         }
     }
 
@@ -632,7 +632,7 @@ FE_ABF<dim>::interpolate(
           unsigned int k = QProjector<dim>::DataSetDescriptor::face (face, false, false, false, n_face_points);
           for (unsigned int i=0; i<boundary_weights_abf.size(1); ++i)
             local_dofs[start_abf_dofs+i] += n_orient * boundary_weights_abf(k + fp, i)
-                                            * values[k + fp](GeometryInfo<dim>::unit_normal_direction[face]+offset);
+                                            * values[face*n_face_points+fp](GeometryInfo<dim>::unit_normal_direction[face]+offset);
         }
     }
 
@@ -694,7 +694,7 @@ FE_ABF<dim>::interpolate(
           unsigned int k = QProjector<dim>::DataSetDescriptor::face (face, false, false, false, n_face_points);
           for (unsigned int i=0; i<boundary_weights_abf.size(1); ++i)
             local_dofs[start_abf_dofs+i] += n_orient * boundary_weights_abf(k + fp, i)
-                                            * values[GeometryInfo<dim>::unit_normal_direction[face]][k + fp];
+                                            * values[GeometryInfo<dim>::unit_normal_direction[face]][face*n_face_points+fp];
         }
     }
 
diff --git a/tests/fe/abf_03.cc b/tests/fe/abf_03.cc
new file mode 100644 (file)
index 0000000..a49feb0
--- /dev/null
@@ -0,0 +1,80 @@
+// ---------------------------------------------------------------------
+//
+// 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 that convert_generalized_support_point_values_to_nodal_values
+// gives the correct when asking for we the nodal values for a
+// function in the FE_ABF<dim(degree) ansatz space.
+
+#include "../tests.h"
+#include <deal.II/fe/fe_abf.h>
+
+
+
+template<unsigned int dim>
+void test (const unsigned int degree)
+{
+  FE_ABF<dim> fe(degree);
+
+  std::vector<double> dof_values (fe.dofs_per_cell);
+  for (unsigned int i=0; i<dof_values.size(); ++i)
+    dof_values[i] = 1. + 2.*(double)Testing::rand()/double(RAND_MAX);
+
+  const std::vector<Point<dim> > &generalized_support_points = fe.get_generalized_support_points();
+  std::vector<Vector<double> > real_values (generalized_support_points.size(), Vector<double>(dim));
+
+  for (unsigned int i=0; i<generalized_support_points.size(); ++i)
+    for (unsigned int j=0; j<fe.dofs_per_cell; ++j)
+      for (unsigned int c=0; c<dim; ++c)
+        real_values[i][c] += dof_values[j] * fe.shape_value_component(j, generalized_support_points[i], c);
+
+  std::vector<double> compare_values (fe.dofs_per_cell);
+  fe.convert_generalized_support_point_values_to_nodal_values(real_values, compare_values);
+
+  for (unsigned int i=0; i<dof_values.size(); ++i)
+    if (std::abs(dof_values[i]-compare_values[i]) > std::abs(dof_values[i])*1.e-6)
+      {
+        deallog << i << ": " << dof_values[i] << " vs. " << compare_values[i] << std::endl;
+        AssertThrow(false, ExcInternalError());
+      }
+  deallog << "OK" << std::endl;
+}
+
+int main()
+{
+  initlog();
+
+  deallog.push("2d");
+  deallog.push("0");
+  test<2>(0);
+  deallog.pop();
+  deallog.push("1");
+  test<2>(1);
+  deallog.pop();
+  deallog.push("2");
+  test<2>(2);
+  deallog.pop();
+  deallog.pop();
+
+  deallog.push("3d");
+  deallog.push("0");
+  test<3>(0);
+  deallog.pop();
+  //deallog.push("1");
+  //test<3>(1);
+  //deallog.pop();
+  deallog.pop();
+}
diff --git a/tests/fe/abf_03.output b/tests/fe/abf_03.output
new file mode 100644 (file)
index 0000000..b698ee2
--- /dev/null
@@ -0,0 +1,5 @@
+
+DEAL:2d:0::OK
+DEAL:2d:1::OK
+DEAL:2d:2::OK
+DEAL:3d:0::OK

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.