Reference documentation for deal.II version 9.3.3
\(\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\}}\)
mapping_q_cache.cc
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2019 - 2021 by the deal.II authors
4//
5// This file is part of the deal.II library.
6//
7// The deal.II library is free software; you can use it, redistribute
8// it, and/or modify it under the terms of the GNU Lesser General
9// Public License as published by the Free Software Foundation; either
10// version 2.1 of the License, or (at your option) any later version.
11// The full text of the license can be found in the file LICENSE.md at
12// the top level directory of deal.II.
13//
14// ---------------------------------------------------------------------
15
19
21
22#include <deal.II/fe/fe_dgq.h>
24#include <deal.II/fe/fe_q.h>
25#include <deal.II/fe/fe_tools.h>
29
34
35#include <functional>
36
38
39template <int dim, int spacedim>
41 const unsigned int polynomial_degree)
42 : MappingQGeneric<dim, spacedim>(polynomial_degree)
43 , uses_level_info(false)
44{}
45
46
47
48template <int dim, int spacedim>
50 const MappingQCache<dim, spacedim> &mapping)
51 : MappingQGeneric<dim, spacedim>(mapping)
52 , support_point_cache(mapping.support_point_cache)
53 , uses_level_info(mapping.uses_level_info)
54{}
55
56
57
58template <int dim, int spacedim>
60{
61 // When this object goes out of scope, we want the cache to get cleared and
62 // free its memory before the signal is disconnected in order to not work on
63 // invalid memory that has been left back by freeing an object of this
64 // class.
65 support_point_cache.reset();
66 clear_signal.disconnect();
67}
68
69
70
71template <int dim, int spacedim>
72std::unique_ptr<Mapping<dim, spacedim>>
74{
75 return std::make_unique<MappingQCache<dim, spacedim>>(*this);
76}
77
78
79
80template <int dim, int spacedim>
81bool
83{
84 return false;
85}
86
87
88
89template <int dim, int spacedim>
90void
92 const Mapping<dim, spacedim> & mapping,
94{
95 // FE and FEValues in the case they are needed
98 fe_values_all;
99
100 this->initialize(
102 [&](const typename Triangulation<dim, spacedim>::cell_iterator &cell) {
103 const auto mapping_q_generic =
104 dynamic_cast<const MappingQGeneric<dim, spacedim> *>(&mapping);
105 if (mapping_q_generic != nullptr &&
106 this->get_degree() == mapping_q_generic->get_degree())
107 {
108 return mapping_q_generic->compute_mapping_support_points(cell);
109 }
110 else
111 {
112 // get FEValues (thread-safe); in the case that this thread has not
113 // created a an FEValues object yet, this helper-function also
114 // creates one with the right quadrature rule
115 auto &fe_values = fe_values_all.get();
116 if (fe_values.get() == nullptr)
117 {
118 QGaussLobatto<dim> quadrature_gl(this->polynomial_degree + 1);
119
120 std::vector<Point<dim>> quadrature_points;
121 for (const auto i :
122 FETools::hierarchic_to_lexicographic_numbering<dim>(
123 this->polynomial_degree))
124 quadrature_points.push_back(quadrature_gl.point(i));
126
127 fe_values = std::make_unique<FEValues<dim, spacedim>>(
128 mapping, fe, quadrature, update_quadrature_points);
129 }
130
131 fe_values->reinit(cell);
132 return fe_values->get_quadrature_points();
133 }
134 });
135}
136
137
138
139template <int dim, int spacedim>
140void
143 const MappingQGeneric<dim, spacedim> &mapping)
144{
145 this->initialize(mapping, triangulation);
146}
147
148
149
150template <int dim, int spacedim>
151void
154 const std::function<std::vector<Point<spacedim>>(
156 &compute_points_on_cell)
157{
158 clear_signal.disconnect();
159 clear_signal = triangulation.signals.any_change.connect(
160 [&]() -> void { this->support_point_cache.reset(); });
161
162 support_point_cache =
163 std::make_shared<std::vector<std::vector<std::vector<Point<spacedim>>>>>(
164 triangulation.n_levels());
165 for (unsigned int l = 0; l < triangulation.n_levels(); ++l)
166 (*support_point_cache)[l].resize(triangulation.n_raw_cells(l));
167
169 triangulation.begin(),
170 triangulation.end(),
171 [&](const typename Triangulation<dim, spacedim>::cell_iterator &cell,
172 void *,
173 void *) {
174 (*support_point_cache)[cell->level()][cell->index()] =
175 compute_points_on_cell(cell);
176 AssertDimension(
177 (*support_point_cache)[cell->level()][cell->index()].size(),
178 Utilities::pow(this->get_degree() + 1, dim));
179 },
180 /* copier */ std::function<void(void *)>(),
181 /* scratch_data */ nullptr,
182 /* copy_data */ nullptr,
184 /* chunk_size = */ 1);
185
186 uses_level_info = true;
187}
188
189
190
191template <int dim, int spacedim>
192void
194 const Mapping<dim, spacedim> & mapping,
196 const std::function<Point<spacedim>(
198 const Point<spacedim> &)> & transformation_function,
199 const bool function_describes_relative_displacement)
200{
201 // FE and FEValues in the case they are needed
204 fe_values_all;
205
206 this->initialize(
207 tria,
208 [&](const typename Triangulation<dim, spacedim>::cell_iterator &cell) {
209 std::vector<Point<spacedim>> points;
210
211 const auto mapping_q_generic =
212 dynamic_cast<const MappingQGeneric<dim, spacedim> *>(&mapping);
213
214 if (mapping_q_generic != nullptr &&
215 this->get_degree() == mapping_q_generic->get_degree())
216 {
217 points = mapping_q_generic->compute_mapping_support_points(cell);
218 }
219 else
220 {
221 // get FEValues (thread-safe); in the case that this thread has not
222 // created a an FEValues object yet, this helper-function also
223 // creates one with the right quadrature rule
224 auto &fe_values = fe_values_all.get();
225 if (fe_values.get() == nullptr)
226 {
227 QGaussLobatto<dim> quadrature_gl(this->polynomial_degree + 1);
228
229 std::vector<Point<dim>> quadrature_points;
230 for (const auto i :
231 FETools::hierarchic_to_lexicographic_numbering<dim>(
232 this->polynomial_degree))
233 quadrature_points.push_back(quadrature_gl.point(i));
235
236 fe_values = std::make_unique<FEValues<dim, spacedim>>(
237 mapping, fe, quadrature, update_quadrature_points);
238 }
239
240 fe_values->reinit(cell);
241 points = fe_values->get_quadrature_points();
242 }
243
244 for (auto &p : points)
245 if (function_describes_relative_displacement)
246 p += transformation_function(cell, p);
247 else
248 p = transformation_function(cell, p);
249
250 return points;
251 });
252
253 uses_level_info = true;
254}
255
256
257
258template <int dim, int spacedim>
259void
261 const Mapping<dim, spacedim> & mapping,
263 const Function<spacedim> & transformation_function,
264 const bool function_describes_relative_displacement)
265{
266 AssertDimension(transformation_function.n_components, spacedim);
267
268 this->initialize(mapping,
269 tria,
270 [&](const auto &, const auto &point) {
271 Point<spacedim> new_point;
272 for (int c = 0; c < spacedim; ++c)
273 new_point[c] = transformation_function.value(point, c);
274 return new_point;
275 },
276 function_describes_relative_displacement);
277
278 uses_level_info = true;
279}
280
281
282
283namespace
284{
285 template <typename VectorType>
286 void
287 copy_locally_owned_data_from(
288 const VectorType &vector,
290 &vector_ghosted)
291 {
293 temp.reinit(vector.locally_owned_elements());
294 temp.import(vector, VectorOperation::insert);
295 vector_ghosted.import(temp, VectorOperation::insert);
296 }
297} // namespace
298
299
300
301template <int dim, int spacedim>
302template <typename VectorType>
303void
305 const Mapping<dim, spacedim> & mapping,
306 const DoFHandler<dim, spacedim> &dof_handler,
307 const VectorType & vector,
308 const bool vector_describes_relative_displacement)
309{
310 AssertDimension(dof_handler.get_fe_collection().size(), 1);
311 const FiniteElement<dim, spacedim> &fe = dof_handler.get_fe();
313 AssertDimension(fe.element_multiplicity(0), spacedim);
314
315 const unsigned int is_fe_q =
316 dynamic_cast<const FE_Q<dim, spacedim> *>(&fe.base_element(0)) != nullptr;
317 const unsigned int is_fe_dgq =
318 dynamic_cast<const FE_DGQ<dim, spacedim> *>(&fe.base_element(0)) != nullptr;
319
322 FETools::hierarchic_to_lexicographic_numbering<spacedim>(
323 this->get_degree()));
324
325 // Step 1: copy global vector so that the ghost values are such that the
326 // cache can be set up for all ghost cells
328 vector_ghosted;
329 IndexSet locally_relevant_dofs;
330 DoFTools::extract_locally_relevant_dofs(dof_handler, locally_relevant_dofs);
331 vector_ghosted.reinit(dof_handler.locally_owned_dofs(),
332 locally_relevant_dofs,
333 dof_handler.get_communicator());
334 copy_locally_owned_data_from(vector, vector_ghosted);
335 vector_ghosted.update_ghost_values();
336
337 // FE and FEValues in the case they are needed
338 FE_Nothing<dim, spacedim> fe_nothing;
340 fe_values_all;
341
342 // Interpolation of values is needed if we cannot just read off locations
343 // from the solution vectors (as in the case of FE_Q and FE_DGQ with the
344 // same polynomial degree as this class has).
345 const bool interpolation_of_values_is_needed =
346 ((is_fe_q || is_fe_dgq) && fe.degree == this->get_degree()) == false;
347
348 // Step 2: loop over all cells
349 this->initialize(
350 dof_handler.get_triangulation(),
351 [&](const typename Triangulation<dim, spacedim>::cell_iterator &cell_tria)
352 -> std::vector<Point<spacedim>> {
353 const bool is_active_non_artificial_cell =
354 (cell_tria->is_active() == true) &&
355 (cell_tria->is_artificial() == false);
356
357 const typename DoFHandler<dim, spacedim>::cell_iterator cell_dofs(
358 &cell_tria->get_triangulation(),
359 cell_tria->level(),
360 cell_tria->index(),
361 &dof_handler);
362
363 const auto mapping_q_generic =
364 dynamic_cast<const MappingQGeneric<dim, spacedim> *>(&mapping);
365
366 // Step 2a) set up and reinit FEValues (if needed)
367 if (
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()) ==
372 false)) /*condition 1: points need to be computed via FEValues*/
373 ||
374 (is_active_non_artificial_cell && interpolation_of_values_is_needed) /*condition 2: interpolation of values is needed*/)
375 {
376 // get FEValues (thread-safe); in the case that this thread has
377 // not created a an FEValues object yet, this helper-function also
378 // creates one with the right quadrature rule
379 auto &fe_values = fe_values_all.get();
380 if (fe_values.get() == nullptr)
381 {
382 QGaussLobatto<dim> quadrature_gl(this->polynomial_degree + 1);
383
384 std::vector<Point<dim>> quadrature_points;
385 for (const auto i :
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);
390
391 fe_values = std::make_unique<FEValues<dim, spacedim>>(
392 mapping,
393 interpolation_of_values_is_needed ?
394 fe :
395 static_cast<const FiniteElement<dim, spacedim> &>(fe_nothing),
396 quadrature,
397 update_quadrature_points | update_values);
398 }
399
400 if (interpolation_of_values_is_needed)
401 fe_values->reinit(cell_dofs);
402 else
403 fe_values->reinit(cell_tria);
404 }
405
406 std::vector<Point<spacedim>> result;
407
408 // Step 2b) read of quadrature points in the relative displacement case
409 // note: we also take this path for non-active or artificial cells so that
410 // these cells are filled with some useful data
411 if (vector_describes_relative_displacement ||
412 is_active_non_artificial_cell == false)
413 {
414 if (mapping_q_generic != nullptr &&
415 this->get_degree() == mapping_q_generic->get_degree())
416 result =
417 mapping_q_generic->compute_mapping_support_points(cell_tria);
418 else
419 result = fe_values_all.get()->get_quadrature_points();
420
421 // for non-active or artificial cells we are done here and return
422 // the absolute positions, since the provided vector cannot contain
423 // any useful information for these cells
424 if (is_active_non_artificial_cell == false)
425 return result;
426 }
427 else
428 {
429 result.resize(
430 Utilities::pow<unsigned int>(this->get_degree() + 1, dim));
431 }
432
433 // Step 2c) read global vector and adjust points accordingly
434 if (interpolation_of_values_is_needed == false)
435 {
436 // case 1: FE_Q or FE_DGQ with same degree as this class has; this
437 // is the simple case since no interpolation is needed
438 std::vector<types::global_dof_index> dof_indices(
439 fe.n_dofs_per_cell());
440 cell_dofs->get_dof_indices(dof_indices);
441
442 for (unsigned int i = 0; i < dof_indices.size(); ++i)
443 {
444 const auto id = fe.system_to_component_index(i);
445
446 if (is_fe_q)
447 {
448 // case 1a: FE_Q
449 if (vector_describes_relative_displacement)
450 result[id.second][id.first] +=
451 vector_ghosted(dof_indices[i]);
452 else
453 result[id.second][id.first] =
454 vector_ghosted(dof_indices[i]);
455 }
456 else
457 {
458 // case 1b: FE_DGQ
459 if (vector_describes_relative_displacement)
461 [id.first] += vector_ghosted(dof_indices[i]);
462 else
464 [id.first] = vector_ghosted(dof_indices[i]);
465 }
466 }
467 }
468 else
469 {
470 // case 2: general case; interpolation is needed
471 // note: the following code could be optimized for tensor-product
472 // elements via application of sum factorization as is done on
473 // MatrixFree/FEEvaluation
474 auto &fe_values = fe_values_all.get();
475
476 std::vector<Vector<typename VectorType::value_type>> values(
477 fe_values->n_quadrature_points,
479
480 fe_values->get_function_values(vector_ghosted, values);
481
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];
486 else
487 result[q][c] = values[q][c];
488 }
489
490 return result;
491 });
492
493 uses_level_info = false;
494}
495
496
497
498template <int dim, int spacedim>
499template <typename VectorType>
500void
502 const Mapping<dim, spacedim> & mapping,
503 const DoFHandler<dim, spacedim> &dof_handler,
504 const MGLevelObject<VectorType> &vectors,
505 const bool vector_describes_relative_displacement)
506{
507 AssertDimension(dof_handler.get_fe_collection().size(), 1);
508 const FiniteElement<dim, spacedim> &fe = dof_handler.get_fe();
510 AssertDimension(fe.element_multiplicity(0), spacedim);
511 AssertDimension(0, vectors.min_level());
513 vectors.max_level());
514
515 const unsigned int is_fe_q =
516 dynamic_cast<const FE_Q<dim, spacedim> *>(&fe.base_element(0)) != nullptr;
517 const unsigned int is_fe_dgq =
518 dynamic_cast<const FE_DGQ<dim, spacedim> *>(&fe.base_element(0)) != nullptr;
519
522 FETools::hierarchic_to_lexicographic_numbering<spacedim>(
523 this->get_degree()));
524
525 // Step 1: copy global vector so that the ghost values are such that the
526 // cache can be set up for all ghost cells
529 vectors_ghosted(vectors.min_level(), vectors.max_level());
530
531 for (unsigned int l = vectors.min_level(); l <= vectors.max_level(); ++l)
532 {
533 IndexSet locally_relevant_dofs;
535 l,
536 locally_relevant_dofs);
537 vectors_ghosted[l].reinit(dof_handler.locally_owned_mg_dofs(l),
538 locally_relevant_dofs,
539 dof_handler.get_communicator());
540 copy_locally_owned_data_from(vectors[l], vectors_ghosted[l]);
541 vectors_ghosted[l].update_ghost_values();
542 }
543
544 // FE and FEValues in the case they are needed
545 FE_Nothing<dim, spacedim> fe_nothing;
547 fe_values_all;
548
549 // Interpolation of values is needed if we cannot just read off locations
550 // from the solution vectors (as in the case of FE_Q and FE_DGQ with the
551 // same polynomial degree as this class has).
552 const bool interpolation_of_values_is_needed =
553 ((is_fe_q || is_fe_dgq) && fe.degree == this->get_degree()) == false;
554
555 // Step 2: loop over all cells
556 this->initialize(
557 dof_handler.get_triangulation(),
558 [&](const typename Triangulation<dim, spacedim>::cell_iterator &cell_tria)
559 -> std::vector<Point<spacedim>> {
560 const bool is_non_artificial_cell =
561 cell_tria->level_subdomain_id() != numbers::artificial_subdomain_id;
562
563 const typename DoFHandler<dim, spacedim>::level_cell_iterator cell_dofs(
564 &cell_tria->get_triangulation(),
565 cell_tria->level(),
566 cell_tria->index(),
567 &dof_handler);
568
569 const auto mapping_q_generic =
570 dynamic_cast<const MappingQGeneric<dim, spacedim> *>(&mapping);
571
572 // Step 2a) set up and reinit FEValues (if needed)
573 if (
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()) ==
578 false)) /*condition 1: points need to be computed via FEValues*/
579 ||
580 (is_non_artificial_cell == true && interpolation_of_values_is_needed) /*condition 2: interpolation of values is needed*/)
581 {
582 // get FEValues (thread-safe); in the case that this thread has
583 // not created a an FEValues object yet, this helper-function also
584 // creates one with the right quadrature rule
585 auto &fe_values = fe_values_all.get();
586 if (fe_values.get() == nullptr)
587 {
588 QGaussLobatto<dim> quadrature_gl(this->polynomial_degree + 1);
589
590 std::vector<Point<dim>> quadrature_points;
591 for (const auto i :
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);
596
597 fe_values = std::make_unique<FEValues<dim, spacedim>>(
598 mapping,
599 interpolation_of_values_is_needed ?
600 fe :
601 static_cast<const FiniteElement<dim, spacedim> &>(fe_nothing),
602 quadrature,
603 update_quadrature_points | update_values);
604 }
605
606 if (interpolation_of_values_is_needed)
607 fe_values->reinit(cell_dofs);
608 else
609 fe_values->reinit(cell_tria);
610 }
611
612 std::vector<Point<spacedim>> result;
613
614 // Step 2b) read of quadrature points in the relative displacement case
615 // note: we also take this path for non-active or artificial cells so that
616 // these cells are filled with some useful data
617 if (vector_describes_relative_displacement ||
618 (is_non_artificial_cell == false))
619 {
620 if (mapping_q_generic != nullptr &&
621 this->get_degree() == mapping_q_generic->get_degree())
622 result =
623 mapping_q_generic->compute_mapping_support_points(cell_tria);
624 else
625 result = fe_values_all.get()->get_quadrature_points();
626
627 // for non-active or artificial cells we are done here and return
628 // the absolute positions, since the provided vector cannot contain
629 // any useful information for these cells
630 if (is_non_artificial_cell == false)
631 return result;
632 }
633 else
634 {
635 result.resize(
636 Utilities::pow<unsigned int>(this->get_degree() + 1, dim));
637 }
638
639 // Step 2c) read global vector and adjust points accordingly
640 if (interpolation_of_values_is_needed == false)
641 {
642 // case 1: FE_Q or FE_DGQ with same degree as this class has; this
643 // is the simple case since no interpolation is needed
644 std::vector<types::global_dof_index> dof_indices(
645 fe.n_dofs_per_cell());
646 cell_dofs->get_mg_dof_indices(dof_indices);
647
648 for (unsigned int i = 0; i < dof_indices.size(); ++i)
649 {
650 const auto id = fe.system_to_component_index(i);
651
652 if (is_fe_q)
653 {
654 // case 1a: FE_Q
655 if (vector_describes_relative_displacement)
656 result[id.second][id.first] +=
657 vectors_ghosted[cell_tria->level()](dof_indices[i]);
658 else
659 result[id.second][id.first] =
660 vectors_ghosted[cell_tria->level()](dof_indices[i]);
661 }
662 else
663 {
664 // case 1b: FE_DGQ
665 if (vector_describes_relative_displacement)
667 [id.first] +=
668 vectors_ghosted[cell_tria->level()](dof_indices[i]);
669 else
671 [id.first] =
672 vectors_ghosted[cell_tria->level()](dof_indices[i]);
673 }
674 }
675 }
676 else
677 {
678 // case 2: general case; interpolation is needed
679 // note: the following code could be optimized for tensor-product
680 // elements via application of sum factorization as is done on
681 // MatrixFree/FEEvaluation
682 auto &fe_values = fe_values_all.get();
683
684 std::vector<types::global_dof_index> dof_indices(
685 fe.n_dofs_per_cell());
686 cell_dofs->get_mg_dof_indices(dof_indices);
687
688 std::vector<typename VectorType::value_type> dof_values(
689 fe.n_dofs_per_cell());
690
691 for (unsigned int i = 0; i < fe.n_dofs_per_cell(); ++i)
692 dof_values[i] = vectors_ghosted[cell_tria->level()](dof_indices[i]);
693
694 for (unsigned int c = 0; c < spacedim; ++c)
695 for (unsigned int i = 0; i < fe.n_dofs_per_cell(); ++i)
696 for (unsigned int q = 0; q < fe_values->n_quadrature_points; ++q)
697 if (vector_describes_relative_displacement == false && i == 0)
698 result[q][c] =
699 dof_values[i] * fe_values->shape_value_component(i, q, c);
700 else
701 result[q][c] +=
702 dof_values[i] * fe_values->shape_value_component(i, q, c);
703 }
704
705 return result;
706 });
707
708 uses_level_info = true;
709}
710
711
712
713template <int dim, int spacedim>
714std::size_t
716{
717 if (support_point_cache.get() != nullptr)
718 return sizeof(*this) +
719 MemoryConsumption::memory_consumption(*support_point_cache);
720 else
721 return sizeof(*this);
722}
723
724
725
726template <int dim, int spacedim>
727std::vector<Point<spacedim>>
729 const typename Triangulation<dim, spacedim>::cell_iterator &cell) const
730{
731 Assert(support_point_cache.get() != nullptr,
732 ExcMessage("Must call MappingQCache::initialize() before "
733 "using it or after mesh has changed!"));
734
735 Assert(uses_level_info || cell->is_active(), ExcInternalError());
736
737 AssertIndexRange(cell->level(), support_point_cache->size());
738 AssertIndexRange(cell->index(), (*support_point_cache)[cell->level()].size());
739 return (*support_point_cache)[cell->level()][cell->index()];
740}
741
742
743
744//--------------------------- Explicit instantiations -----------------------
745#include "mapping_q_cache.inst"
746
747
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
Definition: fe_dgq.h:111
Definition: fe_q.h:549
const unsigned int degree
Definition: fe_base.h:435
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
Definition: function.h:164
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)
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
Definition: collection.h:109
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
@ update_quadrature_points
Transformed quadrature points.
#define Assert(cond, exc)
Definition: exceptions.h:1465
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1622
#define AssertIndexRange(index, range)
Definition: exceptions.h:1690
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={})
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)
void extract_locally_relevant_level_dofs(const DoFHandler< dim, spacedim > &dof_handler, const unsigned int level, IndexSet &dof_set)
Definition: dof_tools.cc:1252
void extract_locally_relevant_dofs(const DoFHandler< dim, spacedim > &dof_handler, IndexSet &dof_set)
Definition: dof_tools.cc:1210
std::vector< unsigned int > lexicographic_to_hierarchic_numbering(unsigned int degree)
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)
Definition: utilities.cc:188
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={})
Definition: generators.cc:451
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)
Definition: utilities.h:1482
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)
Definition: work_stream.h:1252
unsigned int get_degree(const std::vector< BarycentricPolynomial< dim > > &polys)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation