15#ifndef dealii_portable_fe_evaluation_h
16#define dealii_portable_fe_evaluation_h
28#include <deal.II/matrix_free/portable_matrix_free.templates.h>
31#include <Kokkos_Core.hpp>
65 int n_q_points_1d = fe_degree + 1,
66 int n_components_ = 1,
67 typename Number =
double>
74 using value_type = std::conditional_t<(n_components_ == 1),
84 std::conditional_t<n_components_ == dim,
133 const unsigned int dof_index = 0);
244 template <
typename Functor>
317 Kokkos::parallel_for(Kokkos::TeamThreadRange(
data->team_member,
320 for (unsigned int c = 0; c < n_components_; ++c)
321 shared_data->values(i, c) =
322 src[precomputed_data->local_to_global(
323 i + tensor_dofs_per_component * c, cell_id)];
325 data->team_member.team_barrier();
327 for (
unsigned int c = 0; c < n_components_; ++c)
336 Kokkos::subview(
shared_data->values, Kokkos::ALL, c));
351 for (
unsigned int c = 0; c < n_components_; ++c)
360 Kokkos::subview(
shared_data->values, Kokkos::ALL, c));
365 Kokkos::parallel_for(
368 for (unsigned int c = 0; c < n_components_; ++c)
369 dst[precomputed_data->local_to_global(
370 i + tensor_dofs_per_component * c, cell_id)] +=
371 shared_data->values(i, c);
376 Kokkos::parallel_for(
379 for (unsigned int c = 0; c < n_components_; ++c)
380 Kokkos::atomic_add(&dst[precomputed_data->local_to_global(
381 i + (tensor_dofs_per_component)*c, cell_id)],
382 shared_data->values(i, c));
400 if (fe_degree >= 0 && fe_degree + 1 == n_q_points_1d &&
402 ElementType::tensor_symmetric_collocation)
409 else if (fe_degree >= 0 &&
420 ElementType::tensor_symmetric_no_collocation)
427 Kokkos::abort(
"The element type is not yet supported by the portable "
428 "matrix-free module.");
445 if (fe_degree >= 0 && fe_degree + 1 == n_q_points_1d &&
447 ElementType::tensor_symmetric_collocation)
454 else if (fe_degree >= 0 &&
465 ElementType::tensor_symmetric_no_collocation)
472 Kokkos::abort(
"The element type is not yet supported by the portable "
473 "matrix-free module.");
493 if constexpr (n_components_ == 1)
523 if constexpr (n_components_ == 1)
548 if constexpr (n_components_ == 1)
573 if constexpr (n_components_ == 1)
595 Number>::gradient_type
602 if constexpr (n_components_ == 1)
604 for (
unsigned int d_1 = 0; d_1 < dim; ++d_1)
607 for (
unsigned int d_2 = 0; d_2 < dim; ++d_2)
617 for (
unsigned int d_1 = 0; d_1 < dim; ++d_1)
620 for (
unsigned int d_2 = 0; d_2 < dim; ++d_2)
643 if constexpr (n_components_ == 1)
645 for (
unsigned int d_1 = 0; d_1 < dim; ++d_1)
648 for (
unsigned int d_2 = 0; d_2 < dim; ++d_2)
659 for (
unsigned int d_1 = 0; d_1 < dim; ++d_1)
662 for (
unsigned int d_2 = 0; d_2 < dim; ++d_2)
679 template <
typename Functor>
684 Kokkos::parallel_for(Kokkos::TeamThreadRange(
data->team_member,
n_q_points),
685 [&](
const int &i) { func(this, i); });
686 data->team_member.team_barrier();
697 constexpr unsigned int
friend class FEEvaluation
std::conditional_t< n_components_==1, Tensor< 1, dim, Number >, std::conditional_t< n_components_==dim, Tensor< 2, dim, Number >, Tensor< 1, n_components_, Tensor< 1, dim, Number > > > > gradient_type
int get_current_cell_index()
static constexpr unsigned int n_q_points
SharedData< dim, Number > * shared_data
void integrate(const EvaluationFlags::EvaluationFlags integration_flag)
static constexpr unsigned int dimension
const data_type * get_matrix_free_data()
const MatrixFree< dim, Number >::Data * data
static constexpr unsigned int n_components
void evaluate(const EvaluationFlags::EvaluationFlags evaluate_flag)
void read_dof_values(const DeviceVector< Number > &src)
value_type get_value(int q_point) const
void distribute_local_to_global(DeviceVector< Number > &dst) const
static constexpr unsigned int tensor_dofs_per_cell
value_type get_dof_value(int q_point) const
const MatrixFree< dim, Number >::PrecomputedData * precomputed_data
std::conditional_t<(n_components_==1), Number, Tensor< 1, n_components_, Number > > value_type
gradient_type get_gradient(int q_point) const
static constexpr unsigned int tensor_dofs_per_component
typename MatrixFree< dim, Number >::Data data_type
FEEvaluation(const data_type *data, const unsigned int dof_index=0)
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_HOST_DEVICE
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
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, ViewType values)
constexpr bool use_collocation_evaluation(const unsigned int fe_degree, const unsigned int n_q_points_1d)
Kokkos::View< Number *, MemorySpace::Default::kokkos_space > DeviceVector
constexpr T pow(const T base, const int iexp)
static void integrate(const unsigned int n_components, const EvaluationFlags::EvaluationFlags integration_flag, const typename MatrixFree< dim, Number >::Data *data)
static void evaluate(const unsigned int n_components, const EvaluationFlags::EvaluationFlags evaluation_flag, const typename MatrixFree< dim, Number >::Data *data)
static void evaluate(const unsigned int n_components, const EvaluationFlags::EvaluationFlags evaluation_flag, const typename MatrixFree< dim, Number >::Data *data)
static void integrate(const unsigned int n_components, const EvaluationFlags::EvaluationFlags integration_flag, const typename MatrixFree< dim, Number >::Data *data)