deal.II version GIT relicensing-2306-ga6fff1538f 2024-12-29 17:30:00+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
vector_tools_evaluate.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2021 - 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
16#ifndef dealii_vector_tools_evaluation_h
17#define dealii_vector_tools_evaluation_h
18
19#include <deal.II/base/config.h>
20
22
24
26
28#include <deal.II/lac/vector.h>
29
31
33
34namespace VectorTools
35{
40 {
45 {
49 avg = 0,
55 max = 1,
61 min = 2,
65 insert = 3
66 };
67 } // namespace EvaluationFlags
68
139 template <int n_components,
140 template <int, int>
141 class MeshType,
142 int dim,
143 int spacedim,
144 typename VectorType>
147 concepts::is_triangulation_or_dof_handler<MeshType<dim, spacedim>>))
148 std::vector<
149 typename FEPointEvaluation<n_components,
150 dim,
151 spacedim,
152 typename VectorType::value_type>::
153 value_type> point_values(const Mapping<dim> &mapping,
154 const MeshType<dim, spacedim> &mesh,
155 const VectorType &vector,
156 const std::vector<Point<spacedim>>
157 &evaluation_points,
158 Utilities::MPI::RemotePointEvaluation<dim,
159 spacedim>
160 &cache,
161 const EvaluationFlags::EvaluationFlags flags =
162 EvaluationFlags::avg,
163 const unsigned int first_selected_component = 0);
164
180 template <int n_components,
181 template <int, int>
182 class MeshType,
183 int dim,
184 int spacedim,
185 typename VectorType>
187 (concepts::is_dealii_vector_type<VectorType> &&
188 concepts::is_triangulation_or_dof_handler<MeshType<dim, spacedim>>))
189 std::vector<
190 typename FEPointEvaluation<n_components,
191 dim,
192 spacedim,
193 typename VectorType::value_type>::
194 value_type> point_values(const Utilities::MPI::
195 RemotePointEvaluation<dim, spacedim> &cache,
196 const MeshType<dim, spacedim> &mesh,
197 const VectorType &vector,
198 const EvaluationFlags::EvaluationFlags flags =
199 EvaluationFlags::avg,
200 const unsigned int first_selected_component = 0);
201
215 template <int n_components,
216 template <int, int>
217 class MeshType,
218 int dim,
219 int spacedim,
220 typename VectorType>
222 (concepts::is_dealii_vector_type<VectorType> &&
223 concepts::is_triangulation_or_dof_handler<MeshType<dim, spacedim>>))
224 std::vector<
225 typename FEPointEvaluation<n_components,
226 dim,
227 spacedim,
228 typename VectorType::value_type>::
229 gradient_type> point_gradients(const Mapping<dim> &mapping,
230 const MeshType<dim, spacedim> &mesh,
231 const VectorType &vector,
232 const std::vector<Point<spacedim>>
233 &evaluation_points,
234 Utilities::MPI::RemotePointEvaluation<
235 dim,
236 spacedim> &cache,
238 flags = EvaluationFlags::avg,
239 const unsigned int
240 first_selected_component = 0);
241
256 template <int n_components,
257 template <int, int>
258 class MeshType,
259 int dim,
260 int spacedim,
261 typename VectorType>
263 (concepts::is_dealii_vector_type<VectorType> &&
264 concepts::is_triangulation_or_dof_handler<MeshType<dim, spacedim>>))
265 std::vector<
266 typename FEPointEvaluation<n_components,
267 dim,
268 spacedim,
269 typename VectorType::value_type>::
270 gradient_type> point_gradients(const Utilities::MPI::
271 RemotePointEvaluation<dim, spacedim>
272 &cache,
273 const MeshType<dim, spacedim> &mesh,
274 const VectorType &vector,
276 flags = EvaluationFlags::avg,
277 const unsigned int
278 first_selected_component = 0);
279
280
281
282 // inlined functions
283
284
285#ifndef DOXYGEN
286 template <int n_components,
287 template <int, int>
288 class MeshType,
289 int dim,
290 int spacedim,
291 typename VectorType>
294 concepts::is_triangulation_or_dof_handler<MeshType<dim, spacedim>>))
295 inline std::vector<
296 typename FEPointEvaluation<n_components,
297 dim,
298 spacedim,
299 typename VectorType::value_type>::
300 value_type> point_values(const Mapping<dim> &mapping,
301 const MeshType<dim, spacedim> &mesh,
302 const VectorType &vector,
303 const std::vector<Point<spacedim>>
304 &evaluation_points,
305 Utilities::MPI::RemotePointEvaluation<dim,
306 spacedim>
307 &cache,
308 const EvaluationFlags::EvaluationFlags flags,
309 const unsigned int first_selected_component)
310 {
311 cache.reinit(evaluation_points, mesh.get_triangulation(), mapping);
312
313 return point_values<n_components>(
314 cache, mesh, vector, flags, first_selected_component);
315 }
316
317
318
319 template <int n_components,
320 template <int, int>
321 class MeshType,
322 int dim,
323 int spacedim,
324 typename VectorType>
327 concepts::is_triangulation_or_dof_handler<MeshType<dim, spacedim>>))
328 inline std::vector<
329 typename FEPointEvaluation<n_components,
330 dim,
331 spacedim,
332 typename VectorType::value_type>::
333 gradient_type> point_gradients(const Mapping<dim> &mapping,
334 const MeshType<dim, spacedim> &mesh,
335 const VectorType &vector,
336 const std::vector<Point<spacedim>>
337 &evaluation_points,
338 Utilities::MPI::RemotePointEvaluation<
339 dim,
340 spacedim> &cache,
342 flags,
343 const unsigned int
344 first_selected_component)
345 {
346 cache.reinit(evaluation_points, mesh.get_triangulation(), mapping);
347
348 return point_gradients<n_components>(
349 cache, mesh, vector, flags, first_selected_component);
350 }
351
352
353
354 namespace internal
355 {
359 template <typename T>
360 T
362 const ArrayView<const T> &values)
363 {
364 switch (flags)
365 {
366 case EvaluationFlags::avg:
367 {
368 return std::accumulate(values.begin(), values.end(), T{}) /
369 (T(1.0) * values.size());
370 }
371 case EvaluationFlags::max:
372 return *std::max_element(values.begin(), values.end());
373 case EvaluationFlags::min:
374 return *std::min_element(values.begin(), values.end());
375 case EvaluationFlags::insert:
376 return values[0];
377 default:
379 return values[0];
380 }
381 }
382
383
384
388 template <int rank, int dim, typename Number>
391 const ArrayView<const Tensor<rank, dim, Number>> &values)
392 {
393 switch (flags)
394 {
395 case EvaluationFlags::avg:
396 {
397 return std::accumulate(values.begin(),
398 values.end(),
400 (Number(1.0) * values.size());
401 }
402 case EvaluationFlags::insert:
403 return values[0];
404 default:
406 return values[0];
407 }
408 }
409
410
411
416 template <int n_components, int rank, int dim, typename Number>
418 reduce(
420 const ArrayView<const Tensor<1, n_components, Tensor<rank, dim, Number>>>
421 &values)
422 {
423 switch (flags)
424 {
425 case EvaluationFlags::avg:
426 {
428
429 for (unsigned int j = 0; j < values.size(); ++j)
430 for (unsigned int i = 0; i < n_components; ++i)
431 temp[i] = temp[i] + values[j][i];
432
433 for (unsigned int i = 0; i < n_components; ++i)
434 temp[i] /= Number(values.size());
435
436 return temp;
437 }
438 case EvaluationFlags::insert:
439 return values[0];
440 default:
442 return values[0];
443 }
444 }
445
446
447
448 template <int n_components,
449 int dim,
450 int spacedim,
451 typename VectorType,
452 typename value_type>
453 void
454 process_cell(
455 const unsigned int i,
456 const typename Utilities::MPI::RemotePointEvaluation<dim,
457 spacedim>::CellData
458 &cell_data,
460 const DoFHandler<dim, spacedim> &dof_handler,
461 const VectorType &vector,
462 const UpdateFlags update_flags,
463 const ::EvaluationFlags::EvaluationFlags evaluation_flags,
464 const unsigned int first_selected_component,
465 const std::function<
466 value_type(const FEPointEvaluation<n_components,
467 dim,
468 spacedim,
469 typename VectorType::value_type> &,
470 const unsigned int &)> process_quadrature_point,
471 const ArrayView<value_type> &values,
472 std::vector<typename VectorType::value_type> &solution_values,
473 std::vector<
474 std::unique_ptr<FEPointEvaluation<n_components,
475 dim,
476 spacedim,
477 typename VectorType::value_type>>>
478 &evaluators)
479 {
480 if (evaluators.empty())
481 evaluators.resize(dof_handler.get_fe_collection().size());
482
484 &cache.get_triangulation(),
485 cell_data.cells[i].first,
486 cell_data.cells[i].second,
487 &dof_handler};
488
489 const ArrayView<const Point<dim>> unit_points(
490 cell_data.reference_point_values.data() +
491 cell_data.reference_point_ptrs[i],
492 cell_data.reference_point_ptrs[i + 1] -
493 cell_data.reference_point_ptrs[i]);
494
495 solution_values.resize(
496 dof_handler.get_fe(cell->active_fe_index()).n_dofs_per_cell());
497 cell->get_dof_values(vector,
498 solution_values.begin(),
499 solution_values.end());
500
501 if (evaluators[cell->active_fe_index()] == nullptr)
502 evaluators[cell->active_fe_index()] =
503 std::make_unique<FEPointEvaluation<n_components,
504 dim,
505 spacedim,
506 typename VectorType::value_type>>(
507 cache.get_mapping(),
508 cell->get_fe(),
509 update_flags,
510 first_selected_component);
511 auto &evaluator = *evaluators[cell->active_fe_index()];
512
513 evaluator.reinit(cell, unit_points);
514 evaluator.evaluate(solution_values, evaluation_flags);
515
516 for (unsigned int q = 0; q < unit_points.size(); ++q)
517 values[q + cell_data.reference_point_ptrs[i]] =
518 process_quadrature_point(evaluator, q);
519 }
520
521
522
523 template <int dim, int spacedim, typename Number>
524 Number
525 get_value(
527 const Vector<Number> &vector,
529 {
530 (void)tria;
531 AssertDimension(tria.n_active_cells(), vector.size());
532 return vector[cell->active_cell_index()];
533 }
534
535
536
537 template <int dim, int spacedim, typename Number>
538 Number
539 get_value(
543 {
544 const auto distributed_tria =
545 dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(&tria);
546
547 const bool use_distributed_path =
548 (distributed_tria == nullptr) ?
549 false :
550 (vector.get_partitioner().get() ==
551 distributed_tria->global_active_cell_index_partitioner()
552 .lock()
553 .get());
554
555 if (use_distributed_path)
556 {
557 return vector[cell->global_active_cell_index()];
558 }
559 else
560 {
562 return vector[cell->active_cell_index()];
563 }
564 }
565
566
567
568 template <typename Number, typename Number2>
569 void
570 set_value(Number &dst, const Number2 &src)
571 {
572 dst = src;
573 }
574
575
576
577 template <typename Number, int rank, int dim, typename Number2>
578 void
579 set_value(Tensor<rank, dim, Number> &, const Number2 &)
580 {
581 Assert(false,
583 "A cell-data vector can only have a single component."));
584 }
585
586
587
588 template <int n_components,
589 int dim,
590 int spacedim,
591 typename VectorType,
592 typename value_type>
593 void
594 process_cell(
595 const unsigned int i,
596 const typename Utilities::MPI::RemotePointEvaluation<dim,
597 spacedim>::CellData
598 &cell_data,
601 const VectorType &vector,
602 const UpdateFlags,
603 const ::EvaluationFlags::EvaluationFlags evaluation_flags,
604 const unsigned int first_selected_component,
605 const std::function<
606 value_type(const FEPointEvaluation<n_components,
607 dim,
608 spacedim,
609 typename VectorType::value_type> &,
610 const unsigned int &)>,
611 const ArrayView<value_type> &values,
612 std::vector<typename VectorType::value_type> &,
613 std::vector<
614 std::unique_ptr<FEPointEvaluation<n_components,
615 dim,
616 spacedim,
617 typename VectorType::value_type>>> &)
618 {
619 (void)evaluation_flags;
620 (void)first_selected_component;
621
622 Assert(n_components == 1 && first_selected_component == 0,
624 "A cell-data vector can only have a single component."));
625
626 Assert(evaluation_flags ==
628 ExcMessage("For cell-data vectors, only values can be queried."));
629
631 &triangulation, cell_data.cells[i].first, cell_data.cells[i].second};
632
633 const auto value = get_value(triangulation, vector, cell);
634
635 for (unsigned int q = cell_data.reference_point_ptrs[i];
636 q < cell_data.reference_point_ptrs[i + 1];
637 ++q)
638 set_value(values[q], value);
639 }
640
641
642
643 template <int n_components,
644 int dim,
645 int spacedim,
646 typename MeshType,
647 typename VectorType,
648 typename value_type>
652 inline std::vector<value_type> evaluate_at_points(
654 const MeshType &mesh,
655 const VectorType &vector,
657 const unsigned int first_selected_component,
658 const UpdateFlags update_flags,
659 const ::EvaluationFlags::EvaluationFlags evaluation_flags,
660 const std::function<
661 value_type(const FEPointEvaluation<n_components,
662 dim,
663 spacedim,
664 typename VectorType::value_type> &,
665 const unsigned int &)> process_quadrature_point)
666 {
667 Assert(cache.is_ready(),
669 "Utilities::MPI::RemotePointEvaluation is not ready yet! "
670 "Please call Utilities::MPI::RemotePointEvaluation::reinit() "
671 "yourself or another function that does this for you."));
672
673 Assert(
674 &mesh.get_triangulation() == &cache.get_triangulation(),
676 "The provided Utilities::MPI::RemotePointEvaluation and DoFHandler "
677 "object have been set up with different Triangulation objects, "
678 "a scenario not supported!"));
679
680 // evaluate values at points if possible
681 const auto evaluation_point_results = [&]() {
682 // helper function for accessing the global vector and interpolating
683 // the results onto the points
684 const auto evaluation_function = [&](auto &values,
685 const auto &cell_data) {
686 std::vector<typename VectorType::value_type> solution_values;
687
688 std::vector<
689 std::unique_ptr<FEPointEvaluation<n_components,
690 dim,
691 spacedim,
692 typename VectorType::value_type>>>
693 evaluators;
694
695 for (unsigned int i = 0; i < cell_data.cells.size(); ++i)
696 process_cell<n_components, dim, spacedim, VectorType, value_type>(
697 i,
698 cell_data,
699 cache,
700 mesh,
701 vector,
702 update_flags,
703 evaluation_flags,
704 first_selected_component,
705 process_quadrature_point,
706 values,
707 solution_values,
708 evaluators);
709 };
710
711 std::vector<value_type> evaluation_point_results;
712 std::vector<value_type> buffer;
713
714 cache.template evaluate_and_process<value_type>(
715 evaluation_point_results, buffer, evaluation_function);
716
717 return evaluation_point_results;
718 }();
719
720 if (cache.is_map_unique())
721 {
722 // each point has exactly one result (unique map)
723 return evaluation_point_results;
724 }
725 else
726 {
727 // map is not unique (multiple or no results): postprocessing is
728 // needed
729 std::vector<value_type> unique_evaluation_point_results(
730 cache.get_point_ptrs().size() - 1);
731
732 const auto &ptr = cache.get_point_ptrs();
733
734 for (unsigned int i = 0; i < ptr.size() - 1; ++i)
735 {
736 const auto n_entries = ptr[i + 1] - ptr[i];
737 if (n_entries == 0)
738 continue;
739
740 unique_evaluation_point_results[i] =
741 reduce(flags,
743 evaluation_point_results.data() + ptr[i], n_entries));
744 }
745
746 return unique_evaluation_point_results;
747 }
748 }
749 } // namespace internal
750
751 template <int n_components,
752 template <int, int>
753 class MeshType,
754 int dim,
755 int spacedim,
756 typename VectorType>
759 concepts::is_triangulation_or_dof_handler<MeshType<dim, spacedim>>))
760 inline std::vector<
761 typename FEPointEvaluation<n_components,
762 dim,
763 spacedim,
764 typename VectorType::value_type>::
765 value_type> point_values(const Utilities::MPI::
766 RemotePointEvaluation<dim, spacedim> &cache,
767 const MeshType<dim, spacedim> &mesh,
768 const VectorType &vector,
769 const EvaluationFlags::EvaluationFlags flags,
770 const unsigned int first_selected_component)
771 {
772 return internal::evaluate_at_points<
773 n_components,
774 dim,
775 spacedim,
776 MeshType<dim, spacedim>,
778 typename FEPointEvaluation<n_components,
779 dim,
780 spacedim,
781 typename VectorType::value_type>::value_type>(
782 cache,
783 mesh,
784 vector,
785 flags,
786 first_selected_component,
789 [](const auto &evaluator, const auto &q) {
790 return evaluator.get_value(q);
791 });
792 }
793
794
795
796 template <int n_components,
797 template <int, int>
798 class MeshType,
799 int dim,
800 int spacedim,
801 typename VectorType>
804 concepts::is_triangulation_or_dof_handler<MeshType<dim, spacedim>>))
805 inline std::vector<
806 typename FEPointEvaluation<n_components,
807 dim,
808 spacedim,
809 typename VectorType::value_type>::
810 gradient_type> point_gradients(const Utilities::MPI::
811 RemotePointEvaluation<dim, spacedim>
812 &cache,
813 const MeshType<dim, spacedim> &mesh,
814 const VectorType &vector,
816 flags,
817 const unsigned int
818 first_selected_component)
819 {
820 return internal::evaluate_at_points<
821 n_components,
822 dim,
823 spacedim,
824 MeshType<dim, spacedim>,
826 typename FEPointEvaluation<
827 n_components,
828 dim,
829 spacedim,
830 typename VectorType::value_type>::gradient_type>(
831 cache,
832 mesh,
833 vector,
834 flags,
835 first_selected_component,
838 [](const auto &evaluator, const unsigned &q) {
839 return evaluator.get_gradient(q);
840 });
841 }
842
843#endif
844} // namespace VectorTools
845
847
848#endif // dealii_vector_tools_boundary_h
const hp::FECollection< dim, spacedim > & get_fe_collection() const
const FiniteElement< dim, spacedim > & get_fe(const types::fe_index index=0) const
void reinit(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const ArrayView< const Point< dim > > &unit_points)
unsigned int n_dofs_per_cell() const
size_type locally_owned_size() const
Abstract base class for mapping classes.
Definition mapping.h:318
Definition point.h:111
unsigned int n_active_cells() const
const Triangulation< dim, spacedim > & get_triangulation() const
const std::vector< unsigned int > & get_point_ptrs() const
const Mapping< dim, spacedim > & get_mapping() const
virtual size_type size() const override
unsigned int size() const
Definition collection.h:315
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:498
#define DEAL_II_CXX20_REQUIRES(condition)
Definition config.h:175
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcMessage(std::string arg1)
typename ActiveSelector::active_cell_iterator active_cell_iterator
UpdateFlags
@ update_values
Shape function values.
@ update_gradients
Shape function gradients.
#define DEAL_II_NOT_IMPLEMENTED()
The namespace for the EvaluationFlags enum.
EvaluationFlags
The EvaluationFlags enum.
Tpetra::Vector< Number, LO, GO, NodeType< MemorySpace > > VectorType
T reduce(const T &local_value, const MPI_Comm comm, const std::function< T(const T &, const T &)> &combiner, const unsigned int root_process=0)
std::vector< typename FEPointEvaluation< n_components, dim, spacedim, typename VectorType::value_type >::value_type > point_values(const Mapping< dim > &mapping, const MeshType< dim, spacedim > &mesh, const VectorType &vector, const std::vector< Point< spacedim > > &evaluation_points, Utilities::MPI::RemotePointEvaluation< dim, spacedim > &cache, const EvaluationFlags::EvaluationFlags flags=EvaluationFlags::avg, const unsigned int first_selected_component=0)
std::vector< typename FEPointEvaluation< n_components, dim, spacedim, typename VectorType::value_type >::gradient_type > point_gradients(const Mapping< dim > &mapping, const MeshType< dim, spacedim > &mesh, const VectorType &vector, const std::vector< Point< spacedim > > &evaluation_points, Utilities::MPI::RemotePointEvaluation< dim, spacedim > &cache, const EvaluationFlags::EvaluationFlags flags=EvaluationFlags::avg, const unsigned int first_selected_component=0)
STL namespace.
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation