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
51 return (dim == 1 ? threadIdx.x % n_points_1d :
52 dim == 2 ? threadIdx.x % n_points_1d + n_points_1d * threadIdx.y :
53 threadIdx.x % n_points_1d +
55 (threadIdx.y + n_points_1d * threadIdx.z));
86 int n_q_points_1d = fe_degree + 1,
87 int n_components_ = 1,
88 typename Number =
double>
165 evaluate(
const bool evaluate_val,
const bool evaluate_grad);
175 integrate(
const bool integrate_val,
const bool integrate_grad);
231 template <
typename Functor>
241 const ::internal::MatrixFreeFunctions::ConstraintKinds
266 : n_cells(data->n_cells)
267 , padding_length(data->padding_length)
268 , mf_object_id(data->id)
269 , constraint_mask(data->constraint_mask[cell_id])
270 , use_coloring(data->use_coloring)
271 , values(shdata->values)
277 for (
unsigned int i = 0; i < dim; ++i)
292 static_assert(n_components_ == 1,
"This function only supports FE with one \
294 const unsigned int idx = internal::compute_index<dim, n_q_points_1d>();
298 values[idx] = __ldg(&src[src_idx]);
302 internal::resolve_hanging_nodes<dim, fe_degree, false>(constraint_mask,
317 static_assert(n_components_ == 1,
"This function only supports FE with one \
319 internal::resolve_hanging_nodes<dim, fe_degree, true>(constraint_mask,
322 const unsigned int idx = internal::compute_index<dim, n_q_points_1d>();
327 dst[destination_idx] += values[idx];
329 atomicAdd(&dst[destination_idx], values[idx]);
341 const bool evaluate_val,
342 const bool evaluate_grad)
352 evaluator_tensor_product(mf_object_id);
353 if (evaluate_val ==
true && evaluate_grad ==
true)
355 evaluator_tensor_product.value_and_gradient_at_quad_pts(values,
359 else if (evaluate_grad ==
true)
361 evaluator_tensor_product.gradient_at_quad_pts(values, gradients);
364 else if (evaluate_val ==
true)
366 evaluator_tensor_product.value_at_quad_pts(values);
380 const bool integrate_val,
381 const bool integrate_grad)
389 evaluator_tensor_product(mf_object_id);
390 if (integrate_val ==
true && integrate_grad ==
true)
392 evaluator_tensor_product.integrate_value_and_gradient(values,
395 else if (integrate_val ==
true)
397 evaluator_tensor_product.integrate_value(values);
400 else if (integrate_grad ==
true)
402 evaluator_tensor_product.integrate_gradient<
false>(values, gradients);
422 const unsigned int q_point = internal::compute_index<dim, n_q_points_1d>();
423 return values[q_point];
441 const unsigned int dof = internal::compute_index<dim, fe_degree + 1>();
456 const unsigned int q_point = internal::compute_index<dim, n_q_points_1d>();
457 values[q_point] = val_in * JxW[q_point];
471 const unsigned int dof = internal::compute_index<dim, fe_degree + 1>();
472 values[dof] = val_in;
486 Number>::gradient_type
490 static_assert(n_components_ == 1,
"This function only supports FE with one \
494 const unsigned int q_point = internal::compute_index<dim, n_q_points_1d>();
495 const Number * inv_jacobian = &inv_jac[q_point];
497 for (
unsigned int d_1 = 0; d_1 < dim; ++d_1)
500 for (
unsigned int d_2 = 0; d_2 < dim; ++d_2)
501 tmp += inv_jacobian[padding_length * n_cells * (dim * d_2 + d_1)] *
502 gradients[d_2][q_point];
521 const unsigned int q_point = internal::compute_index<dim, n_q_points_1d>();
522 const Number * inv_jacobian = &inv_jac[q_point];
523 for (
unsigned int d_1 = 0; d_1 < dim; ++d_1)
526 for (
unsigned int d_2 = 0; d_2 < dim; ++d_2)
527 tmp += inv_jacobian[n_cells * padding_length * (dim * d_1 + d_2)] *
529 gradients[d_1][q_point] = tmp * JxW[q_point];
540 template <
typename Functor>
static constexpr unsigned int tensor_dofs_per_cell
void integrate(const bool integrate_val, const bool integrate_grad)
void submit_dof_value(const value_type &val_in)
void submit_gradient(const gradient_type &grad_in)
const unsigned int mf_object_id
void submit_value(const value_type &val_in)
types::global_dof_index * local_to_global
void distribute_local_to_global(Number *dst) const
static constexpr unsigned int dimension
const ::internal::MatrixFreeFunctions::ConstraintKinds constraint_mask
void evaluate(const bool evaluate_val, const bool evaluate_grad)
value_type get_value() const
static constexpr unsigned int n_q_points
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)
friend class FEEvaluation
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
unsigned int compute_index()
constexpr T pow(const T base, const int iexp)