Skip to content

Commit

Permalink
Merge pull request #14383 from kronbichler/fix_shape_data
Browse files Browse the repository at this point in the history
ShapeInfo: Implement evaluation of derivatives on faces
  • Loading branch information
peterrum committed Nov 1, 2022
2 parents 0326ce6 + c0f6e6a commit bc07b82
Show file tree
Hide file tree
Showing 2 changed files with 11 additions and 7 deletions.
9 changes: 4 additions & 5 deletions include/deal.II/matrix_free/shape_info.h
Original file line number Diff line number Diff line change
Expand Up @@ -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;

Expand All @@ -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;

Expand Down
9 changes: 7 additions & 2 deletions include/deal.II/matrix_free/shape_info.templates.h
Original file line number Diff line number Diff line change
Expand Up @@ -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];
}
}

Expand Down

0 comments on commit bc07b82

Please sign in to comment.