deal.II version 9.7.0
\(\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\}}\)
Loading...
Searching...
No Matches
fe_poly_tensor.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2005 - 2025 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15
25
29
31#include <deal.II/grid/tria.h>
33
35
36namespace internal
37{
38 namespace FE_PolyTensor
39 {
40 namespace
41 {
42 template <int spacedim>
43 void
44 get_dof_sign_change_h_div(
45 const typename ::Triangulation<1, spacedim>::cell_iterator &,
47 const std::vector<MappingKind> &,
48 std::vector<double> &)
49 {
50 // Nothing to do in 1d.
51 }
52
53
54
55 // TODO: This function is not a consistent fix of the orientation issue
56 // like in 3d. It is rather kept not to break legacy behavior in 2d but
57 // should be replaced. See also the implementation of
58 // FE_RaviartThomas<dim>::initialize_quad_dof_index_permutation_and_sign_change()
59 // or other H(div) conforming elements such as FE_ABF<dim> and
60 // FE_BDM<dim>.
61 template <int spacedim>
62 void
63 get_dof_sign_change_h_div(
64 const typename ::Triangulation<2, spacedim>::cell_iterator &cell,
66 const std::vector<MappingKind> &mapping_kind,
67 std::vector<double> &face_sign)
68 {
69 const unsigned int dim = 2;
70 // const unsigned int spacedim = 2;
71
72 const CellId this_cell_id = cell->id();
73
74 for (unsigned int f = 0; f < GeometryInfo<dim>::faces_per_cell; ++f)
75 {
76 typename ::Triangulation<dim, spacedim>::face_iterator face =
77 cell->face(f);
78 if (!face->at_boundary())
79 {
80 const unsigned int nn = cell->neighbor_face_no(f);
81 const typename ::Triangulation<dim,
82 spacedim>::cell_iterator
83 neighbor_cell_at_face = cell->neighbor(f);
84 const CellId neighbor_cell_id = neighbor_cell_at_face->id();
85
86 // Only fix sign if the orientation is opposite and only do so
87 // on the face dofs on the cell with smaller cell_id.
88 if (((nn + f) % 2 == 0) && this_cell_id < neighbor_cell_id)
89 for (unsigned int j = 0; j < fe.n_dofs_per_face(f); ++j)
90 {
91 const unsigned int cell_j = fe.face_to_cell_index(j, f);
92
93 Assert(f * fe.n_dofs_per_face(f) + j < face_sign.size(),
95 Assert(mapping_kind.size() == 1 ||
96 cell_j < mapping_kind.size(),
98
99 // TODO: This is probably only going to work for those
100 // elements for which all dofs are face dofs
101 if ((mapping_kind.size() > 1 ?
102 mapping_kind[cell_j] :
103 mapping_kind[0]) == mapping_raviart_thomas)
104 face_sign[f * fe.n_dofs_per_face(f) + j] = -1.0;
105 }
106 }
107 }
108 }
109
110
111
112 template <int spacedim>
113 void
114 get_dof_sign_change_h_div(
115 const typename ::Triangulation<3, spacedim>::cell_iterator
116 & /*cell*/,
117 const FiniteElement<3, spacedim> & /*fe*/,
118 const std::vector<MappingKind> & /*mapping_kind*/,
119 std::vector<double> & /*face_sign*/)
120 {
121 // Nothing to do. In 3d we take care of it through the
122 // adjust_quad_dof_sign_for_face_orientation_table
123 }
124
125 template <int spacedim>
126 void
127 get_dof_sign_change_nedelec(
128 const typename ::Triangulation<1, spacedim>::cell_iterator
129 & /*cell*/,
130 const FiniteElement<1, spacedim> & /*fe*/,
131 const std::vector<MappingKind> & /*mapping_kind*/,
132 std::vector<double> & /*line_dof_sign*/)
133 {
134 // nothing to do in 1d
135 }
136
137 template <int spacedim>
138 void
139 get_dof_sign_change_nedelec(
140 const typename ::Triangulation<2, spacedim>::cell_iterator &cell,
142 const std::vector<MappingKind> &mapping_kind,
143 std::vector<double> &line_dof_sign)
144 {
145 // The Nedelec finite elements in two spacial dimensions have two types
146 // of dofs: the line dofs and the quad dofs. The line dofs are
147 // associated with the edges. They are shared between the neighbouring
148 // cells and, consequently, need to be adjusted to compensate for a
149 // possible mismatch in the edge orientation. The quad dofs, they are
150 // associated with the interiors of the cells, are not shared between
151 // the cells in two spacial dimensions and need no adjustments.
152 //
153 // The Nedelec finite elements in two spacial dimensions have
154 // 2*(k+1)*(k+2) dofs per cell. All dofs are distributed between the
155 // line and quad dofs as the following:
156 //
157 // 4*(k+1) line dofs; (k+1) dofs per line.
158 // 2*k*(k+1) quad dofs.
159 //
160 // The dofs are indexed in the following order: first all line dofs,
161 // then all quad dofs.
162 //
163 // Here we adjust the line dofs. The sign of a line dof needs to be
164 // changed if the edge on which the dof resides points in the opposite
165 // direction.
166 const unsigned int k = fe.tensor_degree() - 1;
167
168 for (unsigned int l = 0; l < GeometryInfo<2>::lines_per_cell; ++l)
169 if (cell->line_orientation(l) !=
171 mapping_kind[0] == mapping_nedelec)
172 {
173 if (k == 0)
174 {
175 // The lowest order element (k=0) is straightforward, because
176 // there is a single dof per edge, which needs to be flipped:
177 line_dof_sign[l] = -1.0;
178 }
179 else
180 {
181 // The case k > 0 is a bit more complicated. As we adjust
182 // only edge dofs in this function, we need to concern
183 // ourselves with the first 4*(k+1) entries in line_dof_sign
184 // vector ignoring the rest. There are (k+1) dofs per edge.
185 // Let us consider the local dof indices on one edge,
186 // local_line_dof = 0...k. The shape functions with even
187 // indices are asymmetric. The corresponding dofs need sign
188 // adjustment if the edge points in the opposite direction.
189 // The shape functions with odd indices are symmetric. The
190 // corresponding dofs need no sign adjustment even if the edge
191 // points in the opposite direction. In the current context
192 // the notion of symmetry of a shape function means that a
193 // shape function looks exactly the same if it is looked upon
194 // from the centers of the two neighbouring cells that share
195 // it.
196 for (unsigned int local_line_dof = 0;
197 local_line_dof < (k + 1);
198 local_line_dof++)
199 if (local_line_dof % 2 == 0)
200 line_dof_sign[local_line_dof + l * (k + 1)] = -1.0;
201 }
202 }
203 }
204
205 template <int spacedim>
206 void
207 get_dof_sign_change_nedelec(
208 const typename ::Triangulation<3, spacedim>::cell_iterator &cell,
210 const std::vector<MappingKind> &mapping_kind,
211 std::vector<double> &line_dof_sign)
212 {
213 // This function does half of the job - it adjusts the sign of the
214 // line (edge) dofs. In the three-dimensional space the quad (face) dofs
215 // need to be adjusted as well. The quad dofs are treated by
216 // FE_Nedelec<dim>::initialize_quad_dof_index_permutation_and_sign_change()
217 // in fe_nedelec.cc. The dofs associated with the interior of the cells,
218 // the hex dofs, need no adjustments as they are not shared between the
219 // neighboring cells.
220
221 const unsigned int k = fe.tensor_degree() - 1;
222 // The order of the Nedelec elements equals the tensor degree minus one,
223 // k = n - 1. In the three-dimensional space the Nedelec elements of the
224 // lowermost order, k = 0, have only 12 line (edge) dofs. The Nedelec
225 // elements of the higher orders, k > 0, have 3*(k+1)*(k+2)^2 dofs in
226 // total if dim=3. The dofs in a cell are distributed between lines
227 // (edges), quads (faces), and the hex (the interior of the cell) as the
228 // following:
229 //
230 // 12*(k+1) line dofs; (k+1) dofs per line.
231 // 2*6*k*(k+1) quad dofs; 2*k*(k+1) dofs per quad.
232 // 3*(k+2)^2*(k+1) hex dofs.
233 //
234 // The dofs are indexed in the following order: first all line dofs,
235 // then all quad dofs, and then all hex dofs.
236 //
237 // Here we adjust only the line (edge) dofs. The line dofs need only
238 // sign adjustment. That is, no permutation of the line dofs is needed.
239 for (unsigned int l = 0; l < GeometryInfo<3>::lines_per_cell; ++l)
240 if (cell->line_orientation(l) !=
242 mapping_kind[0] == mapping_nedelec)
243 {
244 if (k == 0)
245 {
246 // The lowest order element (k=0) is straightforward, because
247 // there is a single dof per edge, which needs to be flipped:
248 line_dof_sign[l] = -1.0;
249 }
250 else
251 {
252 // The case k > 0 is a bit more complicated. As we adjust
253 // only edge dofs in this function, we need to concern
254 // ourselves with the first 12*(k+1) entries in line_dof_sign
255 // vector ignoring the rest. There are (k+1) dofs per edge.
256 // Let us consider the local dof indices on one edge,
257 // local_line_dof = 0...k. The shape functions with even
258 // indices are asymmetric. The corresponding dofs need sign
259 // adjustment if the edge points in the opposite direction.
260 // The shape functions with odd indices are symmetric. The
261 // corresponding dofs need no sign adjustment even if the edge
262 // points in the opposite direction. In the current context
263 // the notion of symmetry of a shape function means that a
264 // shape function looks exactly the same if it is looked upon
265 // from the centers of the two neighbouring cells that share
266 // it.
267 for (unsigned int local_line_dof = 0;
268 local_line_dof < (k + 1);
269 local_line_dof++)
270 if (local_line_dof % 2 == 0)
271 line_dof_sign[local_line_dof + l * (k + 1)] = -1.0;
272 }
273 }
274 }
275 } // namespace
276 } // namespace FE_PolyTensor
277} // namespace internal
278
279#ifndef DOXYGEN
280
281template <int dim, int spacedim>
283 const TensorPolynomialsBase<dim> &polynomials,
284 const FiniteElementData<dim> &fe_data,
285 const std::vector<bool> &restriction_is_additive_flags,
286 const std::vector<ComponentMask> &nonzero_components)
287 : FiniteElement<dim, spacedim>(fe_data,
288 restriction_is_additive_flags,
289 nonzero_components)
290 , mapping_kind({MappingKind::mapping_none})
291 , poly_space(polynomials.clone())
292{
293 cached_point[0] = -1;
294 // Set up the table converting
295 // components to base
296 // components. Since we have only
297 // one base element, everything
298 // remains zero except the
299 // component in the base, which is
300 // the component itself
301 for (unsigned int comp = 0; comp < this->n_components(); ++comp)
302 this->component_to_base_table[comp].first.second = comp;
303
304 if (dim == 3)
305 {
306 adjust_quad_dof_sign_for_face_orientation_table.resize(
307 this->n_unique_2d_subobjects());
308
309 for (unsigned int f = 0; f < this->n_unique_2d_subobjects(); ++f)
310 {
311 adjust_quad_dof_sign_for_face_orientation_table[f] =
312 Table<2, bool>(this->n_dofs_per_quad(f),
313 this->reference_cell().n_face_orientations(f));
314 adjust_quad_dof_sign_for_face_orientation_table[f].fill(false);
315 }
316 }
317}
318
319
320
321template <int dim, int spacedim>
323 : FiniteElement<dim, spacedim>(fe)
324 , mapping_kind(fe.mapping_kind)
325 , adjust_quad_dof_sign_for_face_orientation_table(
326 fe.adjust_quad_dof_sign_for_face_orientation_table)
327 , poly_space(fe.poly_space->clone())
328 , inverse_node_matrix(fe.inverse_node_matrix)
329{}
330
331
332
333template <int dim, int spacedim>
334bool
336{
337 return mapping_kind.size() == 1;
338}
339
340
341template <int dim, int spacedim>
342bool
344 const unsigned int index,
345 const unsigned int face,
346 const types::geometric_orientation combined_orientation) const
347{
348 // do nothing in 1d and 2d
349 if (dim < 3)
350 return false;
351
352 // The exception are discontinuous
353 // elements for which there should be no
354 // face dofs anyway (i.e. dofs_per_quad==0
355 // in 3d), so we don't need the table, but
356 // the function should also not have been
357 // called
358 AssertIndexRange(index, this->n_dofs_per_quad(face));
359 Assert(adjust_quad_dof_sign_for_face_orientation_table
360 [this->n_unique_2d_subobjects() == 1 ? 0 : face]
361 .n_elements() ==
362 this->reference_cell().n_face_orientations(face) *
363 this->n_dofs_per_quad(face),
365
366 return adjust_quad_dof_sign_for_face_orientation_table
367 [this->n_unique_2d_subobjects() == 1 ? 0 : face](index,
368 combined_orientation);
369}
370
371
372template <int dim, int spacedim>
374FE_PolyTensor<dim, spacedim>::get_mapping_kind(const unsigned int i) const
375{
376 if (single_mapping_kind())
377 return mapping_kind[0];
378
379 AssertIndexRange(i, mapping_kind.size());
380 return mapping_kind[i];
381}
382
383
384
385template <int dim, int spacedim>
386double
388 const Point<dim> &) const
389
390{
392 return 0.;
393}
394
395
396
397template <int dim, int spacedim>
398double
400 const unsigned int i,
401 const Point<dim> &p,
402 const unsigned int component) const
403{
404 AssertIndexRange(i, this->n_dofs_per_cell());
405 AssertIndexRange(component, dim);
406
407 std::lock_guard<std::mutex> lock(cache_mutex);
408
409 if (cached_point != p || cached_values.empty())
410 {
411 cached_point = p;
412 cached_values.resize(poly_space->n());
413
414 std::vector<Tensor<4, dim>> dummy1;
415 std::vector<Tensor<5, dim>> dummy2;
416 poly_space->evaluate(
417 p, cached_values, cached_grads, cached_grad_grads, dummy1, dummy2);
418 }
419
420 double s = 0;
421 if (inverse_node_matrix.n_cols() == 0)
422 return cached_values[i][component];
423 else
424 for (unsigned int j = 0; j < inverse_node_matrix.n_cols(); ++j)
425 s += inverse_node_matrix(j, i) * cached_values[j][component];
426 return s;
427}
428
429
430
431template <int dim, int spacedim>
432Tensor<1, dim>
434 const Point<dim> &) const
435{
437 return Tensor<1, dim>();
438}
439
440
441
442template <int dim, int spacedim>
443Tensor<1, dim>
445 const unsigned int i,
446 const Point<dim> &p,
447 const unsigned int component) const
448{
449 AssertIndexRange(i, this->n_dofs_per_cell());
450 AssertIndexRange(component, dim);
451
452 std::lock_guard<std::mutex> lock(cache_mutex);
453
454 if (cached_point != p || cached_grads.empty())
455 {
456 cached_point = p;
457 cached_grads.resize(poly_space->n());
458
459 std::vector<Tensor<4, dim>> dummy1;
460 std::vector<Tensor<5, dim>> dummy2;
461 poly_space->evaluate(
462 p, cached_values, cached_grads, cached_grad_grads, dummy1, dummy2);
463 }
464
465 Tensor<1, dim> s;
466 if (inverse_node_matrix.n_cols() == 0)
467 return cached_grads[i][component];
468 else
469 for (unsigned int j = 0; j < inverse_node_matrix.n_cols(); ++j)
470 s += inverse_node_matrix(j, i) * cached_grads[j][component];
471
472 return s;
473}
474
475
476
477template <int dim, int spacedim>
478Tensor<2, dim>
480 const Point<dim> &) const
481{
483 return Tensor<2, dim>();
484}
485
486
487
488template <int dim, int spacedim>
489Tensor<2, dim>
491 const unsigned int i,
492 const Point<dim> &p,
493 const unsigned int component) const
494{
495 AssertIndexRange(i, this->n_dofs_per_cell());
496 AssertIndexRange(component, dim);
497
498 std::lock_guard<std::mutex> lock(cache_mutex);
499
500 if (cached_point != p || cached_grad_grads.empty())
501 {
502 cached_point = p;
503 cached_grad_grads.resize(poly_space->n());
504
505 std::vector<Tensor<4, dim>> dummy1;
506 std::vector<Tensor<5, dim>> dummy2;
507 poly_space->evaluate(
508 p, cached_values, cached_grads, cached_grad_grads, dummy1, dummy2);
509 }
510
511 Tensor<2, dim> s;
512 if (inverse_node_matrix.n_cols() == 0)
513 return cached_grad_grads[i][component];
514 else
515 for (unsigned int j = 0; j < inverse_node_matrix.n_cols(); ++j)
516 s += inverse_node_matrix(i, j) * cached_grad_grads[j][component];
517
518 return s;
519}
520
521
522//---------------------------------------------------------------------------
523// Fill data of FEValues
524//---------------------------------------------------------------------------
525
526template <int dim, int spacedim>
527void
530 const CellSimilarity::Similarity /*cell_similarity*/,
531 const Quadrature<dim> &quadrature,
532 const Mapping<dim, spacedim> &mapping,
533 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
534 const internal::FEValuesImplementation::MappingRelatedData<dim, spacedim>
535 &mapping_data,
536 const typename FiniteElement<dim, spacedim>::InternalDataBase &fe_internal,
537 ::internal::FEValuesImplementation::FiniteElementRelatedData<dim,
538 spacedim>
539 &output_data) const
540{
541 Assert(dynamic_cast<const InternalData *>(&fe_internal) != nullptr,
543
544 compute_fill(cell,
546 quadrature.size(),
547 mapping,
548 mapping_internal,
549 mapping_data,
550 static_cast<const InternalData &>(fe_internal),
551 output_data);
552}
553
554
555
556template <int dim, int spacedim>
557void
560 const unsigned int face_no,
561 const hp::QCollection<dim - 1> &quadrature,
562 const Mapping<dim, spacedim> &mapping,
563 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
564 const internal::FEValuesImplementation::MappingRelatedData<dim, spacedim>
565 &mapping_data,
566 const typename FiniteElement<dim, spacedim>::InternalDataBase &fe_internal,
567 ::internal::FEValuesImplementation::FiniteElementRelatedData<dim,
568 spacedim>
569 &output_data) const
570{
571 AssertDimension(quadrature.size(), 1);
572 Assert(dynamic_cast<const InternalData *>(&fe_internal) != nullptr,
574
575 compute_fill(cell,
577 this->reference_cell(),
578 face_no,
579 // TODO: fix the custom implementation of orientation
581 cell->combined_face_orientation(face_no),
582 quadrature[0].size()),
583 quadrature[0].size(),
584 mapping,
585 mapping_internal,
586 mapping_data,
587 static_cast<const InternalData &>(fe_internal),
588 output_data);
589}
590
591
592
593template <int dim, int spacedim>
594void
597 const unsigned int face_no,
598 const unsigned int sub_no,
599 const Quadrature<dim - 1> &quadrature,
600 const Mapping<dim, spacedim> &mapping,
601 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
602 const internal::FEValuesImplementation::MappingRelatedData<dim, spacedim>
603 &mapping_data,
604 const typename FiniteElement<dim, spacedim>::InternalDataBase &fe_internal,
605 ::internal::FEValuesImplementation::FiniteElementRelatedData<dim,
606 spacedim>
607 &output_data) const
608{
609 Assert(dynamic_cast<const InternalData *>(&fe_internal) != nullptr,
611
612 compute_fill(cell,
614 this->reference_cell(),
615 face_no,
616 sub_no,
617 // TODO: fix the custom implementation of orientation
619 cell->combined_face_orientation(face_no),
620 quadrature.size(),
621 cell->subface_case(face_no)),
622 quadrature.size(),
623 mapping,
624 mapping_internal,
625 mapping_data,
626 static_cast<const InternalData &>(fe_internal),
627 output_data);
628}
629
630
631
632template <int dim, int spacedim>
633void
636 const typename QProjector<dim>::DataSetDescriptor &offset,
637 const unsigned int n_q_points,
638 const Mapping<dim, spacedim> &mapping,
639 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
640 const internal::FEValuesImplementation::MappingRelatedData<dim, spacedim>
641 &mapping_data,
642 const typename FE_PolyTensor<dim, spacedim>::InternalData &fe_data,
643 internal::FEValuesImplementation::FiniteElementRelatedData<dim, spacedim>
644 &output_data) const
645{
646 Assert(!(fe_data.update_each & update_values) ||
647 fe_data.shape_values.size()[0] == this->n_dofs_per_cell(),
648 ExcDimensionMismatch(fe_data.shape_values.size()[0],
649 this->n_dofs_per_cell()));
650
651 // TODO: The dof_sign_change only affects Nedelec elements and is not the
652 // correct thing on complicated meshes for higher order Nedelec elements.
653 // Something similar to FE_Q should be done to permute dofs and to change the
654 // dof signs. A static way using tables (as done in the RaviartThomas<dim>
655 // class) is preferable.
656 std::fill(fe_data.dof_sign_change.begin(),
657 fe_data.dof_sign_change.end(),
658 1.0);
659 if (fe_data.update_each & update_values)
660 internal::FE_PolyTensor::get_dof_sign_change_nedelec(
661 cell, *this, this->mapping_kind, fe_data.dof_sign_change);
662
663 // TODO: This, similarly to the Nedelec case, is just a legacy function in 2d
664 // and affects only face_dofs of H(div) conformal FEs. It does nothing in 1d.
665 // Also nothing in 3d since we take care of it by using the
666 // adjust_quad_dof_sign_for_face_orientation_table.
667 if (fe_data.update_each & update_values)
668 internal::FE_PolyTensor::get_dof_sign_change_h_div(cell,
669 *this,
670 this->mapping_kind,
671 fe_data.dof_sign_change);
672
673 // What is the first dof_index on a quad?
674 const unsigned int first_quad_index = this->get_first_quad_index();
675 // How many dofs per quad and how many quad dofs do we have at all?
676 const unsigned int n_dofs_per_quad = this->n_dofs_per_quad();
677 const unsigned int n_quad_dofs =
678 n_dofs_per_quad * GeometryInfo<dim>::faces_per_cell;
679
680 for (unsigned int dof_index = 0; dof_index < this->n_dofs_per_cell();
681 ++dof_index)
682 {
683 /*
684 * This assumes that the dofs are ordered by first vertices, lines, quads
685 * and volume dofs. Note that in 2d this always gives false.
686 */
687 const bool is_quad_dof =
688 (dim == 2 ? false :
689 (first_quad_index <= dof_index) &&
690 (dof_index < first_quad_index + n_quad_dofs));
691
692 // TODO: This hack is not pretty and it is only here to handle the 2d
693 // case and the Nedelec legacy case. In 2d dof_sign of a face_dof is never
694 // handled by the
695 // >>if(is_quad_dof){...}<< but still a possible dof sign change must be
696 // handled, also for line_dofs in 3d such as in Nedelec. In these cases
697 // this is encoded in the array fe_data.dof_sign_change[dof_index]. In 3d
698 // it is handles with a table. This array is allocated in
699 // fe_poly_tensor.h.
700 double dof_sign = 1.0;
701 // under some circumstances fe_data.dof_sign_change is not allocated
702 if (fe_data.update_each & update_values)
703 dof_sign = fe_data.dof_sign_change[dof_index];
704
705 if (is_quad_dof)
706 {
707 /*
708 * Find the face belonging to this dof_index. This is integer
709 * division.
710 */
711 unsigned int face_index_from_dof_index =
712 (dof_index - first_quad_index) / (n_dofs_per_quad);
713
714 unsigned int local_quad_dof_index = dof_index % n_dofs_per_quad;
715
716 // Correct the dof_sign if necessary
717 if (adjust_quad_dof_sign_for_face_orientation(
718 local_quad_dof_index,
719 face_index_from_dof_index,
720 cell->combined_face_orientation(face_index_from_dof_index)))
721 dof_sign = -1.0;
722 }
723
724 const MappingKind mapping_kind = get_mapping_kind(dof_index);
725
726 const unsigned int first =
728 [dof_index * this->n_components() +
729 this->get_nonzero_components(dof_index).first_selected_component()];
730
731 if (fe_data.update_each & update_values)
732 {
733 switch (mapping_kind)
734 {
735 case mapping_none:
736 {
737 for (unsigned int k = 0; k < n_q_points; ++k)
738 for (unsigned int d = 0; d < dim; ++d)
739 output_data.shape_values(first + d, k) =
740 fe_data.shape_values[dof_index][k + offset][d];
741 break;
742 }
743
746 {
747 const ArrayView<Tensor<1, spacedim>>
748 transformed_shape_values =
750 offset,
751 n_q_points);
753 dof_index,
754 offset,
755 n_q_points),
756 mapping_kind,
757 mapping_internal,
758 transformed_shape_values);
759
760 for (unsigned int k = 0; k < n_q_points; ++k)
761 for (unsigned int d = 0; d < dim; ++d)
762 output_data.shape_values(first + d, k) =
763 transformed_shape_values[k][d];
764
765 break;
766 }
768 case mapping_piola:
769 {
770 const ArrayView<Tensor<1, spacedim>>
771 transformed_shape_values =
773 offset,
774 n_q_points);
776 dof_index,
777 offset,
778 n_q_points),
780 mapping_internal,
781 transformed_shape_values);
782 for (unsigned int k = 0; k < n_q_points; ++k)
783 for (unsigned int d = 0; d < dim; ++d)
784 output_data.shape_values(first + d, k) =
785 dof_sign * transformed_shape_values[k][d];
786 break;
787 }
788
789 case mapping_nedelec:
790 {
791 const ArrayView<Tensor<1, spacedim>>
792 transformed_shape_values =
794 offset,
795 n_q_points);
797 dof_index,
798 offset,
799 n_q_points),
801 mapping_internal,
802 transformed_shape_values);
803
804 for (unsigned int k = 0; k < n_q_points; ++k)
805 for (unsigned int d = 0; d < dim; ++d)
806 output_data.shape_values(first + d, k) =
807 dof_sign * transformed_shape_values[k][d];
808
809 break;
810 }
811
812 default:
814 }
815 }
816
817 if (fe_data.update_each & update_gradients)
818 {
819 const ArrayView<Tensor<2, spacedim>> transformed_shape_grads =
821 offset,
822 n_q_points);
823 switch (mapping_kind)
824 {
825 case mapping_none:
826 {
828 dof_index,
829 offset,
830 n_q_points),
832 mapping_internal,
833 transformed_shape_grads);
834 for (unsigned int k = 0; k < n_q_points; ++k)
835 for (unsigned int d = 0; d < dim; ++d)
836 output_data.shape_gradients[first + d][k] =
837 transformed_shape_grads[k][d];
838 break;
839 }
840
842 {
844 dof_index,
845 offset,
846 n_q_points),
848 mapping_internal,
849 transformed_shape_grads);
850
851 for (unsigned int k = 0; k < n_q_points; ++k)
852 for (unsigned int d = 0; d < spacedim; ++d)
853 for (unsigned int n = 0; n < spacedim; ++n)
854 transformed_shape_grads[k][d] -=
855 output_data.shape_values(first + n, k) *
856 mapping_data.jacobian_pushed_forward_grads[k][n][d];
857
858 for (unsigned int k = 0; k < n_q_points; ++k)
859 for (unsigned int d = 0; d < dim; ++d)
860 output_data.shape_gradients[first + d][k] =
861 transformed_shape_grads[k][d];
862 break;
863 }
864
866 {
867 for (unsigned int k = 0; k < n_q_points; ++k)
868 fe_data.untransformed_shape_grads[k + offset] =
869 fe_data.shape_grads[dof_index][k + offset];
872 offset,
873 n_q_points),
875 mapping_internal,
876 transformed_shape_grads);
877
878 for (unsigned int k = 0; k < n_q_points; ++k)
879 for (unsigned int d = 0; d < spacedim; ++d)
880 for (unsigned int n = 0; n < spacedim; ++n)
881 transformed_shape_grads[k][d] +=
882 output_data.shape_values(first + n, k) *
883 mapping_data.jacobian_pushed_forward_grads[k][d][n];
884
885 for (unsigned int k = 0; k < n_q_points; ++k)
886 for (unsigned int d = 0; d < dim; ++d)
887 output_data.shape_gradients[first + d][k] =
888 transformed_shape_grads[k][d];
889
890 break;
891 }
892
894 case mapping_piola:
895 {
896 for (unsigned int k = 0; k < n_q_points; ++k)
897 fe_data.untransformed_shape_grads[k + offset] =
898 fe_data.shape_grads[dof_index][k + offset];
901 offset,
902 n_q_points),
904 mapping_internal,
905 transformed_shape_grads);
906
907 for (unsigned int k = 0; k < n_q_points; ++k)
908 for (unsigned int d = 0; d < spacedim; ++d)
909 for (unsigned int n = 0; n < spacedim; ++n)
910 transformed_shape_grads[k][d] +=
911 (output_data.shape_values(first + n, k) *
912 mapping_data
914 (output_data.shape_values(first + d, k) *
915 mapping_data.jacobian_pushed_forward_grads[k][n][n]);
916
917 for (unsigned int k = 0; k < n_q_points; ++k)
918 for (unsigned int d = 0; d < dim; ++d)
919 output_data.shape_gradients[first + d][k] =
920 dof_sign * transformed_shape_grads[k][d];
921
922 break;
923 }
924
925 case mapping_nedelec:
926 {
927 // treat the gradients of
928 // this particular shape
929 // function at all
930 // q-points. if Dv is the
931 // gradient of the shape
932 // function on the unit
933 // cell, then
934 // (J^-T)Dv(J^-1) is the
935 // value we want to have on
936 // the real cell.
937 for (unsigned int k = 0; k < n_q_points; ++k)
938 fe_data.untransformed_shape_grads[k + offset] =
939 fe_data.shape_grads[dof_index][k + offset];
940
943 offset,
944 n_q_points),
946 mapping_internal,
947 transformed_shape_grads);
948
949 for (unsigned int k = 0; k < n_q_points; ++k)
950 for (unsigned int d = 0; d < spacedim; ++d)
951 for (unsigned int n = 0; n < spacedim; ++n)
952 transformed_shape_grads[k][d] -=
953 output_data.shape_values(first + n, k) *
954 mapping_data.jacobian_pushed_forward_grads[k][n][d];
955
956 for (unsigned int k = 0; k < n_q_points; ++k)
957 for (unsigned int d = 0; d < dim; ++d)
958 output_data.shape_gradients[first + d][k] =
959 dof_sign * transformed_shape_grads[k][d];
960
961 break;
962 }
963
964 default:
966 }
967 }
968
969 if (fe_data.update_each & update_hessians)
970 {
971 switch (mapping_kind)
972 {
973 case mapping_none:
974 {
975 const ArrayView<Tensor<3, spacedim>>
976 transformed_shape_hessians =
978 offset,
979 n_q_points);
981 dof_index,
982 offset,
983 n_q_points),
985 mapping_internal,
986 transformed_shape_hessians);
987
988 for (unsigned int k = 0; k < n_q_points; ++k)
989 for (unsigned int d = 0; d < spacedim; ++d)
990 for (unsigned int n = 0; n < spacedim; ++n)
991 transformed_shape_hessians[k][d] -=
992 output_data.shape_gradients[first + d][k][n] *
993 mapping_data.jacobian_pushed_forward_grads[k][n];
994
995 for (unsigned int k = 0; k < n_q_points; ++k)
996 for (unsigned int d = 0; d < dim; ++d)
997 output_data.shape_hessians[first + d][k] =
998 transformed_shape_hessians[k][d];
999
1000 break;
1001 }
1002 case mapping_covariant:
1003 {
1004 for (unsigned int k = 0; k < n_q_points; ++k)
1005 fe_data.untransformed_shape_hessian_tensors[k + offset] =
1006 fe_data.shape_grad_grads[dof_index][k + offset];
1007
1008 const ArrayView<Tensor<3, spacedim>>
1009 transformed_shape_hessians =
1011 offset,
1012 n_q_points);
1015 offset,
1016 n_q_points),
1018 mapping_internal,
1019 transformed_shape_hessians);
1020
1021 for (unsigned int k = 0; k < n_q_points; ++k)
1022 for (unsigned int d = 0; d < spacedim; ++d)
1023 for (unsigned int n = 0; n < spacedim; ++n)
1024 for (unsigned int i = 0; i < spacedim; ++i)
1025 for (unsigned int j = 0; j < spacedim; ++j)
1026 {
1027 transformed_shape_hessians[k][d][i][j] -=
1028 (output_data.shape_values(first + n, k) *
1029 mapping_data
1030 .jacobian_pushed_forward_2nd_derivatives
1031 [k][n][d][i][j]) +
1032 (output_data.shape_gradients[first + d][k][n] *
1033 mapping_data
1034 .jacobian_pushed_forward_grads[k][n][i][j]) +
1035 (output_data.shape_gradients[first + n][k][i] *
1036 mapping_data
1037 .jacobian_pushed_forward_grads[k][n][d][j]) +
1038 (output_data.shape_gradients[first + n][k][j] *
1039 mapping_data
1040 .jacobian_pushed_forward_grads[k][n][i][d]);
1041 }
1042
1043 for (unsigned int k = 0; k < n_q_points; ++k)
1044 for (unsigned int d = 0; d < dim; ++d)
1045 output_data.shape_hessians[first + d][k] =
1046 transformed_shape_hessians[k][d];
1047
1048 break;
1049 }
1050
1052 {
1053 for (unsigned int k = 0; k < n_q_points; ++k)
1054 fe_data.untransformed_shape_hessian_tensors[k + offset] =
1055 fe_data.shape_grad_grads[dof_index][k + offset];
1056
1057 const ArrayView<Tensor<3, spacedim>>
1058 transformed_shape_hessians =
1060 offset,
1061 n_q_points);
1064 offset,
1065 n_q_points),
1067 mapping_internal,
1068 transformed_shape_hessians);
1069
1070 for (unsigned int k = 0; k < n_q_points; ++k)
1071 for (unsigned int d = 0; d < spacedim; ++d)
1072 for (unsigned int n = 0; n < spacedim; ++n)
1073 for (unsigned int i = 0; i < spacedim; ++i)
1074 for (unsigned int j = 0; j < spacedim; ++j)
1075 {
1076 transformed_shape_hessians[k][d][i][j] +=
1077 (output_data.shape_values(first + n, k) *
1078 mapping_data
1079 .jacobian_pushed_forward_2nd_derivatives
1080 [k][d][n][i][j]) +
1081 (output_data.shape_gradients[first + n][k][i] *
1082 mapping_data
1083 .jacobian_pushed_forward_grads[k][d][n][j]) +
1084 (output_data.shape_gradients[first + n][k][j] *
1085 mapping_data
1086 .jacobian_pushed_forward_grads[k][d][i][n]) -
1087 (output_data.shape_gradients[first + d][k][n] *
1088 mapping_data
1089 .jacobian_pushed_forward_grads[k][n][i][j]);
1090 for (unsigned int m = 0; m < spacedim; ++m)
1091 transformed_shape_hessians[k][d][i][j] -=
1092 (mapping_data
1093 .jacobian_pushed_forward_grads[k][d][i]
1094 [m] *
1095 mapping_data
1096 .jacobian_pushed_forward_grads[k][m][n]
1097 [j] *
1098 output_data.shape_values(first + n, k)) +
1099 (mapping_data
1101 [j] *
1102 mapping_data
1104 [n] *
1105 output_data.shape_values(first + n, k));
1106 }
1107
1108 for (unsigned int k = 0; k < n_q_points; ++k)
1109 for (unsigned int d = 0; d < dim; ++d)
1110 output_data.shape_hessians[first + d][k] =
1111 transformed_shape_hessians[k][d];
1112
1113 break;
1114 }
1115
1117 case mapping_piola:
1118 {
1119 for (unsigned int k = 0; k < n_q_points; ++k)
1120 fe_data.untransformed_shape_hessian_tensors[k + offset] =
1121 fe_data.shape_grad_grads[dof_index][k + offset];
1122
1123 const ArrayView<Tensor<3, spacedim>>
1124 transformed_shape_hessians =
1126 offset,
1127 n_q_points);
1130 offset,
1131 n_q_points),
1133 mapping_internal,
1134 transformed_shape_hessians);
1135
1136 for (unsigned int k = 0; k < n_q_points; ++k)
1137 for (unsigned int d = 0; d < spacedim; ++d)
1138 for (unsigned int n = 0; n < spacedim; ++n)
1139 for (unsigned int i = 0; i < spacedim; ++i)
1140 for (unsigned int j = 0; j < spacedim; ++j)
1141 {
1142 transformed_shape_hessians[k][d][i][j] +=
1143 (output_data.shape_values(first + n, k) *
1144 mapping_data
1145 .jacobian_pushed_forward_2nd_derivatives
1146 [k][d][n][i][j]) +
1147 (output_data.shape_gradients[first + n][k][i] *
1148 mapping_data
1149 .jacobian_pushed_forward_grads[k][d][n][j]) +
1150 (output_data.shape_gradients[first + n][k][j] *
1151 mapping_data
1152 .jacobian_pushed_forward_grads[k][d][i][n]) -
1153 (output_data.shape_gradients[first + d][k][n] *
1154 mapping_data
1155 .jacobian_pushed_forward_grads[k][n][i][j]);
1156
1157 transformed_shape_hessians[k][d][i][j] -=
1158 (output_data.shape_values(first + d, k) *
1159 mapping_data
1160 .jacobian_pushed_forward_2nd_derivatives
1161 [k][n][n][i][j]) +
1162 (output_data.shape_gradients[first + d][k][i] *
1163 mapping_data
1164 .jacobian_pushed_forward_grads[k][n][n][j]) +
1165 (output_data.shape_gradients[first + d][k][j] *
1166 mapping_data
1167 .jacobian_pushed_forward_grads[k][n][n][i]);
1168
1169 for (unsigned int m = 0; m < spacedim; ++m)
1170 {
1171 transformed_shape_hessians[k][d][i][j] -=
1172 (mapping_data
1174 [m] *
1175 mapping_data
1177 [j] *
1178 output_data.shape_values(first + n, k)) +
1179 (mapping_data
1181 [j] *
1182 mapping_data
1184 [n] *
1185 output_data.shape_values(first + n, k));
1186
1187 transformed_shape_hessians[k][d][i][j] +=
1188 (mapping_data
1190 [m] *
1191 mapping_data
1193 [j] *
1194 output_data.shape_values(first + d, k)) +
1195 (mapping_data
1197 [j] *
1198 mapping_data
1200 [n] *
1201 output_data.shape_values(first + d, k));
1202 }
1203 }
1204
1205 for (unsigned int k = 0; k < n_q_points; ++k)
1206 for (unsigned int d = 0; d < dim; ++d)
1207 output_data.shape_hessians[first + d][k] =
1208 dof_sign * transformed_shape_hessians[k][d];
1209
1210 break;
1211 }
1212
1213 case mapping_nedelec:
1214 {
1215 for (unsigned int k = 0; k < n_q_points; ++k)
1216 fe_data.untransformed_shape_hessian_tensors[k + offset] =
1217 fe_data.shape_grad_grads[dof_index][k + offset];
1218
1219 const ArrayView<Tensor<3, spacedim>>
1220 transformed_shape_hessians =
1222 offset,
1223 n_q_points);
1226 offset,
1227 n_q_points),
1229 mapping_internal,
1230 transformed_shape_hessians);
1231
1232 for (unsigned int k = 0; k < n_q_points; ++k)
1233 for (unsigned int d = 0; d < spacedim; ++d)
1234 for (unsigned int n = 0; n < spacedim; ++n)
1235 for (unsigned int i = 0; i < spacedim; ++i)
1236 for (unsigned int j = 0; j < spacedim; ++j)
1237 {
1238 transformed_shape_hessians[k][d][i][j] -=
1239 (output_data.shape_values(first + n, k) *
1240 mapping_data
1241 .jacobian_pushed_forward_2nd_derivatives
1242 [k][n][d][i][j]) +
1243 (output_data.shape_gradients[first + d][k][n] *
1244 mapping_data
1245 .jacobian_pushed_forward_grads[k][n][i][j]) +
1246 (output_data.shape_gradients[first + n][k][i] *
1247 mapping_data
1248 .jacobian_pushed_forward_grads[k][n][d][j]) +
1249 (output_data.shape_gradients[first + n][k][j] *
1250 mapping_data
1251 .jacobian_pushed_forward_grads[k][n][i][d]);
1252 }
1253
1254 for (unsigned int k = 0; k < n_q_points; ++k)
1255 for (unsigned int d = 0; d < dim; ++d)
1256 output_data.shape_hessians[first + d][k] =
1257 dof_sign * transformed_shape_hessians[k][d];
1258
1259 break;
1260 }
1261
1262 default:
1264 }
1265 }
1266
1267 // third derivatives are not implemented
1268 if (fe_data.update_each & update_3rd_derivatives)
1269 {
1271 }
1272 }
1273}
1274
1275
1276
1277template <int dim, int spacedim>
1280 const UpdateFlags flags) const
1281{
1283
1284 for (unsigned int i = 0; i < this->n_dofs_per_cell(); ++i)
1285 {
1286 const MappingKind mapping_kind = get_mapping_kind(i);
1287
1288 switch (mapping_kind)
1289 {
1290 case mapping_none:
1291 {
1292 if (flags & update_values)
1293 out |= update_values;
1294
1295 if (flags & update_gradients)
1298
1299 if (flags & update_hessians)
1303 break;
1304 }
1306 case mapping_piola:
1307 {
1308 if (flags & update_values)
1309 out |= update_values | update_piola;
1310
1311 if (flags & update_gradients)
1316
1317 if (flags & update_hessians)
1322
1323 break;
1324 }
1325
1326
1328 {
1329 if (flags & update_values)
1330 out |= update_values | update_piola;
1331
1332 if (flags & update_gradients)
1337
1338 if (flags & update_hessians)
1343
1344 break;
1345 }
1346
1347 case mapping_nedelec:
1348 case mapping_covariant:
1349 {
1350 if (flags & update_values)
1352
1353 if (flags & update_gradients)
1357
1358 if (flags & update_hessians)
1363
1364 break;
1365 }
1366
1367 default:
1368 {
1370 }
1371 }
1372 }
1373
1374 return out;
1375}
1376
1377#endif
1378// explicit instantiations
1379#include "fe/fe_poly_tensor.inst"
1380
1381
ArrayView< std::remove_reference_t< typename std::iterator_traits< Iterator >::reference >, MemorySpaceType > make_array_view(const Iterator begin, const Iterator end)
Definition array_view.h:954
internal::SubfaceCase< dim > subface_case(const unsigned int face_no) const
std::vector< Tensor< 3, spacedim > > transformed_shape_hessians
Table< 2, DerivativeForm< 2, dim, spacedim > > shape_grad_grads
Table< 2, DerivativeForm< 1, dim, spacedim > > shape_grads
std::vector< Tensor< 2, spacedim > > transformed_shape_grads
std::vector< Tensor< 3, dim > > untransformed_shape_hessian_tensors
std::vector< Tensor< 2, dim > > untransformed_shape_grads
Table< 2, Tensor< 1, dim > > shape_values
std::vector< Tensor< 1, spacedim > > transformed_shape_values
std::vector< double > dof_sign_change
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
bool adjust_quad_dof_sign_for_face_orientation(const unsigned int index, const unsigned int face_no, const types::geometric_orientation combined_orientation) const
virtual Tensor< 1, dim > shape_grad(const unsigned int i, const Point< dim > &p) const override
virtual Tensor< 2, dim > shape_grad_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
virtual Tensor< 2, dim > shape_grad_grad(const unsigned int i, const Point< dim > &p) const override
void compute_fill(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const typename QProjector< dim >::DataSetDescriptor &offset, const unsigned int n_q_points, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const typename FE_PolyTensor< dim, spacedim >::InternalData &fe_internal, internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const
virtual void fill_fe_face_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const hp::QCollection< dim - 1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const typename FiniteElement< dim, spacedim >::InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
virtual double shape_value(const unsigned int i, const Point< dim > &p) const override
virtual void fill_fe_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< dim > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const typename FiniteElement< dim, spacedim >::InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
virtual void fill_fe_subface_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const Quadrature< dim - 1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const typename FiniteElement< dim, spacedim >::InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
MappingKind get_mapping_kind(const unsigned int i) const
FE_PolyTensor(const TensorPolynomialsBase< dim > &polynomials, const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags, const std::vector< ComponentMask > &nonzero_components)
bool single_mapping_kind() const
virtual Tensor< 1, dim > shape_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
unsigned int n_dofs_per_face(unsigned int face_no=0, unsigned int child=0) const
unsigned int tensor_degree() const
virtual unsigned int face_to_cell_index(const unsigned int face_dof_index, const unsigned int face, const types::geometric_orientation combined_orientation=numbers::default_geometric_orientation) const
friend class InternalDataBase
Definition fe.h:3194
virtual void transform(const ArrayView< const Tensor< 1, dim > > &input, const MappingKind kind, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const ArrayView< Tensor< 1, spacedim > > &output) const override
unsigned int size() const
unsigned int size() const
Definition collection.h:316
std::vector< Tensor< 3, spacedim > > jacobian_pushed_forward_grads
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:40
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:41
#define DEAL_II_NOT_IMPLEMENTED()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcFENotPrimitive()
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
TriaIterator< CellAccessor< dim, spacedim > > cell_iterator
Definition tria.h:1557
UpdateFlags
@ update_jacobian_pushed_forward_2nd_derivatives
@ update_contravariant_transformation
Contravariant transformation.
@ update_jacobian_pushed_forward_grads
@ update_hessians
Second derivatives of shape functions.
@ update_values
Shape function values.
@ update_3rd_derivatives
Third derivatives of shape functions.
@ update_covariant_transformation
Covariant transformation.
@ update_gradients
Shape function gradients.
@ update_default
No update.
@ update_piola
Values needed for Piola transform.
MappingKind
Definition mapping.h:81
@ mapping_piola
Definition mapping.h:116
@ mapping_covariant_gradient
Definition mapping.h:102
@ mapping_covariant
Definition mapping.h:91
@ mapping_nedelec
Definition mapping.h:131
@ mapping_contravariant
Definition mapping.h:96
@ mapping_raviart_thomas
Definition mapping.h:136
@ mapping_contravariant_hessian
Definition mapping.h:158
@ mapping_covariant_hessian
Definition mapping.h:152
@ mapping_contravariant_gradient
Definition mapping.h:108
@ mapping_none
Definition mapping.h:86
@ mapping_piola_gradient
Definition mapping.h:122
@ mapping_piola_hessian
Definition mapping.h:164
types::geometric_orientation combined_face_orientation(const unsigned int face) const
MappingQ< dim, spacedim > StaticMappingQ1< dim, spacedim >::mapping
Definition mapping_q1.h:104
void reference_cell(Triangulation< dim, spacedim > &tria, const ReferenceCell &reference_cell)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
constexpr types::geometric_orientation default_geometric_orientation
Definition types.h:352
unsigned char geometric_orientation
Definition types.h:40
static constexpr unsigned int faces_per_cell