Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Work on interpolation from cell to face quadrature points #10870

Merged
merged 1 commit into from
Sep 4, 2020

Conversation

peterrum
Copy link
Member

@peterrum peterrum commented Sep 1, 2020

Extracted from #10830. In the ECL context, values are available at cell quadrature points. The face quadrature points can be simply computed via a 1D interpolation as proposed here.

Comment on lines 1745 to 1753
template <bool do_evaluate, bool add_into_output>
static void
interpolate(const internal::MatrixFreeFunctions::ShapeInfo<Number> &data,
const Number * input,
Number * output,
const bool do_gradients,
const unsigned int face_no)
{
AssertDimension(n_points_1d, data.data.front().quadrature.size());
Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Since this function is almost the same as the other FEFaceNormalEvaluationImpl, I wonder whether we can unify the calls and let both function call an internal function interpolate_generic. What I'm thinking about is to let interpolate and interpolate_quadrature extract the data array (shape_data_on_face[face_no % 2] or quadrature_data_on_face[face_no % 2]) as well as the in_stride and out_stride and let the interpolate_generic (or whatever name we choose) do the job.

Copy link
Member Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Good idea. What I was thinking about as a follow-up PR is to make the function dimension independent since I want to use it in 6D. What do you think?

Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Sounds good.

Copy link
Member

@kronbichler kronbichler left a comment

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

A few more comments.

Comment on lines 233 to 234
* Collects all data of 1D nodal shape values (that defined by the
* support of the quadrature rule) evaluated at the point 0 and 1
Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Suggested change
* Collects all data of 1D nodal shape values (that defined by the
* support of the quadrature rule) evaluated at the point 0 and 1
* Collects all data of 1D nodal shape values (defined by the Lagrange
* polynomials in the points of the quadrature rule) evaluated at the point 0 and 1

Comment on lines 303 to 318
quadrature_data_on_face[0].resize(quad.size() * 3);
quadrature_data_on_face[1].resize(quad.size() * 3);

dealii::FE_DGQArbitraryNodes<1> fe_quad(quad);

for (unsigned int i = 0; i < quad.size(); ++i)
{
Point<1> q_point;
q_point[0] = 0;
quadrature_data_on_face[0][i] = fe_quad.shape_value(i, q_point);
q_point[0] = 1;
quadrature_data_on_face[1][i] = fe_quad.shape_value(i, q_point);
}
Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I think you need to encapsulate this part in some if for the case the number of quadrature points is very high, say a few hundreds, look at what we do elsewhere in this function for the auxiliary FE. In that case, the FE_DGQArbitraryNodes polynomials would underflow and result in a division by zero. The collocation is not really helpful in that case anyway.

@peterrum peterrum force-pushed the quadrature_data_on_face branch 2 times, most recently from c55b518 to e94bc78 Compare September 3, 2020 17:19
@peterrum
Copy link
Member Author

peterrum commented Sep 3, 2020

Done :)

Copy link
Member

@kronbichler kronbichler left a comment

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Excellent.

Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment
Projects
None yet
Development

Successfully merging this pull request may close these issues.

None yet

2 participants