]> https://gitweb.dealii.org/ - dealii.git/commitdiff
ShapeInfo: Implement evaluation of derivatives on faces 14383/head
authorMartin Kronbichler <martin.kronbichler@uni-a.de>
Tue, 1 Nov 2022 16:31:01 +0000 (17:31 +0100)
committerMartin Kronbichler <martin.kronbichler@uni-a.de>
Tue, 1 Nov 2022 16:31:01 +0000 (17:31 +0100)
include/deal.II/matrix_free/shape_info.h
include/deal.II/matrix_free/shape_info.templates.h

index 36ee133b22bacdf346259569694b514a4377b735..2715745dc504648681caf5a7db39c808669ce7f6 100644 (file)
@@ -247,8 +247,8 @@ namespace internal
 
       /**
        * Collects all data of 1D shape values evaluated at the point 0 and 1
-       * (the vertices) in one data structure. Sorting is first the values,
-       * then gradients, then second derivatives.
+       * (the vertices) in one data structure. The sorting of data is to
+       * start with the values, then gradients, then second derivatives.
        */
       std::array<AlignedVector<Number>, 2> shape_data_on_face;
 
@@ -258,9 +258,8 @@ namespace internal
        * point 0 and 1 (the vertices) in one data structure.
        *
        * This data structure can be used to interpolate from the cell to the
-       * face quadrature points.
-       *
-       * @note In contrast to shape_data_on_face, only the vales are evaluated.
+       * face quadrature points. The sorting of data is to start with the
+       * values, then gradients, then second derivatives.
        */
       std::array<AlignedVector<Number>, 2> quadrature_data_on_face;
 
index b9a7ba44199f6d6dd30bbbf3bee82ac34ff9ba7c..dd289ccd1af1100034708e9d721cd3d970dc0bee 100644 (file)
@@ -748,8 +748,13 @@ namespace internal
 
           for (unsigned int i = 0; i < quad.size(); ++i)
             {
-              quadrature_data_on_face[0][i] = poly_coll[i].value(0.0);
-              quadrature_data_on_face[1][i] = poly_coll[i].value(1.0);
+              std::array<double, 3> values;
+              poly_coll[i].value(0.0, 2, values.data());
+              for (unsigned int d = 0; d < 3; ++d)
+                quadrature_data_on_face[0][i + d * quad.size()] = values[d];
+              poly_coll[i].value(1.0, 2, values.data());
+              for (unsigned int d = 0; d < 3; ++d)
+                quadrature_data_on_face[1][i + d * quad.size()] = values[d];
             }
         }
 

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.