15#ifndef dealii_portable_fe_evaluation_h
16#define dealii_portable_fe_evaluation_h
27#include <deal.II/matrix_free/portable_matrix_free.templates.h>
30#include <Kokkos_Core.hpp>
66 int n_q_points_1d = fe_degree + 1,
67 int n_components_ = 1,
68 typename Number =
double>
152 DEAL_II_DEPRECATED_EARLY_WITH_COMMENT(
153 "Use the version taking EvaluationFlags.")
156 evaluate(const
bool evaluate_val, const
bool evaluate_grad);
175 DEAL_II_DEPRECATED_EARLY_WITH_COMMENT(
179 integrate(const
bool integrate_val, const
bool integrate_grad);
235 template <typename Functor>
271 static_assert(n_components_ == 1,
"This function only supports FE with one \
274 Kokkos::parallel_for(Kokkos::TeamThreadRange(
shared_data->team_member,
277 shared_data->values(i) =
278 src[data->local_to_global(cell_id, i)];
284 data->constraint_weights,
300 static_assert(n_components_ == 1,
"This function only supports FE with one \
305 data->constraint_weights,
309 if (
data->use_coloring)
311 Kokkos::parallel_for(Kokkos::TeamThreadRange(
shared_data->team_member,
314 dst[data->local_to_global(cell_id, i)] +=
315 shared_data->values(i);
320 Kokkos::parallel_for(
323 Kokkos::atomic_add(&dst[data->local_to_global(cell_id, i)],
324 shared_data->values(i));
350 data->shape_gradients,
351 data->co_shape_gradients);
356 evaluator_tensor_product.evaluate_values_and_gradients(
362 evaluator_tensor_product.evaluate_gradients(
shared_data->values,
368 evaluator_tensor_product.evaluate_values(
shared_data->values);
382 const bool evaluate_val,
383 const bool evaluate_grad)
409 data->shape_gradients,
410 data->co_shape_gradients);
415 evaluator_tensor_product.integrate_values_and_gradients(
420 evaluator_tensor_product.integrate_values(
shared_data->values);
425 evaluator_tensor_product.template integrate_gradients<false>(
440 const bool integrate_val,
441 const bool integrate_grad)
523 Number>::gradient_type
527 static_assert(n_components_ == 1,
"This function only supports FE with one \
531 for (
unsigned int d_1 = 0; d_1 < dim; ++d_1)
534 for (
unsigned int d_2 = 0; d_2 < dim; ++d_2)
535 tmp +=
data->inv_jacobian(
cell_id, q_point, d_2, d_1) *
554 for (
unsigned int d_1 = 0; d_1 < dim; ++d_1)
557 for (
unsigned int d_2 = 0; d_2 < dim; ++d_2)
558 tmp +=
data->inv_jacobian(
cell_id, q_point, d_1, d_2) * grad_in[d_2];
571 template <
typename Functor>
576 Kokkos::parallel_for(Kokkos::TeamThreadRange(
shared_data->team_member,
578 [&](
const int &i) { func(this, i); });
590 constexpr unsigned int
static constexpr unsigned int n_q_points
SharedData< dim, Number > * shared_data
void integrate(const EvaluationFlags::EvaluationFlags integration_flag)
static constexpr unsigned int dimension
static constexpr unsigned int n_components
FEEvaluation(const data_type *data, SharedData< dim, Number > *shdata)
value_type get_value(int q_point) const
static constexpr unsigned int tensor_dofs_per_cell
value_type get_dof_value(int q_point) const
void evaluate(const EvaluationFlags::EvaluationFlags evaluate_flag)
void distribute_local_to_global(Number *dst) const
gradient_type get_gradient(int q_point) const
typename MatrixFree< dim, Number >::Data data_type
void read_dof_values(const Number *src)
void submit_dof_value(const value_type &val_in, int q_point)
void submit_value(const value_type &val_in, int q_point)
void submit_gradient(const gradient_type &grad_in, int q_point)
void apply_for_each_quad_point(const Functor &func)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
The namespace for the EvaluationFlags enum.
EvaluationFlags
The EvaluationFlags enum.
void resolve_hanging_nodes(const Kokkos::TeamPolicy< MemorySpace::Default::kokkos_space::execution_space >::member_type &team_member, Kokkos::View< Number *, MemorySpace::Default::kokkos_space > constraint_weights, const ::internal::MatrixFreeFunctions::ConstraintKinds constraint_mask, Kokkos::View< Number *, MemorySpace::Default::kokkos_space::execution_space::scratch_memory_space, Kokkos::MemoryTraits< Kokkos::Unmanaged > > values)
constexpr T pow(const T base, const int iexp)
#define DEAL_II_HOST_DEVICE