16#ifndef dealii_cuda_fe_evaluation_h
17#define dealii_cuda_fe_evaluation_h
21#ifdef DEAL_II_COMPILER_CUDA_AWARE
31# include <deal.II/matrix_free/cuda_matrix_free.templates.h>
47 template <
int dim,
int n_po
ints_1d>
48 __device__
inline unsigned int
52 threadIdx.x % n_points_1d :
54 threadIdx.x % n_points_1d + n_points_1d * threadIdx.y :
55 threadIdx.x % n_points_1d +
56 n_points_1d * (threadIdx.y + n_points_1d * threadIdx.z));
87 int n_q_points_1d = fe_degree + 1,
88 int n_components_ = 1,
89 typename Number =
double>
166 evaluate(
const bool evaluate_val,
const bool evaluate_grad);
176 integrate(
const bool integrate_val,
const bool integrate_grad);
185 get_value(
const unsigned int q_point)
const;
293 template <
typename Functor>
309 template <
typename Functor>
344 , padding_length(data->padding_length)
345 , mf_object_id(data->id)
346 , constraint_mask(data->constraint_mask[cell_id])
347 , use_coloring(data->use_coloring)
354 for (
unsigned int i = 0; i < dim; ++i)
369 static_assert(n_components_ == 1,
"This function only supports FE with one \
371 const unsigned int idx = internal::compute_index<dim, n_q_points_1d>();
375 values[idx] = __ldg(&src[src_idx]);
379 internal::resolve_hanging_nodes<dim, fe_degree, false>(constraint_mask,
394 static_assert(n_components_ == 1,
"This function only supports FE with one \
396 internal::resolve_hanging_nodes<dim, fe_degree, true>(constraint_mask,
399 const unsigned int idx = internal::compute_index<dim, n_q_points_1d>();
404 dst[destination_idx] +=
values[idx];
406 atomicAdd(&dst[destination_idx],
values[idx]);
418 const bool evaluate_val,
419 const bool evaluate_grad)
429 evaluator_tensor_product(mf_object_id);
430 if (evaluate_val ==
true && evaluate_grad ==
true)
432 evaluator_tensor_product.value_and_gradient_at_quad_pts(
values,
436 else if (evaluate_grad ==
true)
441 else if (evaluate_val ==
true)
443 evaluator_tensor_product.value_at_quad_pts(
values);
457 const bool integrate_val,
458 const bool integrate_grad)
466 evaluator_tensor_product(mf_object_id);
467 if (integrate_val ==
true && integrate_grad ==
true)
469 evaluator_tensor_product.integrate_value_and_gradient(
values,
472 else if (integrate_val ==
true)
474 evaluator_tensor_product.integrate_value(
values);
477 else if (integrate_grad ==
true)
497 const unsigned int q_point)
const
517 const unsigned int q_point = internal::compute_index<dim, n_q_points_1d>();
554 const unsigned int dof = internal::compute_index<dim, fe_degree + 1>();
569 values[q_point] = val_in * JxW[q_point];
583 const unsigned int q_point = internal::compute_index<dim, n_q_points_1d>();
584 values[q_point] = val_in * JxW[q_point];
612 const unsigned int dof = internal::compute_index<dim, fe_degree + 1>();
627 Number>::gradient_type
631 static_assert(n_components_ == 1,
"This function only supports FE with one \
634 const Number *inv_jacobian = &inv_jac[q_point];
636 for (
int d_1 = 0; d_1 < dim; ++d_1)
639 for (
int d_2 = 0; d_2 < dim; ++d_2)
640 tmp += inv_jacobian[padding_length *
n_cells * (dim * d_2 + d_1)] *
659 Number>::gradient_type
663 static_assert(n_components_ == 1,
"This function only supports FE with one \
667 const unsigned int q_point = internal::compute_index<dim, n_q_points_1d>();
668 const Number * inv_jacobian = &inv_jac[q_point];
670 for (
int d_1 = 0; d_1 < dim; ++d_1)
673 for (
int d_2 = 0; d_2 < dim; ++d_2)
674 tmp += inv_jacobian[padding_length *
n_cells * (dim * d_2 + d_1)] *
694 const Number *inv_jacobian = &inv_jac[q_point];
695 for (
int d_1 = 0; d_1 < dim; ++d_1)
698 for (
int d_2 = 0; d_2 < dim; ++d_2)
699 tmp += inv_jacobian[
n_cells * padding_length * (dim * d_1 + d_2)] *
701 gradients[d_1][q_point] = tmp * JxW[q_point];
717 const unsigned int q_point = internal::compute_index<dim, n_q_points_1d>();
718 const Number * inv_jacobian = &inv_jac[q_point];
719 for (
int d_1 = 0; d_1 < dim; ++d_1)
722 for (
int d_2 = 0; d_2 < dim; ++d_2)
723 tmp += inv_jacobian[
n_cells * padding_length * (dim * d_1 + d_2)] *
725 gradients[d_1][q_point] = tmp * JxW[q_point];
736 template <
typename Functor>
741 func(
this, internal::compute_index<dim, n_q_points_1d>());
753 template <
typename Functor>
FEEvaluation(const unsigned int cell_id, const data_type *data, SharedData< dim, Number > *shdata)
static constexpr unsigned int tensor_dofs_per_cell
void integrate(const bool integrate_val, const bool integrate_grad)
void apply_quad_point_operations(const Functor &func)
void submit_value(const value_type &val_in, const unsigned int q_point)
const unsigned int mf_object_id
types::global_dof_index * local_to_global
void distribute_local_to_global(Number *dst) const
static constexpr unsigned int dimension
void submit_gradient(const gradient_type &grad_in, const unsigned int q_point)
void evaluate(const bool evaluate_val, const bool evaluate_grad)
value_type get_value() const
static constexpr unsigned int n_q_points
void submit_dof_value(const value_type &val_in, const unsigned int dof)
const unsigned int constraint_mask
static constexpr unsigned int n_components
value_type get_dof_value() const
typename MatrixFree< dim, Number >::Data data_type
void read_dof_values(const Number *src)
unsigned int padding_length
gradient_type get_gradient() const
void apply_for_each_quad_point(const Functor &func)
#define DEAL_II_DEPRECATED
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
unsigned int compute_index()
constexpr T pow(const T base, const int iexp)
unsigned int n_cells(const internal::TriangulationImplementation::NumberCache< 1 > &c)