39template <
int dim,
int spacedim>
41 const unsigned int 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_generic =
105 if (mapping_q_generic !=
nullptr &&
106 this->
get_degree() == mapping_q_generic->get_degree())
115 auto &fe_values = fe_values_all.
get();
116 if (fe_values.get() ==
nullptr)
122 FETools::hierarchic_to_lexicographic_numbering<dim>(
123 this->polynomial_degree))
127 fe_values = std::make_unique<FEValues<dim, spacedim>>(
131 fe_values->reinit(cell);
132 return fe_values->get_quadrature_points();
139template <
int dim,
int spacedim>
150template <
int dim,
int spacedim>
156 &compute_points_on_cell)
158 clear_signal.disconnect();
160 [&]() ->
void { this->support_point_cache.reset(); });
162 support_point_cache =
163 std::make_shared<std::vector<std::vector<std::vector<Point<spacedim>>>>>(
174 (*support_point_cache)[cell->level()][cell->index()] =
175 compute_points_on_cell(cell);
177 (*support_point_cache)[cell->level()][cell->index()].size(),
178 Utilities::pow(this->get_degree() + 1, dim));
180 std::function<
void(
void *)>(),
186 uses_level_info =
true;
191template <
int dim,
int spacedim>
199 const bool function_describes_relative_displacement)
209 std::vector<Point<spacedim>> points;
211 const auto mapping_q_generic =
214 if (mapping_q_generic !=
nullptr &&
215 this->
get_degree() == mapping_q_generic->get_degree())
224 auto &fe_values = fe_values_all.
get();
225 if (fe_values.get() ==
nullptr)
231 FETools::hierarchic_to_lexicographic_numbering<dim>(
232 this->polynomial_degree))
236 fe_values = std::make_unique<FEValues<dim, spacedim>>(
240 fe_values->reinit(cell);
241 points = fe_values->get_quadrature_points();
244 for (
auto &p : points)
245 if (function_describes_relative_displacement)
246 p += transformation_function(cell, p);
248 p = transformation_function(cell, p);
253 uses_level_info =
true;
258template <
int dim,
int spacedim>
264 const bool function_describes_relative_displacement)
268 this->initialize(mapping,
270 [&](
const auto &,
const auto &
point) {
272 for (
int c = 0; c < spacedim; ++c)
273 new_point[c] = transformation_function.
value(
point, c);
276 function_describes_relative_displacement);
278 uses_level_info =
true;
285 template <
typename VectorType>
287 copy_locally_owned_data_from(
288 const VectorType &vector,
293 temp.
reinit(vector.locally_owned_elements());
301template <
int dim,
int spacedim>
302template <
typename VectorType>
307 const VectorType & vector,
308 const bool vector_describes_relative_displacement)
315 const unsigned int is_fe_q =
317 const unsigned int is_fe_dgq =
322 FETools::hierarchic_to_lexicographic_numbering<spacedim>(
332 locally_relevant_dofs,
334 copy_locally_owned_data_from(vector, vector_ghosted);
345 const bool interpolation_of_values_is_needed =
346 ((is_fe_q || is_fe_dgq) && fe.
degree == this->get_degree()) ==
false;
353 const bool is_active_non_artificial_cell =
354 (cell_tria->is_active() == true) &&
355 (cell_tria->is_artificial() == false);
357 const typename DoFHandler<dim, spacedim>::cell_iterator cell_dofs(
358 &cell_tria->get_triangulation(),
363 const auto mapping_q_generic =
364 dynamic_cast<const MappingQGeneric<dim, spacedim> *>(&mapping);
368 ((vector_describes_relative_displacement ||
369 (is_active_non_artificial_cell == false)) &&
370 ((mapping_q_generic != nullptr &&
371 this->get_degree() == mapping_q_generic->get_degree()) ==
374 (is_active_non_artificial_cell && interpolation_of_values_is_needed) )
379 auto &fe_values = fe_values_all.get();
380 if (fe_values.get() == nullptr)
382 QGaussLobatto<dim> quadrature_gl(this->polynomial_degree + 1);
384 std::vector<Point<dim>> quadrature_points;
386 FETools::hierarchic_to_lexicographic_numbering<dim>(
387 this->polynomial_degree))
388 quadrature_points.push_back(quadrature_gl.point(i));
389 Quadrature<dim> quadrature(quadrature_points);
391 fe_values = std::make_unique<FEValues<dim, spacedim>>(
393 interpolation_of_values_is_needed ?
395 static_cast<const FiniteElement<dim, spacedim> &>(fe_nothing),
397 update_quadrature_points | update_values);
400 if (interpolation_of_values_is_needed)
401 fe_values->reinit(cell_dofs);
403 fe_values->reinit(cell_tria);
411 if (vector_describes_relative_displacement ||
412 is_active_non_artificial_cell ==
false)
414 if (mapping_q_generic !=
nullptr &&
415 this->
get_degree() == mapping_q_generic->get_degree())
417 mapping_q_generic->compute_mapping_support_points(cell_tria);
419 result = fe_values_all.
get()->get_quadrature_points();
424 if (is_active_non_artificial_cell ==
false)
430 Utilities::pow<unsigned int>(this->
get_degree() + 1, dim));
434 if (interpolation_of_values_is_needed ==
false)
438 std::vector<types::global_dof_index> dof_indices(
440 cell_dofs->get_dof_indices(dof_indices);
442 for (
unsigned int i = 0; i < dof_indices.size(); ++i)
449 if (vector_describes_relative_displacement)
450 result[
id.second][
id.first] +=
451 vector_ghosted(dof_indices[i]);
453 result[
id.second][
id.first] =
454 vector_ghosted(dof_indices[i]);
459 if (vector_describes_relative_displacement)
461 [
id.first] += vector_ghosted(dof_indices[i]);
464 [
id.first] = vector_ghosted(dof_indices[i]);
474 auto &fe_values = fe_values_all.
get();
476 std::vector<Vector<typename VectorType::value_type>>
values(
477 fe_values->n_quadrature_points,
480 fe_values->get_function_values(vector_ghosted,
values);
482 for (
unsigned int q = 0; q < fe_values->n_quadrature_points; ++q)
483 for (
unsigned int c = 0; c < spacedim; ++c)
484 if (vector_describes_relative_displacement)
485 result[q][c] +=
values[q][c];
487 result[q][c] =
values[q][c];
493 uses_level_info =
false;
498template <
int dim,
int spacedim>
499template <
typename VectorType>
505 const bool vector_describes_relative_displacement)
515 const unsigned int is_fe_q =
517 const unsigned int is_fe_dgq =
522 FETools::hierarchic_to_lexicographic_numbering<spacedim>(
536 locally_relevant_dofs);
538 locally_relevant_dofs,
540 copy_locally_owned_data_from(vectors[
l], vectors_ghosted[
l]);
541 vectors_ghosted[
l].update_ghost_values();
552 const bool interpolation_of_values_is_needed =
553 ((is_fe_q || is_fe_dgq) && fe.
degree == this->get_degree()) ==
false;
560 const bool is_non_artificial_cell =
561 cell_tria->level_subdomain_id() != numbers::artificial_subdomain_id;
563 const typename DoFHandler<dim, spacedim>::level_cell_iterator cell_dofs(
564 &cell_tria->get_triangulation(),
569 const auto mapping_q_generic =
570 dynamic_cast<const MappingQGeneric<dim, spacedim> *>(&mapping);
574 ((vector_describes_relative_displacement ||
575 (is_non_artificial_cell == false)) &&
576 ((mapping_q_generic != nullptr &&
577 this->get_degree() == mapping_q_generic->get_degree()) ==
580 (is_non_artificial_cell == true && interpolation_of_values_is_needed) )
585 auto &fe_values = fe_values_all.get();
586 if (fe_values.get() == nullptr)
588 QGaussLobatto<dim> quadrature_gl(this->polynomial_degree + 1);
590 std::vector<Point<dim>> quadrature_points;
592 FETools::hierarchic_to_lexicographic_numbering<dim>(
593 this->polynomial_degree))
594 quadrature_points.push_back(quadrature_gl.point(i));
595 Quadrature<dim> quadrature(quadrature_points);
597 fe_values = std::make_unique<FEValues<dim, spacedim>>(
599 interpolation_of_values_is_needed ?
601 static_cast<const FiniteElement<dim, spacedim> &>(fe_nothing),
603 update_quadrature_points | update_values);
606 if (interpolation_of_values_is_needed)
607 fe_values->reinit(cell_dofs);
609 fe_values->reinit(cell_tria);
617 if (vector_describes_relative_displacement ||
618 (is_non_artificial_cell ==
false))
620 if (mapping_q_generic !=
nullptr &&
621 this->
get_degree() == mapping_q_generic->get_degree())
623 mapping_q_generic->compute_mapping_support_points(cell_tria);
625 result = fe_values_all.
get()->get_quadrature_points();
630 if (is_non_artificial_cell ==
false)
636 Utilities::pow<unsigned int>(this->
get_degree() + 1, dim));
640 if (interpolation_of_values_is_needed ==
false)
644 std::vector<types::global_dof_index> dof_indices(
646 cell_dofs->get_mg_dof_indices(dof_indices);
648 for (
unsigned int i = 0; i < dof_indices.size(); ++i)
655 if (vector_describes_relative_displacement)
656 result[
id.second][
id.first] +=
657 vectors_ghosted[cell_tria->level()](dof_indices[i]);
659 result[
id.second][
id.first] =
660 vectors_ghosted[cell_tria->level()](dof_indices[i]);
665 if (vector_describes_relative_displacement)
668 vectors_ghosted[cell_tria->level()](dof_indices[i]);
672 vectors_ghosted[cell_tria->level()](dof_indices[i]);
682 auto &fe_values = fe_values_all.
get();
684 std::vector<types::global_dof_index> dof_indices(
686 cell_dofs->get_mg_dof_indices(dof_indices);
688 std::vector<typename VectorType::value_type> dof_values(
692 dof_values[i] = vectors_ghosted[cell_tria->level()](dof_indices[i]);
694 for (
unsigned int c = 0; c < spacedim; ++c)
696 for (
unsigned int q = 0; q < fe_values->n_quadrature_points; ++q)
697 if (vector_describes_relative_displacement ==
false && i == 0)
699 dof_values[i] * fe_values->shape_value_component(i, q, c);
702 dof_values[i] * fe_values->shape_value_component(i, q, c);
708 uses_level_info =
true;
713template <
int dim,
int spacedim>
717 if (support_point_cache.get() !=
nullptr)
718 return sizeof(*this) +
721 return sizeof(*this);
726template <
int dim,
int spacedim>
727std::vector<Point<spacedim>>
731 Assert(support_point_cache.get() !=
nullptr,
732 ExcMessage(
"Must call MappingQCache::initialize() before "
733 "using it or after mesh has changed!"));
738 AssertIndexRange(cell->index(), (*support_point_cache)[cell->level()].size());
739 return (*support_point_cache)[cell->level()][cell->index()];
745#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 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)
friend class MappingQCache
virtual std::vector< Point< spacedim > > compute_mapping_support_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
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)
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
void quadrature_points(const Triangulation< dim, spacedim > &triangulation, const Quadrature< dim > &quadrature, const std::vector< std::vector< BoundingBox< spacedim > > > &global_bounding_boxes, ParticleHandler< dim, spacedim > &particle_handler, const Mapping< dim, spacedim > &mapping=(ReferenceCells::get_hypercube< dim >() .template get_default_linear_mapping< dim, spacedim >()), const std::vector< std::vector< double > > &properties={})
Tensor< 2, dim, Number > l(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
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)
unsigned int get_degree(const std::vector< BarycentricPolynomial< dim > > &polys)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation