39template <
int dim,
int spacedim>
41 const unsigned int polynomial_degree)
42 :
MappingQ<dim, spacedim>(polynomial_degree)
43 , uses_level_info(false)
48template <
int dim,
int spacedim>
52 , support_point_cache(mapping.support_point_cache)
53 , uses_level_info(mapping.uses_level_info)
58template <
int dim,
int spacedim>
65 support_point_cache.reset();
66 clear_signal.disconnect();
71template <
int dim,
int spacedim>
72std::unique_ptr<Mapping<dim, spacedim>>
75 return std::make_unique<MappingQCache<dim, spacedim>>(*this);
80template <
int dim,
int spacedim>
89template <
int dim,
int spacedim>
103 const auto mapping_q =
105 if (mapping_q !=
nullptr && this->get_degree() == mapping_q->get_degree())
114 auto &fe_values = fe_values_all.
get();
115 if (fe_values.get() ==
nullptr)
119 std::vector<Point<dim>> quadrature_points;
121 FETools::hierarchic_to_lexicographic_numbering<dim>(
122 this->polynomial_degree))
123 quadrature_points.push_back(quadrature_gl.
point(i));
126 fe_values = std::make_unique<FEValues<dim, spacedim>>(
130 fe_values->reinit(cell);
131 return fe_values->get_quadrature_points();
138template <
int dim,
int spacedim>
149template <
int dim,
int spacedim>
155 &compute_points_on_cell)
157 clear_signal.disconnect();
159 [&]() ->
void { this->support_point_cache.reset(); });
161 support_point_cache =
162 std::make_shared<std::vector<std::vector<std::vector<Point<spacedim>>>>>(
165 (*support_point_cache)[l].resize(
triangulation.n_raw_cells(l));
173 (*support_point_cache)[cell->level()][cell->index()] =
174 compute_points_on_cell(cell);
176 (*support_point_cache)[cell->level()][cell->index()].size(),
177 Utilities::pow(this->get_degree() + 1, dim));
179 std::function<
void(
void *)>(),
185 uses_level_info =
true;
190template <
int dim,
int spacedim>
198 const bool function_describes_relative_displacement)
208 std::vector<Point<spacedim>> points;
210 const auto mapping_q =
213 if (mapping_q !=
nullptr && this->get_degree() == mapping_q->get_degree())
222 auto &fe_values = fe_values_all.
get();
223 if (fe_values.get() ==
nullptr)
227 std::vector<Point<dim>> quadrature_points;
229 FETools::hierarchic_to_lexicographic_numbering<dim>(
230 this->polynomial_degree))
231 quadrature_points.push_back(quadrature_gl.
point(i));
234 fe_values = std::make_unique<FEValues<dim, spacedim>>(
238 fe_values->reinit(cell);
239 points = fe_values->get_quadrature_points();
242 for (
auto &p : points)
243 if (function_describes_relative_displacement)
244 p += transformation_function(cell, p);
246 p = transformation_function(cell, p);
251 uses_level_info =
true;
256template <
int dim,
int spacedim>
262 const bool function_describes_relative_displacement)
269 [&](
const auto &,
const auto &point) {
271 for (
unsigned int c = 0; c < spacedim; ++c)
272 new_point[c] = transformation_function.
value(point, c);
275 function_describes_relative_displacement);
277 uses_level_info =
true;
284 template <
typename VectorType>
286 copy_locally_owned_data_from(
287 const VectorType &vector,
292 temp.
reinit(vector.locally_owned_elements());
300template <
int dim,
int spacedim>
301template <
typename VectorType>
306 const VectorType & vector,
307 const bool vector_describes_relative_displacement)
314 const unsigned int is_fe_q =
316 const unsigned int is_fe_dgq =
319 const auto lexicographic_to_hierarchic_numbering =
321 FETools::hierarchic_to_lexicographic_numbering<spacedim>(
322 this->get_degree()));
331 locally_relevant_dofs,
333 copy_locally_owned_data_from(vector, vector_ghosted);
344 const bool interpolation_of_values_is_needed =
345 ((is_fe_q || is_fe_dgq) && fe.
degree == this->get_degree()) ==
false;
352 const bool is_active_non_artificial_cell =
353 (cell_tria->is_active() == true) &&
354 (cell_tria->is_artificial() == false);
356 const typename DoFHandler<dim, spacedim>::cell_iterator cell_dofs(
357 &cell_tria->get_triangulation(),
362 const auto mapping_q =
363 dynamic_cast<const MappingQ<dim, spacedim> *>(&mapping);
367 ((vector_describes_relative_displacement ||
368 (is_active_non_artificial_cell == false)) &&
369 ((mapping_q != nullptr &&
370 this->get_degree() == mapping_q->get_degree()) ==
373 (is_active_non_artificial_cell && interpolation_of_values_is_needed) )
378 auto &fe_values = fe_values_all.get();
379 if (fe_values.get() == nullptr)
381 QGaussLobatto<dim> quadrature_gl(this->polynomial_degree + 1);
383 std::vector<Point<dim>> quadrature_points;
385 FETools::hierarchic_to_lexicographic_numbering<dim>(
386 this->polynomial_degree))
387 quadrature_points.push_back(quadrature_gl.point(i));
388 Quadrature<dim> quadrature(quadrature_points);
390 fe_values = std::make_unique<FEValues<dim, spacedim>>(
392 interpolation_of_values_is_needed ?
394 static_cast<const FiniteElement<dim, spacedim> &>(fe_nothing),
396 update_quadrature_points | update_values);
399 if (interpolation_of_values_is_needed)
400 fe_values->reinit(cell_dofs);
402 fe_values->reinit(cell_tria);
410 if (vector_describes_relative_displacement ||
411 is_active_non_artificial_cell ==
false)
413 if (mapping_q !=
nullptr &&
414 this->get_degree() == mapping_q->get_degree())
415 result = mapping_q->compute_mapping_support_points(cell_tria);
417 result = fe_values_all.
get()->get_quadrature_points();
422 if (is_active_non_artificial_cell ==
false)
428 Utilities::pow<unsigned int>(this->get_degree() + 1, dim));
432 if (interpolation_of_values_is_needed ==
false)
436 std::vector<types::global_dof_index> dof_indices(
438 cell_dofs->get_dof_indices(dof_indices);
440 for (
unsigned int i = 0; i < dof_indices.size(); ++i)
447 if (vector_describes_relative_displacement)
448 result[
id.second][
id.first] +=
449 vector_ghosted(dof_indices[i]);
451 result[
id.second][
id.first] =
452 vector_ghosted(dof_indices[i]);
457 if (vector_describes_relative_displacement)
458 result[lexicographic_to_hierarchic_numbering[
id.second]]
459 [
id.first] += vector_ghosted(dof_indices[i]);
461 result[lexicographic_to_hierarchic_numbering[
id.second]]
462 [
id.first] = vector_ghosted(dof_indices[i]);
472 auto &fe_values = fe_values_all.
get();
474 std::vector<Vector<typename VectorType::value_type>> values(
475 fe_values->n_quadrature_points,
478 fe_values->get_function_values(vector_ghosted, values);
480 for (
unsigned int q = 0; q < fe_values->n_quadrature_points; ++q)
481 for (
unsigned int c = 0; c < spacedim; ++c)
482 if (vector_describes_relative_displacement)
483 result[q][c] += values[q][c];
485 result[q][c] = values[q][c];
491 uses_level_info =
false;
496template <
int dim,
int spacedim>
497template <
typename VectorType>
503 const bool vector_describes_relative_displacement)
513 const unsigned int is_fe_q =
515 const unsigned int is_fe_dgq =
518 const auto lexicographic_to_hierarchic_numbering =
520 FETools::hierarchic_to_lexicographic_numbering<spacedim>(
521 this->get_degree()));
534 locally_relevant_dofs);
536 locally_relevant_dofs,
538 copy_locally_owned_data_from(vectors[l], vectors_ghosted[l]);
539 vectors_ghosted[l].update_ghost_values();
550 const bool interpolation_of_values_is_needed =
551 ((is_fe_q || is_fe_dgq) && fe.
degree == this->get_degree()) ==
false;
558 const bool is_non_artificial_cell =
559 cell_tria->level_subdomain_id() != numbers::artificial_subdomain_id;
561 const typename DoFHandler<dim, spacedim>::level_cell_iterator cell_dofs(
562 &cell_tria->get_triangulation(),
567 const auto mapping_q =
568 dynamic_cast<const MappingQ<dim, spacedim> *>(&mapping);
572 ((vector_describes_relative_displacement ||
573 (is_non_artificial_cell == false)) &&
574 ((mapping_q != nullptr &&
575 this->get_degree() == mapping_q->get_degree()) ==
578 (is_non_artificial_cell == true && interpolation_of_values_is_needed) )
583 auto &fe_values = fe_values_all.get();
584 if (fe_values.get() == nullptr)
586 QGaussLobatto<dim> quadrature_gl(this->polynomial_degree + 1);
588 std::vector<Point<dim>> quadrature_points;
590 FETools::hierarchic_to_lexicographic_numbering<dim>(
591 this->polynomial_degree))
592 quadrature_points.push_back(quadrature_gl.point(i));
593 Quadrature<dim> quadrature(quadrature_points);
595 fe_values = std::make_unique<FEValues<dim, spacedim>>(
597 interpolation_of_values_is_needed ?
599 static_cast<const FiniteElement<dim, spacedim> &>(fe_nothing),
601 update_quadrature_points | update_values);
604 if (interpolation_of_values_is_needed)
605 fe_values->reinit(cell_dofs);
607 fe_values->reinit(cell_tria);
615 if (vector_describes_relative_displacement ||
616 (is_non_artificial_cell ==
false))
618 if (mapping_q !=
nullptr &&
619 this->get_degree() == mapping_q->get_degree())
620 result = mapping_q->compute_mapping_support_points(cell_tria);
622 result = fe_values_all.
get()->get_quadrature_points();
627 if (is_non_artificial_cell ==
false)
633 Utilities::pow<unsigned int>(this->get_degree() + 1, dim));
637 if (interpolation_of_values_is_needed ==
false)
641 std::vector<types::global_dof_index> dof_indices(
643 cell_dofs->get_mg_dof_indices(dof_indices);
645 for (
unsigned int i = 0; i < dof_indices.size(); ++i)
652 if (vector_describes_relative_displacement)
653 result[
id.second][
id.first] +=
654 vectors_ghosted[cell_tria->level()](dof_indices[i]);
656 result[
id.second][
id.first] =
657 vectors_ghosted[cell_tria->level()](dof_indices[i]);
662 if (vector_describes_relative_displacement)
663 result[lexicographic_to_hierarchic_numbering[
id.second]]
665 vectors_ghosted[cell_tria->level()](dof_indices[i]);
667 result[lexicographic_to_hierarchic_numbering[
id.second]]
669 vectors_ghosted[cell_tria->level()](dof_indices[i]);
679 auto &fe_values = fe_values_all.
get();
681 std::vector<types::global_dof_index> dof_indices(
683 cell_dofs->get_mg_dof_indices(dof_indices);
685 std::vector<typename VectorType::value_type> dof_values(
689 dof_values[i] = vectors_ghosted[cell_tria->level()](dof_indices[i]);
691 for (
unsigned int c = 0; c < spacedim; ++c)
693 for (
unsigned int q = 0; q < fe_values->n_quadrature_points; ++q)
694 if (vector_describes_relative_displacement ==
false && i == 0)
696 dof_values[i] * fe_values->shape_value_component(i, q, c);
699 dof_values[i] * fe_values->shape_value_component(i, q, c);
705 uses_level_info =
true;
710template <
int dim,
int spacedim>
714 if (support_point_cache.get() !=
nullptr)
715 return sizeof(*this) +
718 return sizeof(*this);
723template <
int dim,
int spacedim>
724std::vector<Point<spacedim>>
728 Assert(support_point_cache.get() !=
nullptr,
729 ExcMessage(
"Must call MappingQCache::initialize() before "
730 "using it or after mesh has changed!"));
735 AssertIndexRange(cell->index(), (*support_point_cache)[cell->level()].size());
736 return (*support_point_cache)[cell->level()][cell->index()];
741template <
int dim,
int spacedim>
742boost::container::small_vector<Point<spacedim>,
747 Assert(support_point_cache.get() !=
nullptr,
748 ExcMessage(
"Must call MappingQCache::initialize() before "
749 "using it or after mesh has changed!"));
754 AssertIndexRange(cell->index(), (*support_point_cache)[cell->level()].size());
755 const auto ptr = (*support_point_cache)[cell->level()][cell->index()].begin();
756 return boost::container::small_vector<Point<spacedim>,
758 ptr, ptr + cell->n_vertices());
764#include "mapping_q_cache.inst"
const FiniteElement< dim, spacedim > & get_fe(const unsigned int index=0) const
const hp::FECollection< dim, spacedim > & get_fe_collection() const
const IndexSet & locally_owned_mg_dofs(const unsigned int level) const
const Triangulation< dim, spacedim > & get_triangulation() const
const IndexSet & locally_owned_dofs() const
MPI_Comm get_communicator() const
const unsigned int degree
unsigned int n_dofs_per_cell() const
virtual const FiniteElement< dim, spacedim > & base_element(const unsigned int index) const
std::pair< unsigned int, unsigned int > system_to_component_index(const unsigned int index) const
unsigned int element_multiplicity(const unsigned int index) const
unsigned int n_base_elements() const
const unsigned int n_components
virtual RangeNumberType value(const Point< dim > &p, const unsigned int component=0) const
unsigned int max_level() const
unsigned int min_level() const
std::size_t memory_consumption() const
virtual boost::container::small_vector< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > get_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const override
virtual bool preserves_vertex_locations() const override
virtual std::vector< Point< spacedim > > compute_mapping_support_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const override
virtual std::unique_ptr< Mapping< dim, spacedim > > clone() const override
void initialize(const Mapping< dim, spacedim > &mapping, const Triangulation< dim, spacedim > &triangulation)
virtual std::vector< Point< spacedim > > compute_mapping_support_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
friend class MappingQCache
Abstract base class for mapping classes.
static unsigned int n_threads()
const Point< dim > & point(const unsigned int i) const
A class that provides a separate storage location on each thread that accesses the object.
virtual unsigned int n_global_levels() const
unsigned int size() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
@ update_quadrature_points
Transformed quadrature points.
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
void import(const ::Vector< Number > &vec, VectorOperation::values operation, const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &communication_pattern={})
void update_ghost_values() const
virtual void reinit(const size_type size, const bool omit_zeroing_entries=false)
void reinit(const size_type size, const bool omit_zeroing_entries=false)
void import(const Vector< Number, MemorySpace2 > &src, VectorOperation::values operation)
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
std::vector< Integer > invert_permutation(const std::vector< Integer > &permutation)
void run(const std::vector< std::vector< Iterator > > &colored_iterators, Worker worker, Copier copier, const ScratchData &sample_scratch_data, const CopyData &sample_copy_data, const unsigned int queue_length=2 *MultithreadInfo::n_threads(), const unsigned int chunk_size=8)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
const ::Triangulation< dim, spacedim > & tria