16#ifndef dealii_vector_tools_evaluation_h
17#define dealii_vector_tools_evaluation_h
139 template <
int n_components,
152 typename VectorType::value_type>::
154 const MeshType<dim, spacedim> &mesh,
155 const VectorType &vector,
163 const
unsigned int first_selected_component = 0);
180 template <
int n_components,
187 (
concepts::is_dealii_vector_type<VectorType> &&
188 concepts::is_triangulation_or_dof_handler<MeshType<dim, spacedim>>))
193 typename VectorType::value_type>::
195 RemotePointEvaluation<dim, spacedim> &cache,
196 const MeshType<dim, spacedim> &mesh,
197 const VectorType &vector,
200 const
unsigned int first_selected_component = 0);
215 template <
int n_components,
222 (
concepts::is_dealii_vector_type<VectorType> &&
223 concepts::is_triangulation_or_dof_handler<MeshType<dim, spacedim>>))
228 typename VectorType::value_type>::
230 const MeshType<dim, spacedim> &mesh,
231 const VectorType &vector,
240 first_selected_component = 0);
256 template <
int n_components,
263 (
concepts::is_dealii_vector_type<VectorType> &&
264 concepts::is_triangulation_or_dof_handler<MeshType<dim, spacedim>>))
269 typename VectorType::value_type>::
271 RemotePointEvaluation<dim, spacedim>
273 const MeshType<dim, spacedim> &mesh,
274 const VectorType &vector,
278 first_selected_component = 0);
286 template <
int n_components,
299 typename VectorType::value_type>::
301 const MeshType<dim, spacedim> &mesh,
302 const VectorType &vector,
309 const
unsigned int first_selected_component)
311 cache.
reinit(evaluation_points, mesh.get_triangulation(), mapping);
314 cache, mesh, vector, flags, first_selected_component);
319 template <
int n_components,
332 typename VectorType::value_type>::
334 const MeshType<dim, spacedim> &mesh,
335 const VectorType &vector,
344 first_selected_component)
346 cache.
reinit(evaluation_points, mesh.get_triangulation(), mapping);
349 cache, mesh, vector, flags, first_selected_component);
359 template <
typename T>
366 case EvaluationFlags::avg:
368 return std::accumulate(values.begin(), values.end(), T{}) /
369 (T(1.0) * values.size());
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:
388 template <
int rank,
int dim,
typename Number>
395 case EvaluationFlags::avg:
397 return std::accumulate(values.begin(),
400 (Number(1.0) * values.size());
402 case EvaluationFlags::insert:
416 template <
int n_components,
int rank,
int dim,
typename Number>
425 case EvaluationFlags::avg:
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];
433 for (
unsigned int i = 0; i < n_components; ++i)
434 temp[i] /= Number(values.size());
438 case EvaluationFlags::insert:
448 template <
int n_components,
455 const unsigned int i,
461 const VectorType &vector,
463 const ::EvaluationFlags::EvaluationFlags evaluation_flags,
464 const unsigned int first_selected_component,
469 typename VectorType::value_type> &,
470 const unsigned int &)> process_quadrature_point,
472 std::vector<typename VectorType::value_type> &solution_values,
477 typename VectorType::value_type>>>
480 if (evaluators.empty())
485 cell_data.cells[i].first,
486 cell_data.cells[i].second,
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]);
495 solution_values.resize(
497 cell->get_dof_values(vector,
498 solution_values.begin(),
499 solution_values.end());
501 if (evaluators[cell->active_fe_index()] ==
nullptr)
502 evaluators[cell->active_fe_index()] =
506 typename VectorType::value_type>>(
510 first_selected_component);
511 auto &evaluator = *evaluators[cell->active_fe_index()];
513 evaluator.reinit(cell, unit_points);
514 evaluator.evaluate(solution_values, evaluation_flags);
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);
523 template <
int dim,
int spacedim,
typename Number>
532 return vector[cell->active_cell_index()];
537 template <
int dim,
int spacedim,
typename Number>
544 const auto distributed_tria =
547 const bool use_distributed_path =
548 (distributed_tria ==
nullptr) ?
550 (vector.get_partitioner().get() ==
551 distributed_tria->global_active_cell_index_partitioner()
555 if (use_distributed_path)
557 return vector[cell->global_active_cell_index()];
562 return vector[cell->active_cell_index()];
568 template <
typename Number,
typename Number2>
570 set_value(Number &dst,
const Number2 &src)
577 template <
typename Number,
int rank,
int dim,
typename Number2>
583 "A cell-data vector can only have a single component."));
588 template <
int n_components,
595 const unsigned int i,
601 const VectorType &vector,
603 const ::EvaluationFlags::EvaluationFlags evaluation_flags,
604 const unsigned int first_selected_component,
609 typename VectorType::value_type> &,
610 const unsigned int &)>,
612 std::vector<typename VectorType::value_type> &,
617 typename VectorType::value_type>>> &)
619 (void)evaluation_flags;
620 (void)first_selected_component;
622 Assert(n_components == 1 && first_selected_component == 0,
624 "A cell-data vector can only have a single component."));
626 Assert(evaluation_flags ==
628 ExcMessage(
"For cell-data vectors, only values can be queried."));
631 &
triangulation, cell_data.cells[i].first, cell_data.cells[i].second};
635 for (
unsigned int q = cell_data.reference_point_ptrs[i];
636 q < cell_data.reference_point_ptrs[i + 1];
638 set_value(values[q], value);
643 template <
int n_components,
652 inline std::vector<value_type> evaluate_at_points(
654 const MeshType &mesh,
655 const VectorType &vector,
657 const unsigned int first_selected_component,
659 const ::EvaluationFlags::EvaluationFlags evaluation_flags,
664 typename VectorType::value_type> &,
665 const unsigned int &)> process_quadrature_point)
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."));
676 "The provided Utilities::MPI::RemotePointEvaluation and DoFHandler "
677 "object have been set up with different Triangulation objects, "
678 "a scenario not supported!"));
681 const auto evaluation_point_results = [&]() {
684 const auto evaluation_function = [&](
auto &values,
685 const auto &cell_data) {
686 std::vector<typename VectorType::value_type> solution_values;
692 typename VectorType::value_type>>>
695 for (
unsigned int i = 0; i < cell_data.cells.size(); ++i)
696 process_cell<n_components, dim, spacedim, VectorType, value_type>(
704 first_selected_component,
705 process_quadrature_point,
711 std::vector<value_type> evaluation_point_results;
712 std::vector<value_type> buffer;
714 cache.template evaluate_and_process<value_type>(
715 evaluation_point_results, buffer, evaluation_function);
717 return evaluation_point_results;
723 return evaluation_point_results;
729 std::vector<value_type> unique_evaluation_point_results(
734 for (
unsigned int i = 0; i < ptr.size() - 1; ++i)
736 const auto n_entries = ptr[i + 1] - ptr[i];
740 unique_evaluation_point_results[i] =
743 evaluation_point_results.data() + ptr[i], n_entries));
746 return unique_evaluation_point_results;
751 template <
int n_components,
764 typename VectorType::value_type>::
766 RemotePointEvaluation<dim, spacedim> &cache,
767 const MeshType<dim, spacedim> &mesh,
768 const VectorType &vector,
770 const
unsigned int first_selected_component)
772 return internal::evaluate_at_points<
776 MeshType<dim, spacedim>,
781 typename VectorType::value_type>::value_type>(
786 first_selected_component,
789 [](
const auto &evaluator,
const auto &q) {
790 return evaluator.get_value(q);
796 template <
int n_components,
809 typename VectorType::value_type>::
811 RemotePointEvaluation<dim, spacedim>
813 const MeshType<dim, spacedim> &mesh,
814 const VectorType &vector,
818 first_selected_component)
820 return internal::evaluate_at_points<
824 MeshType<dim, spacedim>,
830 typename VectorType::value_type>::gradient_type>(
835 first_selected_component,
838 [](
const auto &evaluator,
const unsigned &q) {
839 return evaluator.get_gradient(q);
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.
unsigned int n_active_cells() const
const Triangulation< dim, spacedim > & get_triangulation() const
const std::vector< unsigned int > & get_point_ptrs() const
bool is_map_unique() const
const Mapping< dim, spacedim > & get_mapping() const
virtual size_type size() const override
unsigned int size() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_CXX20_REQUIRES(condition)
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcMessage(std::string arg1)
typename ActiveSelector::active_cell_iterator active_cell_iterator
@ 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.
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)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation