Reference documentation for deal.II version GIT relicensing-487-ge9eb5ab491 2024-04-25 07:20:02+00:00
\(\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{
65 template <int dim,
66 int fe_degree,
67 int n_q_points_1d = fe_degree + 1,
68 int n_components_ = 1,
69 typename Number = double>
71 {
72 public:
76 using value_type = Number;
77
82
87
91 static constexpr unsigned int dimension = dim;
92
96 static constexpr unsigned int n_components = n_components_;
97
101 static constexpr unsigned int n_q_points =
102 Utilities::pow(n_q_points_1d, dim);
103
107 static constexpr unsigned int tensor_dofs_per_cell =
108 Utilities::pow(fe_degree + 1, dim);
109
115
125 read_dof_values(const Number *src);
126
134 distribute_local_to_global(Number *dst) const;
135
144 evaluate(const EvaluationFlags::EvaluationFlags evaluate_flag);
145
153 DEAL_II_DEPRECATED_EARLY_WITH_COMMENT(
154 "Use the version taking EvaluationFlags.")
156 void
157 evaluate(const bool evaluate_val, const bool evaluate_grad);
158
167 integrate(const EvaluationFlags::EvaluationFlags integration_flag);
168
176 DEAL_II_DEPRECATED_EARLY_WITH_COMMENT(
177 "Use the version taking EvaluationFlags.")
179 void
180 integrate(const bool integrate_val, const bool integrate_grad);
181
187 get_value(int q_point) const;
188
194 get_dof_value(int q_point) const;
195
201 submit_value(const value_type &val_in, int q_point);
202
208 submit_dof_value(const value_type &val_in, int q_point);
209
215 get_gradient(int q_point) const;
216
222 submit_gradient(const gradient_type &grad_in, int q_point);
223
224 // clang-format off
235 // clang-format on
236 template <typename Functor>
238 apply_for_each_quad_point(const Functor &func);
239
240 private:
242 SharedData<dim, Number> *shared_data;
244 };
245
246
247
248 template <int dim,
249 int fe_degree,
250 int n_q_points_1d,
251 int n_components_,
252 typename Number>
254 FEEvaluation<dim, fe_degree, n_q_points_1d, n_components_, Number>::
255 FEEvaluation(const data_type *data, SharedData<dim, Number> *shdata)
256 : data(data)
257 , shared_data(shdata)
258 , cell_id(shared_data->team_member.league_rank())
259 {}
260
261
262
263 template <int dim,
264 int fe_degree,
265 int n_q_points_1d,
266 int n_components_,
267 typename Number>
270 read_dof_values(const Number *src)
271 {
272 static_assert(n_components_ == 1, "This function only supports FE with one \
273 components");
274 // Populate the scratch memory
275 Kokkos::parallel_for(Kokkos::TeamThreadRange(shared_data->team_member,
276 n_q_points),
277 [&](const int &i) {
278 shared_data->values(i) =
279 src[data->local_to_global(cell_id, i)];
280 });
281 shared_data->team_member.team_barrier();
282
283 internal::resolve_hanging_nodes<dim, fe_degree, false>(
284 shared_data->team_member,
285 data->constraint_weights,
286 data->constraint_mask(cell_id),
287 shared_data->values);
288 }
289
290
291
292 template <int dim,
293 int fe_degree,
294 int n_q_points_1d,
295 int n_components_,
296 typename Number>
299 distribute_local_to_global(Number *dst) const
300 {
301 static_assert(n_components_ == 1, "This function only supports FE with one \
302 components");
303
304 internal::resolve_hanging_nodes<dim, fe_degree, true>(
305 shared_data->team_member,
306 data->constraint_weights,
307 data->constraint_mask(cell_id),
308 shared_data->values);
309
310 if (data->use_coloring)
311 {
312 Kokkos::parallel_for(Kokkos::TeamThreadRange(shared_data->team_member,
313 n_q_points),
314 [&](const int &i) {
315 dst[data->local_to_global(cell_id, i)] +=
316 shared_data->values(i);
317 });
318 }
319 else
320 {
321 Kokkos::parallel_for(
322 Kokkos::TeamThreadRange(shared_data->team_member, n_q_points),
323 [&](const int &i) {
324 Kokkos::atomic_add(&dst[data->local_to_global(cell_id, i)],
325 shared_data->values(i));
326 });
327 }
328 }
329
330
331
332 template <int dim,
333 int fe_degree,
334 int n_q_points_1d,
335 int n_components_,
336 typename Number>
339 const EvaluationFlags::EvaluationFlags evaluate_flag)
340 {
341 // First evaluate the gradients because it requires values that will be
342 // changed if evaluate_val is true
345 dim,
346 fe_degree,
347 n_q_points_1d,
348 Number>
349 evaluator_tensor_product(shared_data->team_member,
350 data->shape_values,
351 data->shape_gradients,
352 data->co_shape_gradients);
353
354 if ((evaluate_flag & EvaluationFlags::values) &&
355 (evaluate_flag & EvaluationFlags::gradients))
356 {
357 evaluator_tensor_product.evaluate_values_and_gradients(
358 shared_data->values, shared_data->gradients);
359 shared_data->team_member.team_barrier();
360 }
361 else if (evaluate_flag & EvaluationFlags::gradients)
362 {
363 evaluator_tensor_product.evaluate_gradients(shared_data->values,
364 shared_data->gradients);
365 shared_data->team_member.team_barrier();
366 }
367 else if (evaluate_flag & EvaluationFlags::values)
368 {
369 evaluator_tensor_product.evaluate_values(shared_data->values);
370 shared_data->team_member.team_barrier();
371 }
372 }
373
374
375
376 template <int dim,
377 int fe_degree,
378 int n_q_points_1d,
379 int n_components_,
380 typename Number>
383 const bool evaluate_val,
384 const bool evaluate_grad)
385 {
386 evaluate(
389 }
390
391
392
393 template <int dim,
394 int fe_degree,
395 int n_q_points_1d,
396 int n_components_,
397 typename Number>
400 const EvaluationFlags::EvaluationFlags integration_flag)
401 {
404 dim,
405 fe_degree,
406 n_q_points_1d,
407 Number>
408 evaluator_tensor_product(shared_data->team_member,
409 data->shape_values,
410 data->shape_gradients,
411 data->co_shape_gradients);
412
413 if ((integration_flag & EvaluationFlags::values) &&
414 (integration_flag & EvaluationFlags::gradients))
415 {
416 evaluator_tensor_product.integrate_values_and_gradients(
417 shared_data->values, shared_data->gradients);
418 }
419 else if (integration_flag & EvaluationFlags::values)
420 {
421 evaluator_tensor_product.integrate_values(shared_data->values);
422 shared_data->team_member.team_barrier();
423 }
424 else if (integration_flag & EvaluationFlags::gradients)
425 {
426 evaluator_tensor_product.template integrate_gradients<false>(
427 shared_data->values, shared_data->gradients);
428 shared_data->team_member.team_barrier();
429 }
430 }
431
432
433
434 template <int dim,
435 int fe_degree,
436 int n_q_points_1d,
437 int n_components_,
438 typename Number>
441 const bool integrate_val,
442 const bool integrate_grad)
443 {
444 integrate(
447 }
448
449
450
451 template <int dim,
452 int fe_degree,
453 int n_q_points_1d,
454 int n_components_,
455 typename Number>
457 fe_degree,
458 n_q_points_1d,
459 n_components_,
460 Number>::value_type
462 int q_point) const
463 {
464 return shared_data->values(q_point);
465 }
466
467
468
469 template <int dim,
470 int fe_degree,
471 int n_q_points_1d,
472 int n_components_,
473 typename Number>
475 fe_degree,
476 n_q_points_1d,
477 n_components_,
478 Number>::value_type
480 get_dof_value(int q_point) const
481 {
482 return shared_data->values(q_point);
483 }
484
485
486
487 template <int dim,
488 int fe_degree,
489 int n_q_points_1d,
490 int n_components_,
491 typename Number>
494 submit_value(const value_type &val_in, int q_point)
495 {
496 shared_data->values(q_point) = val_in * data->JxW(cell_id, q_point);
497 }
498
499
500
501 template <int dim,
502 int fe_degree,
503 int n_q_points_1d,
504 int n_components_,
505 typename Number>
508 submit_dof_value(const value_type &val_in, int q_point)
509 {
510 shared_data->values(q_point) = val_in;
511 }
512
513
514
515 template <int dim,
516 int fe_degree,
517 int n_q_points_1d,
518 int n_components_,
519 typename Number>
521 fe_degree,
522 n_q_points_1d,
523 n_components_,
524 Number>::gradient_type
526 get_gradient(int q_point) const
527 {
528 static_assert(n_components_ == 1, "This function only supports FE with one \
529 components");
530
531 gradient_type grad;
532 for (unsigned int d_1 = 0; d_1 < dim; ++d_1)
533 {
534 Number tmp = 0.;
535 for (unsigned int d_2 = 0; d_2 < dim; ++d_2)
536 tmp += data->inv_jacobian(cell_id, q_point, d_2, d_1) *
537 shared_data->gradients(q_point, d_2);
538 grad[d_1] = tmp;
539 }
540
541 return grad;
542 }
543
544
545
546 template <int dim,
547 int fe_degree,
548 int n_q_points_1d,
549 int n_components_,
550 typename Number>
553 submit_gradient(const gradient_type &grad_in, int q_point)
554 {
555 for (unsigned int d_1 = 0; d_1 < dim; ++d_1)
556 {
557 Number tmp = 0.;
558 for (unsigned int d_2 = 0; d_2 < dim; ++d_2)
559 tmp += data->inv_jacobian(cell_id, q_point, d_1, d_2) * grad_in[d_2];
560 shared_data->gradients(q_point, d_1) =
561 tmp * data->JxW(cell_id, q_point);
562 }
563 }
564
565
566
567 template <int dim,
568 int fe_degree,
569 int n_q_points_1d,
570 int n_components_,
571 typename Number>
572 template <typename Functor>
575 apply_for_each_quad_point(const Functor &func)
576 {
577 Kokkos::parallel_for(Kokkos::TeamThreadRange(shared_data->team_member,
578 n_q_points),
579 [&](const int &i) { func(this, i); });
580 shared_data->team_member.team_barrier();
581 }
582
583
584
585#ifndef DOXYGEN
586 template <int dim,
587 int fe_degree,
588 int n_q_points_1d,
589 int n_components_,
590 typename Number>
591 constexpr unsigned int
594#endif
595} // namespace Portable
596
598
599#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
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
Tensor< 1, dim, Number > gradient_type
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:502
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:503
The namespace for the EvaluationFlags enum.
EvaluationFlags
The EvaluationFlags enum.
constexpr T pow(const T base, const int iexp)
Definition utilities.h:966
#define DEAL_II_HOST_DEVICE
Definition numbers.h:34