38template <
int dim,
int spacedim>
40 const unsigned int polynomial_degree)
41 :
MappingQ<dim, spacedim>(polynomial_degree)
42 , uses_level_info(false)
47template <
int dim,
int spacedim>
51 , support_point_cache(mapping.support_point_cache)
52 , uses_level_info(mapping.uses_level_info)
57template <
int dim,
int spacedim>
64 support_point_cache.reset();
65 clear_signal.disconnect();
70template <
int dim,
int spacedim>
71std::unique_ptr<Mapping<dim, spacedim>>
74 return std::make_unique<MappingQCache<dim, spacedim>>(*this);
79template <
int dim,
int spacedim>
88template <
int dim,
int spacedim>
102 const auto mapping_q =
104 if (mapping_q !=
nullptr && this->get_degree() == mapping_q->get_degree())
113 auto &fe_values = fe_values_all.
get();
114 if (fe_values.get() ==
nullptr)
118 std::vector<Point<dim>> quadrature_points;
120 FETools::hierarchic_to_lexicographic_numbering<dim>(
121 this->polynomial_degree))
122 quadrature_points.push_back(quadrature_gl.
point(i));
125 fe_values = std::make_unique<FEValues<dim, spacedim>>(
129 fe_values->reinit(cell);
130 return fe_values->get_quadrature_points();
137template <
int dim,
int spacedim>
148template <
int dim,
int spacedim>
154 &compute_points_on_cell)
156 clear_signal.disconnect();
158 [&]() ->
void { this->support_point_cache.reset(); });
160 support_point_cache =
161 std::make_shared<std::vector<std::vector<std::vector<Point<spacedim>>>>>(
164 (*support_point_cache)[l].resize(
triangulation.n_raw_cells(l));
172 (*support_point_cache)[cell->level()][cell->index()] =
173 compute_points_on_cell(cell);
175 (*support_point_cache)[cell->level()][cell->index()].size(),
176 Utilities::pow(this->get_degree() + 1, dim));
178 std::function<
void(
void *)>(),
184 uses_level_info =
true;
189template <
int dim,
int spacedim>
197 const bool function_describes_relative_displacement)
207 std::vector<Point<spacedim>> points;
209 const auto mapping_q =
212 if (mapping_q !=
nullptr && this->get_degree() == mapping_q->get_degree())
221 auto &fe_values = fe_values_all.
get();
222 if (fe_values.get() ==
nullptr)
226 std::vector<Point<dim>> quadrature_points;
228 FETools::hierarchic_to_lexicographic_numbering<dim>(
229 this->polynomial_degree))
230 quadrature_points.push_back(quadrature_gl.
point(i));
233 fe_values = std::make_unique<FEValues<dim, spacedim>>(
237 fe_values->reinit(cell);
238 points = fe_values->get_quadrature_points();
241 for (
auto &p : points)
242 if (function_describes_relative_displacement)
243 p += transformation_function(cell, p);
245 p = transformation_function(cell, p);
250 uses_level_info =
true;
255template <
int dim,
int spacedim>
261 const bool function_describes_relative_displacement)
268 [&](
const auto &,
const auto &point) {
270 for (
unsigned int c = 0; c < spacedim; ++c)
271 new_point[c] = transformation_function.
value(point, c);
274 function_describes_relative_displacement);
276 uses_level_info =
true;
283 template <
typename VectorType>
285 copy_locally_owned_data_from(
286 const VectorType &vector,
291 temp.
reinit(vector.locally_owned_elements());
299template <
int dim,
int spacedim>
300template <
typename VectorType>
305 const VectorType & vector,
306 const bool vector_describes_relative_displacement)
313 const unsigned int is_fe_q =
315 const unsigned int is_fe_dgq =
318 const auto lexicographic_to_hierarchic_numbering =
320 FETools::hierarchic_to_lexicographic_numbering<spacedim>(
321 this->get_degree()));
330 locally_relevant_dofs,
332 copy_locally_owned_data_from(vector, vector_ghosted);
343 const bool interpolation_of_values_is_needed =
344 ((is_fe_q || is_fe_dgq) && fe.
degree == this->get_degree()) ==
false;
351 const bool is_active_non_artificial_cell =
352 (cell_tria->is_active() == true) &&
353 (cell_tria->is_artificial() == false);
355 const typename DoFHandler<dim, spacedim>::cell_iterator cell_dofs(
356 &cell_tria->get_triangulation(),
361 const auto mapping_q =
362 dynamic_cast<const MappingQ<dim, spacedim> *>(&mapping);
366 ((vector_describes_relative_displacement ||
367 (is_active_non_artificial_cell == false)) &&
368 ((mapping_q != nullptr &&
369 this->get_degree() == mapping_q->get_degree()) ==
372 (is_active_non_artificial_cell && interpolation_of_values_is_needed) )
377 auto &fe_values = fe_values_all.get();
378 if (fe_values.get() == nullptr)
380 QGaussLobatto<dim> quadrature_gl(this->polynomial_degree + 1);
382 std::vector<Point<dim>> quadrature_points;
384 FETools::hierarchic_to_lexicographic_numbering<dim>(
385 this->polynomial_degree))
386 quadrature_points.push_back(quadrature_gl.point(i));
387 Quadrature<dim> quadrature(quadrature_points);
389 fe_values = std::make_unique<FEValues<dim, spacedim>>(
391 interpolation_of_values_is_needed ?
393 static_cast<const FiniteElement<dim, spacedim> &>(fe_nothing),
395 update_quadrature_points | update_values);
398 if (interpolation_of_values_is_needed)
399 fe_values->reinit(cell_dofs);
401 fe_values->reinit(cell_tria);
409 if (vector_describes_relative_displacement ||
410 is_active_non_artificial_cell ==
false)
412 if (mapping_q !=
nullptr &&
413 this->get_degree() == mapping_q->get_degree())
414 result = mapping_q->compute_mapping_support_points(cell_tria);
416 result = fe_values_all.
get()->get_quadrature_points();
421 if (is_active_non_artificial_cell ==
false)
427 Utilities::pow<unsigned int>(this->get_degree() + 1, dim));
431 if (interpolation_of_values_is_needed ==
false)
435 std::vector<types::global_dof_index> dof_indices(
437 cell_dofs->get_dof_indices(dof_indices);
439 for (
unsigned int i = 0; i < dof_indices.size(); ++i)
446 if (vector_describes_relative_displacement)
447 result[
id.second][
id.first] +=
448 vector_ghosted(dof_indices[i]);
450 result[
id.second][
id.first] =
451 vector_ghosted(dof_indices[i]);
456 if (vector_describes_relative_displacement)
457 result[lexicographic_to_hierarchic_numbering[
id.second]]
458 [
id.first] += vector_ghosted(dof_indices[i]);
460 result[lexicographic_to_hierarchic_numbering[
id.second]]
461 [
id.first] = vector_ghosted(dof_indices[i]);
471 auto &fe_values = fe_values_all.
get();
473 std::vector<Vector<typename VectorType::value_type>> values(
474 fe_values->n_quadrature_points,
477 fe_values->get_function_values(vector_ghosted, values);
479 for (
unsigned int q = 0; q < fe_values->n_quadrature_points; ++q)
480 for (
unsigned int c = 0; c < spacedim; ++c)
481 if (vector_describes_relative_displacement)
482 result[q][c] += values[q][c];
484 result[q][c] = values[q][c];
490 uses_level_info =
false;
495template <
int dim,
int spacedim>
496template <
typename VectorType>
502 const bool vector_describes_relative_displacement)
512 const unsigned int is_fe_q =
514 const unsigned int is_fe_dgq =
517 const auto lexicographic_to_hierarchic_numbering =
519 FETools::hierarchic_to_lexicographic_numbering<spacedim>(
520 this->get_degree()));
533 locally_relevant_dofs);
535 locally_relevant_dofs,
537 copy_locally_owned_data_from(vectors[l], vectors_ghosted[l]);
538 vectors_ghosted[l].update_ghost_values();
549 const bool interpolation_of_values_is_needed =
550 ((is_fe_q || is_fe_dgq) && fe.
degree == this->get_degree()) ==
false;
557 const bool is_non_artificial_cell =
558 cell_tria->level_subdomain_id() != numbers::artificial_subdomain_id;
560 const typename DoFHandler<dim, spacedim>::level_cell_iterator cell_dofs(
561 &cell_tria->get_triangulation(),
566 const auto mapping_q =
567 dynamic_cast<const MappingQ<dim, spacedim> *>(&mapping);
571 ((vector_describes_relative_displacement ||
572 (is_non_artificial_cell == false)) &&
573 ((mapping_q != nullptr &&
574 this->get_degree() == mapping_q->get_degree()) ==
577 (is_non_artificial_cell == true && interpolation_of_values_is_needed) )
582 auto &fe_values = fe_values_all.get();
583 if (fe_values.get() == nullptr)
585 QGaussLobatto<dim> quadrature_gl(this->polynomial_degree + 1);
587 std::vector<Point<dim>> quadrature_points;
589 FETools::hierarchic_to_lexicographic_numbering<dim>(
590 this->polynomial_degree))
591 quadrature_points.push_back(quadrature_gl.point(i));
592 Quadrature<dim> quadrature(quadrature_points);
594 fe_values = std::make_unique<FEValues<dim, spacedim>>(
596 interpolation_of_values_is_needed ?
598 static_cast<const FiniteElement<dim, spacedim> &>(fe_nothing),
600 update_quadrature_points | update_values);
603 if (interpolation_of_values_is_needed)
604 fe_values->reinit(cell_dofs);
606 fe_values->reinit(cell_tria);
614 if (vector_describes_relative_displacement ||
615 (is_non_artificial_cell ==
false))
617 if (mapping_q !=
nullptr &&
618 this->get_degree() == mapping_q->get_degree())
619 result = mapping_q->compute_mapping_support_points(cell_tria);
621 result = fe_values_all.
get()->get_quadrature_points();
626 if (is_non_artificial_cell ==
false)
632 Utilities::pow<unsigned int>(this->get_degree() + 1, dim));
636 if (interpolation_of_values_is_needed ==
false)
640 std::vector<types::global_dof_index> dof_indices(
642 cell_dofs->get_mg_dof_indices(dof_indices);
644 for (
unsigned int i = 0; i < dof_indices.size(); ++i)
651 if (vector_describes_relative_displacement)
652 result[
id.second][
id.first] +=
653 vectors_ghosted[cell_tria->level()](dof_indices[i]);
655 result[
id.second][
id.first] =
656 vectors_ghosted[cell_tria->level()](dof_indices[i]);
661 if (vector_describes_relative_displacement)
662 result[lexicographic_to_hierarchic_numbering[
id.second]]
664 vectors_ghosted[cell_tria->level()](dof_indices[i]);
666 result[lexicographic_to_hierarchic_numbering[
id.second]]
668 vectors_ghosted[cell_tria->level()](dof_indices[i]);
678 auto &fe_values = fe_values_all.
get();
680 std::vector<types::global_dof_index> dof_indices(
682 cell_dofs->get_mg_dof_indices(dof_indices);
684 std::vector<typename VectorType::value_type> dof_values(
688 dof_values[i] = vectors_ghosted[cell_tria->level()](dof_indices[i]);
690 for (
unsigned int c = 0; c < spacedim; ++c)
692 for (
unsigned int q = 0; q < fe_values->n_quadrature_points; ++q)
693 if (vector_describes_relative_displacement ==
false && i == 0)
695 dof_values[i] * fe_values->shape_value_component(i, q, c);
698 dof_values[i] * fe_values->shape_value_component(i, q, c);
704 uses_level_info =
true;
709template <
int dim,
int spacedim>
713 if (support_point_cache.get() !=
nullptr)
714 return sizeof(*this) +
717 return sizeof(*this);
722template <
int dim,
int spacedim>
723std::vector<Point<spacedim>>
727 Assert(support_point_cache.get() !=
nullptr,
728 ExcMessage(
"Must call MappingQCache::initialize() before "
729 "using it or after mesh has changed!"));
734 AssertIndexRange(cell->index(), (*support_point_cache)[cell->level()].size());
735 return (*support_point_cache)[cell->level()][cell->index()];
740template <
int dim,
int spacedim>
741boost::container::small_vector<Point<spacedim>,
746 Assert(support_point_cache.get() !=
nullptr,
747 ExcMessage(
"Must call MappingQCache::initialize() before "
748 "using it or after mesh has changed!"));
753 AssertIndexRange(cell->index(), (*support_point_cache)[cell->level()].size());
754 const auto ptr = (*support_point_cache)[cell->level()][cell->index()].begin();
755 return boost::container::small_vector<Point<spacedim>,
757 ptr, ptr + cell->n_vertices());
763#include "mapping_q_cache.inst"
const hp::FECollection< dim, spacedim > & get_fe_collection() const
const FiniteElement< dim, spacedim > & get_fe(const types::fe_index index=0) 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
void import_elements(const ::Vector< Number > &vec, VectorOperation::values operation, const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &communication_pattern={})
virtual void reinit(const size_type size, const bool omit_zeroing_entries=false)
void update_ghost_values() const
void import_elements(const Vector< Number, MemorySpace2 > &src, VectorOperation::values operation)
void reinit(const size_type size, const bool omit_zeroing_entries=false)
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
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
@ update_quadrature_points
Transformed quadrature points.
std::enable_if_t< std::is_fundamental< T >::value, std::size_t > 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