From: Peter Munch Date: Sun, 7 Jan 2024 13:06:29 +0000 (+0100) Subject: Rename some functions in CUDAWrappers::EvaluatorTensorProduct X-Git-Tag: relicensing~183^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=refs%2Fpull%2F16424%2Fhead;p=dealii.git Rename some functions in CUDAWrappers::EvaluatorTensorProduct --- diff --git a/include/deal.II/matrix_free/cuda_fe_evaluation.h b/include/deal.II/matrix_free/cuda_fe_evaluation.h index 282fcaf298..a5f26290fa 100644 --- a/include/deal.II/matrix_free/cuda_fe_evaluation.h +++ b/include/deal.II/matrix_free/cuda_fe_evaluation.h @@ -357,19 +357,19 @@ namespace CUDAWrappers if ((evaluate_flag & EvaluationFlags::values) && (evaluate_flag & EvaluationFlags::gradients)) { - evaluator_tensor_product.value_and_gradient_at_quad_pts( + evaluator_tensor_product.evaluate_values_and_gradients( shared_data->values, shared_data->gradients); shared_data->team_member.team_barrier(); } else if (evaluate_flag & EvaluationFlags::gradients) { - evaluator_tensor_product.gradient_at_quad_pts(shared_data->values, - shared_data->gradients); + evaluator_tensor_product.evaluate_gradients(shared_data->values, + shared_data->gradients); shared_data->team_member.team_barrier(); } else if (evaluate_flag & EvaluationFlags::values) { - evaluator_tensor_product.value_at_quad_pts(shared_data->values); + evaluator_tensor_product.evaluate_values(shared_data->values); shared_data->team_member.team_barrier(); } } @@ -416,17 +416,17 @@ namespace CUDAWrappers if ((integration_flag & EvaluationFlags::values) && (integration_flag & EvaluationFlags::gradients)) { - evaluator_tensor_product.integrate_value_and_gradient( + evaluator_tensor_product.integrate_values_and_gradients( shared_data->values, shared_data->gradients); } else if (integration_flag & EvaluationFlags::values) { - evaluator_tensor_product.integrate_value(shared_data->values); + evaluator_tensor_product.integrate_values(shared_data->values); shared_data->team_member.team_barrier(); } else if (integration_flag & EvaluationFlags::gradients) { - evaluator_tensor_product.template integrate_gradient( + evaluator_tensor_product.template integrate_gradients( shared_data->values, shared_data->gradients); shared_data->team_member.team_barrier(); } diff --git a/include/deal.II/matrix_free/cuda_tensor_product_kernels.h b/include/deal.II/matrix_free/cuda_tensor_product_kernels.h index 2a614281cc..8f9f678edf 100644 --- a/include/deal.II/matrix_free/cuda_tensor_product_kernels.h +++ b/include/deal.II/matrix_free/cuda_tensor_product_kernels.h @@ -357,6 +357,7 @@ namespace CUDAWrappers n_q_points_1d, Number> { + public: using TeamHandle = Kokkos::TeamPolicy< MemorySpace::Default::kokkos_space::execution_space>::member_type; @@ -370,89 +371,90 @@ namespace CUDAWrappers co_shape_gradients); /** - * Evaluate the values of a finite element function at the quadrature - * points. + * Evaluate the finite element function at the quadrature points. */ - template + template DEAL_II_HOST_DEVICE void - values(const ViewTypeIn in, ViewTypeOut out) const; + evaluate_values(ViewType u); /** - * Evaluate the gradient of a finite element function at the quadrature - * points for a given @p direction. + * Evaluate the gradients of the finite element function at the quadrature + * points. */ - template + template DEAL_II_HOST_DEVICE void - gradients(const ViewTypeIn in, ViewTypeOut out) const; + evaluate_gradients(const ViewTypeIn u, ViewTypeOut grad_u); /** - * Evaluate the gradient of a finite element function at the quadrature - * points for a given @p direction for collocation methods. + * Evaluate the values and the gradients of the finite element function at + * the quadrature points. */ - template + template DEAL_II_HOST_DEVICE void - co_gradients(const ViewTypeIn in, ViewTypeOut out) const; + evaluate_values_and_gradients(ViewType1 u, ViewType2 grad_u); /** - * Evaluate the finite element function at the quadrature points. + * Helper function for integrate(). Integrate the finite element function. */ template DEAL_II_HOST_DEVICE void - value_at_quad_pts(ViewType u); + integrate_values(ViewType u); /** - * Helper function for integrate(). Integrate the finite element function. + * Helper function for integrate(). Integrate the gradients of the finite + * element function. */ - template + template DEAL_II_HOST_DEVICE void - integrate_value(ViewType u); + integrate_gradients(ViewType1 u, ViewType2 grad_u); /** - * Evaluate the gradients of the finite element function at the quadrature - * points. + * Helper function for integrate(). Integrate the values and the gradients + * of the finite element function. */ - template + template DEAL_II_HOST_DEVICE void - gradient_at_quad_pts(const ViewTypeIn u, ViewTypeOut grad_u); + integrate_values_and_gradients(ViewType1 u, ViewType2 grad_u); /** - * Evaluate the values and the gradients of the finite element function at - * the quadrature points. + * Evaluate/integrate the values of a finite element function at the + * quadrature points for a given @p direction. */ - template + template DEAL_II_HOST_DEVICE void - value_and_gradient_at_quad_pts(ViewType1 u, ViewType2 grad_u); + values(const ViewTypeIn in, ViewTypeOut out) const; /** - * Helper function for integrate(). Integrate the gradients of the finite - * element function. + * Evaluate/integrate the gradient of a finite element function at the + * quadrature points for a given @p direction. */ - template + template DEAL_II_HOST_DEVICE void - integrate_gradient(ViewType1 u, ViewType2 grad_u); + gradients(const ViewTypeIn in, ViewTypeOut out) const; + public: /** - * Helper function for integrate(). Integrate the values and the gradients - * of the finite element function. + * Evaluate the gradient of a finite element function at the quadrature + * points for a given @p direction for collocation methods. */ - template + template DEAL_II_HOST_DEVICE void - integrate_value_and_gradient(ViewType1 u, ViewType2 grad_u); + co_gradients(const ViewTypeIn in, ViewTypeOut out) const; /** * TeamPolicy handle. @@ -571,7 +573,7 @@ namespace CUDAWrappers dim, fe_degree, n_q_points_1d, - Number>::value_at_quad_pts(ViewType u) + Number>::evaluate_values(ViewType u) { if constexpr (dim == 1) values<0, true, false, true>(u, u); @@ -602,7 +604,7 @@ namespace CUDAWrappers dim, fe_degree, n_q_points_1d, - Number>::integrate_value(ViewType u) + Number>::integrate_values(ViewType u) { if constexpr (dim == 1) values<0, false, false, true>(u, u); @@ -633,8 +635,8 @@ namespace CUDAWrappers dim, fe_degree, n_q_points_1d, - Number>::gradient_at_quad_pts(const ViewTypeIn u, - ViewTypeOut grad_u) + Number>::evaluate_gradients(const ViewTypeIn u, + ViewTypeOut grad_u) { if constexpr (dim == 1) { @@ -698,9 +700,9 @@ namespace CUDAWrappers dim, fe_degree, n_q_points_1d, - Number>::value_and_gradient_at_quad_pts(ViewType1 u, - ViewType2 - grad_u) + Number>::evaluate_values_and_gradients(ViewType1 u, + ViewType2 + grad_u) { if constexpr (dim == 1) { @@ -751,8 +753,8 @@ namespace CUDAWrappers dim, fe_degree, n_q_points_1d, - Number>::integrate_gradient(ViewType1 u, - ViewType2 grad_u) + Number>::integrate_gradients(ViewType1 u, + ViewType2 grad_u) { if constexpr (dim == 1) { @@ -829,9 +831,9 @@ namespace CUDAWrappers dim, fe_degree, n_q_points_1d, - Number>::integrate_value_and_gradient(ViewType1 u, - ViewType2 - grad_u) + Number>::integrate_values_and_gradients(ViewType1 u, + ViewType2 + grad_u) { if constexpr (dim == 1) {