Reference documentation for deal.II version 9.6.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
portable_fe_evaluation.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2023 - 2024 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_portable_fe_evaluation_h
16#define dealii_portable_fe_evaluation_h
17
18#include <deal.II/base/config.h>
19
21#include <deal.II/base/tensor.h>
23
27#include <deal.II/matrix_free/portable_matrix_free.templates.h>
29
30#include <Kokkos_Core.hpp>
31
33
37namespace Portable
38{
64 template <int dim,
65 int fe_degree,
66 int n_q_points_1d = fe_degree + 1,
67 int n_components_ = 1,
68 typename Number = double>
70 {
71 public:
75 using value_type = Number;
76
81
86
90 static constexpr unsigned int dimension = dim;
91
95 static constexpr unsigned int n_components = n_components_;
96
100 static constexpr unsigned int n_q_points =
101 Utilities::pow(n_q_points_1d, dim);
102
106 static constexpr unsigned int tensor_dofs_per_cell =
107 Utilities::pow(fe_degree + 1, dim);
108
114
124 read_dof_values(const Number *src);
125
133 distribute_local_to_global(Number *dst) const;
134
143 evaluate(const EvaluationFlags::EvaluationFlags evaluate_flag);
144
152 DEAL_II_DEPRECATED_EARLY_WITH_COMMENT(
153 "Use the version taking EvaluationFlags.")
155 void
156 evaluate(const bool evaluate_val, const bool evaluate_grad);
157
166 integrate(const EvaluationFlags::EvaluationFlags integration_flag);
167
175 DEAL_II_DEPRECATED_EARLY_WITH_COMMENT(
176 "Use the version taking EvaluationFlags.")
178 void
179 integrate(const bool integrate_val, const bool integrate_grad);
180
186 get_value(int q_point) const;
187
193 get_dof_value(int q_point) const;
194
200 submit_value(const value_type &val_in, int q_point);
201
207 submit_dof_value(const value_type &val_in, int q_point);
208
214 get_gradient(int q_point) const;
215
221 submit_gradient(const gradient_type &grad_in, int q_point);
222
223 // clang-format off
234 // clang-format on
235 template <typename Functor>
237 apply_for_each_quad_point(const Functor &func);
238
239 private:
241 SharedData<dim, Number> *shared_data;
243 };
244
245
246
247 template <int dim,
248 int fe_degree,
249 int n_q_points_1d,
250 int n_components_,
251 typename Number>
253 FEEvaluation<dim, fe_degree, n_q_points_1d, n_components_, Number>::
254 FEEvaluation(const data_type *data, SharedData<dim, Number> *shdata)
255 : data(data)
256 , shared_data(shdata)
257 , cell_id(shared_data->team_member.league_rank())
258 {}
259
260
261
262 template <int dim,
263 int fe_degree,
264 int n_q_points_1d,
265 int n_components_,
266 typename Number>
269 read_dof_values(const Number *src)
270 {
271 static_assert(n_components_ == 1, "This function only supports FE with one \
272 components");
273 // Populate the scratch memory
274 Kokkos::parallel_for(Kokkos::TeamThreadRange(shared_data->team_member,
275 n_q_points),
276 [&](const int &i) {
277 shared_data->values(i) =
278 src[data->local_to_global(cell_id, i)];
279 });
280 shared_data->team_member.team_barrier();
281
283 shared_data->team_member,
284 data->constraint_weights,
285 data->constraint_mask(cell_id),
286 shared_data->values);
287 }
288
289
290
291 template <int dim,
292 int fe_degree,
293 int n_q_points_1d,
294 int n_components_,
295 typename Number>
298 distribute_local_to_global(Number *dst) const
299 {
300 static_assert(n_components_ == 1, "This function only supports FE with one \
301 components");
302
304 shared_data->team_member,
305 data->constraint_weights,
306 data->constraint_mask(cell_id),
307 shared_data->values);
308
309 if (data->use_coloring)
310 {
311 Kokkos::parallel_for(Kokkos::TeamThreadRange(shared_data->team_member,
312 n_q_points),
313 [&](const int &i) {
314 dst[data->local_to_global(cell_id, i)] +=
315 shared_data->values(i);
316 });
317 }
318 else
319 {
320 Kokkos::parallel_for(
321 Kokkos::TeamThreadRange(shared_data->team_member, n_q_points),
322 [&](const int &i) {
323 Kokkos::atomic_add(&dst[data->local_to_global(cell_id, i)],
324 shared_data->values(i));
325 });
326 }
327 }
328
329
330
331 template <int dim,
332 int fe_degree,
333 int n_q_points_1d,
334 int n_components_,
335 typename Number>
338 const EvaluationFlags::EvaluationFlags evaluate_flag)
339 {
340 // First evaluate the gradients because it requires values that will be
341 // changed if evaluate_val is true
344 dim,
345 fe_degree,
346 n_q_points_1d,
347 Number>
348 evaluator_tensor_product(shared_data->team_member,
349 data->shape_values,
350 data->shape_gradients,
351 data->co_shape_gradients);
352
353 if ((evaluate_flag & EvaluationFlags::values) &&
354 (evaluate_flag & EvaluationFlags::gradients))
355 {
356 evaluator_tensor_product.evaluate_values_and_gradients(
357 shared_data->values, shared_data->gradients);
358 shared_data->team_member.team_barrier();
359 }
360 else if (evaluate_flag & EvaluationFlags::gradients)
361 {
362 evaluator_tensor_product.evaluate_gradients(shared_data->values,
363 shared_data->gradients);
364 shared_data->team_member.team_barrier();
365 }
366 else if (evaluate_flag & EvaluationFlags::values)
367 {
368 evaluator_tensor_product.evaluate_values(shared_data->values);
369 shared_data->team_member.team_barrier();
370 }
371 }
372
373
374
375 template <int dim,
376 int fe_degree,
377 int n_q_points_1d,
378 int n_components_,
379 typename Number>
382 const bool evaluate_val,
383 const bool evaluate_grad)
384 {
385 evaluate(
388 }
389
390
391
392 template <int dim,
393 int fe_degree,
394 int n_q_points_1d,
395 int n_components_,
396 typename Number>
399 const EvaluationFlags::EvaluationFlags integration_flag)
400 {
403 dim,
404 fe_degree,
405 n_q_points_1d,
406 Number>
407 evaluator_tensor_product(shared_data->team_member,
408 data->shape_values,
409 data->shape_gradients,
410 data->co_shape_gradients);
411
412 if ((integration_flag & EvaluationFlags::values) &&
413 (integration_flag & EvaluationFlags::gradients))
414 {
415 evaluator_tensor_product.integrate_values_and_gradients(
416 shared_data->values, shared_data->gradients);
417 }
418 else if (integration_flag & EvaluationFlags::values)
419 {
420 evaluator_tensor_product.integrate_values(shared_data->values);
421 shared_data->team_member.team_barrier();
422 }
423 else if (integration_flag & EvaluationFlags::gradients)
424 {
425 evaluator_tensor_product.template integrate_gradients<false>(
426 shared_data->values, shared_data->gradients);
427 shared_data->team_member.team_barrier();
428 }
429 }
430
431
432
433 template <int dim,
434 int fe_degree,
435 int n_q_points_1d,
436 int n_components_,
437 typename Number>
440 const bool integrate_val,
441 const bool integrate_grad)
442 {
443 integrate(
446 }
447
448
449
450 template <int dim,
451 int fe_degree,
452 int n_q_points_1d,
453 int n_components_,
454 typename Number>
456 fe_degree,
457 n_q_points_1d,
458 n_components_,
459 Number>::value_type
461 int q_point) const
462 {
463 return shared_data->values(q_point);
464 }
465
466
467
468 template <int dim,
469 int fe_degree,
470 int n_q_points_1d,
471 int n_components_,
472 typename Number>
474 fe_degree,
475 n_q_points_1d,
476 n_components_,
477 Number>::value_type
479 get_dof_value(int q_point) const
480 {
481 return shared_data->values(q_point);
482 }
483
484
485
486 template <int dim,
487 int fe_degree,
488 int n_q_points_1d,
489 int n_components_,
490 typename Number>
493 submit_value(const value_type &val_in, int q_point)
494 {
495 shared_data->values(q_point) = val_in * data->JxW(cell_id, q_point);
496 }
497
498
499
500 template <int dim,
501 int fe_degree,
502 int n_q_points_1d,
503 int n_components_,
504 typename Number>
507 submit_dof_value(const value_type &val_in, int q_point)
508 {
509 shared_data->values(q_point) = val_in;
510 }
511
512
513
514 template <int dim,
515 int fe_degree,
516 int n_q_points_1d,
517 int n_components_,
518 typename Number>
520 fe_degree,
521 n_q_points_1d,
522 n_components_,
523 Number>::gradient_type
525 get_gradient(int q_point) const
526 {
527 static_assert(n_components_ == 1, "This function only supports FE with one \
528 components");
529
530 gradient_type grad;
531 for (unsigned int d_1 = 0; d_1 < dim; ++d_1)
532 {
533 Number tmp = 0.;
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) *
536 shared_data->gradients(q_point, d_2);
537 grad[d_1] = tmp;
538 }
539
540 return grad;
541 }
542
543
544
545 template <int dim,
546 int fe_degree,
547 int n_q_points_1d,
548 int n_components_,
549 typename Number>
552 submit_gradient(const gradient_type &grad_in, int q_point)
553 {
554 for (unsigned int d_1 = 0; d_1 < dim; ++d_1)
555 {
556 Number tmp = 0.;
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];
559 shared_data->gradients(q_point, d_1) =
560 tmp * data->JxW(cell_id, q_point);
561 }
562 }
563
564
565
566 template <int dim,
567 int fe_degree,
568 int n_q_points_1d,
569 int n_components_,
570 typename Number>
571 template <typename Functor>
574 apply_for_each_quad_point(const Functor &func)
575 {
576 Kokkos::parallel_for(Kokkos::TeamThreadRange(shared_data->team_member,
577 n_q_points),
578 [&](const int &i) { func(this, i); });
579 shared_data->team_member.team_barrier();
580 }
581
582
583
584#ifndef DOXYGEN
585 template <int dim,
586 int fe_degree,
587 int n_q_points_1d,
588 int n_components_,
589 typename Number>
590 constexpr unsigned int
593#endif
594} // namespace Portable
595
597
598#endif
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
Definition config.h:503
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:504
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)
Definition utilities.h:966
#define DEAL_II_HOST_DEVICE
Definition numbers.h:34