deal.II version GIT relicensing-2289-g1e5549a87a 2024-12-21 21:30:00+00:00
\(\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.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2004 - 2023 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
16#include <deal.II/base/config.h>
17
26
27#include <deal.II/fe/fe_poly.h>
30
32
33#ifndef DOXYGEN
34
35template <int dim, int spacedim>
37 : FiniteElement<dim, spacedim>(fe)
38 , poly_space(fe.poly_space->clone())
39{}
40
41template <int dim, int spacedim>
43 const ScalarPolynomialsBase<dim> &poly_space,
44 const FiniteElementData<dim> &fe_data,
45 const std::vector<bool> &restriction_is_additive_flags,
46 const std::vector<ComponentMask> &nonzero_components)
47 : FiniteElement<dim, spacedim>(fe_data,
48 restriction_is_additive_flags,
49 nonzero_components)
50 , poly_space(poly_space.clone())
51{}
52
53
54template <int dim, int spacedim>
55unsigned int
57{
58 return this->degree;
59}
60
61
62template <int dim, int spacedim>
63double
64FE_Poly<dim, spacedim>::shape_value(const unsigned int i,
65 const Point<dim> &p) const
66{
67 AssertIndexRange(i, this->n_dofs_per_cell());
68 return poly_space->compute_value(i, p);
69}
70
71
72template <int dim, int spacedim>
73double
75 const unsigned int i,
76 const Point<dim> &p,
77 const unsigned int component) const
78{
79 (void)component;
80 AssertIndexRange(i, this->n_dofs_per_cell());
81 AssertIndexRange(component, 1);
82 return poly_space->compute_value(i, p);
83}
84
85
86
87template <int dim, int spacedim>
89FE_Poly<dim, spacedim>::shape_grad(const unsigned int i,
90 const Point<dim> &p) const
91{
92 AssertIndexRange(i, this->n_dofs_per_cell());
93 return poly_space->template compute_derivative<1>(i, p);
94}
95
96
97
98template <int dim, int spacedim>
101 const Point<dim> &p,
102 const unsigned int component) const
103{
104 (void)component;
105 AssertIndexRange(i, this->n_dofs_per_cell());
106 AssertIndexRange(component, 1);
107 return poly_space->template compute_derivative<1>(i, p);
108}
109
110
111
112template <int dim, int spacedim>
114FE_Poly<dim, spacedim>::shape_grad_grad(const unsigned int i,
115 const Point<dim> &p) const
116{
117 AssertIndexRange(i, this->n_dofs_per_cell());
118 return poly_space->template compute_derivative<2>(i, p);
119}
120
121
122
123template <int dim, int spacedim>
126 const unsigned int i,
127 const Point<dim> &p,
128 const unsigned int component) const
129{
130 (void)component;
131 AssertIndexRange(i, this->n_dofs_per_cell());
132 AssertIndexRange(component, 1);
133 return poly_space->template compute_derivative<2>(i, p);
134}
135
136
137
138template <int dim, int spacedim>
141 const Point<dim> &p) const
142{
143 AssertIndexRange(i, this->n_dofs_per_cell());
144 return poly_space->template compute_derivative<3>(i, p);
145}
146
147
148
149template <int dim, int spacedim>
152 const unsigned int i,
153 const Point<dim> &p,
154 const unsigned int component) const
155{
156 (void)component;
157 AssertIndexRange(i, this->n_dofs_per_cell());
158 AssertIndexRange(component, 1);
159 return poly_space->template compute_derivative<3>(i, p);
160}
161
162
163
164template <int dim, int spacedim>
167 const Point<dim> &p) const
168{
169 AssertIndexRange(i, this->n_dofs_per_cell());
170 return poly_space->template compute_derivative<4>(i, p);
171}
172
173
174
175template <int dim, int spacedim>
178 const unsigned int i,
179 const Point<dim> &p,
180 const unsigned int component) const
181{
182 (void)component;
183 AssertIndexRange(i, this->n_dofs_per_cell());
184 AssertIndexRange(component, 1);
185 return poly_space->template compute_derivative<4>(i, p);
186}
187
188
189
190//---------------------------------------------------------------------------
191// Auxiliary functions
192//---------------------------------------------------------------------------
193
194
195template <int dim, int spacedim>
198{
200
201 if (flags & update_values)
202 out |= update_values;
203 if (flags & update_gradients)
205 if (flags & update_hessians)
208 if (flags & update_3rd_derivatives)
213 if (flags & update_normal_vectors)
215
216 return out;
217}
218
219
220
221//---------------------------------------------------------------------------
222// Fill data of FEValues
223//---------------------------------------------------------------------------
224
225
226
236template <int dim, int spacedim>
237bool
238higher_derivatives_need_correcting(
239 const Mapping<dim, spacedim> &mapping,
241 &mapping_data,
242 const unsigned int n_q_points,
243 const UpdateFlags update_flags)
244{
245 // If higher derivatives weren't requested we don't need to correct them.
246 const bool update_higher_derivatives =
247 (update_flags & update_hessians) || (update_flags & update_3rd_derivatives);
248 if (!update_higher_derivatives)
249 return false;
250
251 // If we have a Cartesian mapping, we know that jacoban_pushed_forward_grads
252 // are identically zero.
253 if (dynamic_cast<const MappingCartesian<dim> *>(&mapping))
254 return false;
255
256 // Here, we should check if jacobian_pushed_forward_grads are zero at the
257 // quadrature points. This is yet to be implemented.
258 (void)mapping_data;
259 (void)n_q_points;
260
261 return true;
262}
263
264
265
266template <int dim, int spacedim>
267void
270 const CellSimilarity::Similarity cell_similarity,
271 const Quadrature<dim> &quadrature,
272 const Mapping<dim, spacedim> &mapping,
273 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
275 &mapping_data,
276 const typename FiniteElement<dim, spacedim>::InternalDataBase &fe_internal,
278 spacedim>
279 &output_data) const
280{
281 // convert data object to internal
282 // data for this class. fails with
283 // an exception if that is not
284 // possible
285 Assert(dynamic_cast<const InternalData *>(&fe_internal) != nullptr,
287 const InternalData &fe_data = static_cast<const InternalData &>(fe_internal);
288
289 const UpdateFlags flags(fe_data.update_each);
290
291 const bool need_to_correct_higher_derivatives =
292 higher_derivatives_need_correcting(mapping,
293 mapping_data,
294 quadrature.size(),
295 flags);
296
297 // transform gradients and higher derivatives. there is nothing to do
298 // for values since we already emplaced them into output_data when
299 // we were in get_data()
300 if ((flags & update_gradients) &&
301 (cell_similarity != CellSimilarity::translation))
302 for (unsigned int k = 0; k < this->n_dofs_per_cell(); ++k)
303 mapping.transform(make_array_view(fe_data.shape_gradients, k),
305 mapping_internal,
306 make_array_view(output_data.shape_gradients, k));
307
308 if ((flags & update_hessians) &&
309 (cell_similarity != CellSimilarity::translation))
310 {
311 for (unsigned int k = 0; k < this->n_dofs_per_cell(); ++k)
312 mapping.transform(make_array_view(fe_data.shape_hessians, k),
314 mapping_internal,
315 make_array_view(output_data.shape_hessians, k));
316
317 if (need_to_correct_higher_derivatives)
318 correct_hessians(output_data, mapping_data, quadrature.size());
319 }
320
321 if ((flags & update_3rd_derivatives) &&
322 (cell_similarity != CellSimilarity::translation))
323 {
324 for (unsigned int k = 0; k < this->n_dofs_per_cell(); ++k)
325 mapping.transform(make_array_view(fe_data.shape_3rd_derivatives, k),
327 mapping_internal,
328 make_array_view(output_data.shape_3rd_derivatives,
329 k));
330
331 if (need_to_correct_higher_derivatives)
332 correct_third_derivatives(output_data, mapping_data, quadrature.size());
333 }
334}
335
336
337
338template <int dim, int spacedim>
339void
342 const unsigned int face_no,
343 const hp::QCollection<dim - 1> &quadrature,
344 const Mapping<dim, spacedim> &mapping,
345 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
347 &mapping_data,
348 const typename FiniteElement<dim, spacedim>::InternalDataBase &fe_internal,
350 spacedim>
351 &output_data) const
352{
353 const unsigned int n_q_points =
354 quadrature[quadrature.size() == 1 ? 0 : face_no].size();
355
356 // convert data object to internal
357 // data for this class. fails with
358 // an exception if that is not
359 // possible
360 Assert(dynamic_cast<const InternalData *>(&fe_internal) != nullptr,
362 const InternalData &fe_data = static_cast<const InternalData &>(fe_internal);
363
364 // offset determines which data set
365 // to take (all data sets for all
366 // faces are stored contiguously)
367
368 const auto offset =
370 face_no,
371 cell->face_orientation(face_no),
372 cell->face_flip(face_no),
373 cell->face_rotation(face_no),
374 quadrature);
375
376 const UpdateFlags flags(fe_data.update_each);
377
378 const bool need_to_correct_higher_derivatives =
379 higher_derivatives_need_correcting(mapping,
380 mapping_data,
381 n_q_points,
382 flags);
383
384 // transform gradients and higher derivatives. we also have to copy
385 // the values (unlike in the case of fill_fe_values()) since
386 // we need to take into account the offsets
387 if (flags & update_values)
388 for (unsigned int k = 0; k < this->n_dofs_per_cell(); ++k)
389 for (unsigned int i = 0; i < n_q_points; ++i)
390 output_data.shape_values(k, i) = fe_data.shape_values[k][i + offset];
391
392 if (flags & update_gradients)
393 for (unsigned int k = 0; k < this->n_dofs_per_cell(); ++k)
394 mapping.transform(
395 make_array_view(fe_data.shape_gradients, k, offset, n_q_points),
397 mapping_internal,
398 make_array_view(output_data.shape_gradients, k));
399
400 if (flags & update_hessians)
401 {
402 for (unsigned int k = 0; k < this->n_dofs_per_cell(); ++k)
403 mapping.transform(
404 make_array_view(fe_data.shape_hessians, k, offset, n_q_points),
406 mapping_internal,
407 make_array_view(output_data.shape_hessians, k));
408
409 if (need_to_correct_higher_derivatives)
410 correct_hessians(output_data, mapping_data, n_q_points);
411 }
412
413 if (flags & update_3rd_derivatives)
414 {
415 for (unsigned int k = 0; k < this->n_dofs_per_cell(); ++k)
416 mapping.transform(
417 make_array_view(fe_data.shape_3rd_derivatives, k, offset, n_q_points),
419 mapping_internal,
420 make_array_view(output_data.shape_3rd_derivatives, k));
421
422 if (need_to_correct_higher_derivatives)
423 correct_third_derivatives(output_data, mapping_data, n_q_points);
424 }
425}
426
427
428
429template <int dim, int spacedim>
430void
433 const unsigned int face_no,
434 const unsigned int sub_no,
435 const Quadrature<dim - 1> &quadrature,
436 const Mapping<dim, spacedim> &mapping,
437 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
439 &mapping_data,
440 const typename FiniteElement<dim, spacedim>::InternalDataBase &fe_internal,
442 spacedim>
443 &output_data) const
444{
445 // convert data object to internal
446 // data for this class. fails with
447 // an exception if that is not
448 // possible
449 Assert(dynamic_cast<const InternalData *>(&fe_internal) != nullptr,
451 const InternalData &fe_data = static_cast<const InternalData &>(fe_internal);
452
453 // offset determines which data set
454 // to take (all data sets for all
455 // sub-faces are stored contiguously)
456
457 const auto offset =
459 face_no,
460 sub_no,
461 cell->face_orientation(face_no),
462 cell->face_flip(face_no),
463 cell->face_rotation(face_no),
464 quadrature.size(),
465 cell->subface_case(face_no));
466
467 const UpdateFlags flags(fe_data.update_each);
468
469 const bool need_to_correct_higher_derivatives =
470 higher_derivatives_need_correcting(mapping,
471 mapping_data,
472 quadrature.size(),
473 flags);
474
475 // transform gradients and higher derivatives. we also have to copy
476 // the values (unlike in the case of fill_fe_values()) since
477 // we need to take into account the offsets
478 if (flags & update_values)
479 for (unsigned int k = 0; k < this->n_dofs_per_cell(); ++k)
480 for (unsigned int i = 0; i < quadrature.size(); ++i)
481 output_data.shape_values(k, i) = fe_data.shape_values[k][i + offset];
482
483 if (flags & update_gradients)
484 for (unsigned int k = 0; k < this->n_dofs_per_cell(); ++k)
485 mapping.transform(
486 make_array_view(fe_data.shape_gradients, k, offset, quadrature.size()),
488 mapping_internal,
489 make_array_view(output_data.shape_gradients, k));
490
491 if (flags & update_hessians)
492 {
493 for (unsigned int k = 0; k < this->n_dofs_per_cell(); ++k)
494 mapping.transform(
495 make_array_view(fe_data.shape_hessians, k, offset, quadrature.size()),
497 mapping_internal,
498 make_array_view(output_data.shape_hessians, k));
499
500 if (need_to_correct_higher_derivatives)
501 correct_hessians(output_data, mapping_data, quadrature.size());
502 }
503
504 if (flags & update_3rd_derivatives)
505 {
506 for (unsigned int k = 0; k < this->n_dofs_per_cell(); ++k)
507 mapping.transform(make_array_view(fe_data.shape_3rd_derivatives,
508 k,
509 offset,
510 quadrature.size()),
512 mapping_internal,
513 make_array_view(output_data.shape_3rd_derivatives,
514 k));
515
516 if (need_to_correct_higher_derivatives)
517 correct_third_derivatives(output_data, mapping_data, quadrature.size());
518 }
519}
520
521
522
523template <int dim, int spacedim>
524inline void
527 &output_data,
529 &mapping_data,
530 const unsigned int n_q_points) const
531{
532 for (unsigned int dof = 0; dof < this->n_dofs_per_cell(); ++dof)
533 for (unsigned int i = 0; i < n_q_points; ++i)
534 for (unsigned int j = 0; j < spacedim; ++j)
535 output_data.shape_hessians[dof][i] -=
536 mapping_data.jacobian_pushed_forward_grads[i][j] *
537 output_data.shape_gradients[dof][i][j];
538}
539
540
541
542template <int dim, int spacedim>
543inline void
546 &output_data,
548 &mapping_data,
549 const unsigned int n_q_points) const
550{
551 for (unsigned int dof = 0; dof < this->n_dofs_per_cell(); ++dof)
552 for (unsigned int i = 0; i < n_q_points; ++i)
553 for (unsigned int j = 0; j < spacedim; ++j)
554 for (unsigned int k = 0; k < spacedim; ++k)
555 for (unsigned int l = 0; l < spacedim; ++l)
556 for (unsigned int m = 0; m < spacedim; ++m)
557 {
558 output_data.shape_3rd_derivatives[dof][i][j][k][l] -=
559 (mapping_data.jacobian_pushed_forward_grads[i][m][j][l] *
560 output_data.shape_hessians[dof][i][k][m]) +
561 (mapping_data.jacobian_pushed_forward_grads[i][m][k][l] *
562 output_data.shape_hessians[dof][i][j][m]) +
563 (mapping_data.jacobian_pushed_forward_grads[i][m][j][k] *
564 output_data.shape_hessians[dof][i][l][m]) +
565 (mapping_data
566 .jacobian_pushed_forward_2nd_derivatives[i][m][j][k][l] *
567 output_data.shape_gradients[dof][i][m]);
568 }
569}
570
571
572
573template <int dim, int spacedim>
574inline const ScalarPolynomialsBase<dim> &
576{
577 return *poly_space;
578}
579
580
581
582template <int dim, int spacedim>
583std::vector<unsigned int>
585{
586 auto *const space_tensor_prod =
587 dynamic_cast<TensorProductPolynomials<dim> *>(this->poly_space.get());
588 if (space_tensor_prod != nullptr)
589 return space_tensor_prod->get_numbering();
590
591 auto *const space_tensor_prod_aniso =
592 dynamic_cast<AnisotropicPolynomials<dim> *>(this->poly_space.get());
593 if (space_tensor_prod_aniso != nullptr)
594 return space_tensor_prod_aniso->get_numbering();
595
596 auto *const space_tensor_prod_piecewise = dynamic_cast<
598 this->poly_space.get());
599 if (space_tensor_prod_piecewise != nullptr)
600 return space_tensor_prod_piecewise->get_numbering();
601
602 auto *const space_tensor_prod_bubbles =
604 this->poly_space.get());
605 if (space_tensor_prod_bubbles != nullptr)
606 return space_tensor_prod_bubbles->get_numbering();
607
608 auto *const space_tensor_prod_const =
609 dynamic_cast<TensorProductPolynomialsConst<dim> *>(this->poly_space.get());
610 if (space_tensor_prod_const != nullptr)
611 return space_tensor_prod_const->get_numbering();
612
614 return std::vector<unsigned int>();
615}
616
617
618
619template <int dim, int spacedim>
620std::vector<unsigned int>
622{
623 return Utilities::invert_permutation(get_poly_space_numbering());
624}
625
626
627
628template <int dim, int spacedim>
629std::size_t
631{
633 poly_space->memory_consumption();
634}
635
636#endif
637
638#include "fe_poly.inst"
639
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:949
const std::vector< unsigned int > & get_numbering() const
virtual Tensor< 2, dim > shape_grad_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 Tensor< 3, dim > shape_3rd_derivative_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
std::vector< unsigned int > get_poly_space_numbering_inverse() const
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
FE_Poly(const ScalarPolynomialsBase< dim > &poly_space, const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags, const std::vector< ComponentMask > &nonzero_components)
void correct_hessians(internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data, const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const unsigned int n_q_points) const
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 Tensor< 1, dim > shape_grad(const unsigned int i, const Point< dim > &p) const override
unsigned int get_degree() const
const ScalarPolynomialsBase< dim > & get_poly_space() const
virtual Tensor< 4, dim > shape_4th_derivative_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
void correct_third_derivatives(internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data, const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const unsigned int n_q_points) const
virtual double shape_value(const unsigned int i, const Point< dim > &p) const override
virtual Tensor< 3, dim > shape_3rd_derivative(const unsigned int i, const Point< dim > &p) const override
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 UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
virtual Tensor< 4, dim > shape_4th_derivative(const unsigned int i, const Point< dim > &p) const override
std::vector< unsigned int > get_poly_space_numbering() const
virtual std::size_t memory_consumption() 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
virtual Tensor< 1, dim > shape_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
virtual std::size_t memory_consumption() const
Abstract base class for mapping classes.
Definition mapping.h:318
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 =0
Definition point.h:111
static DataSetDescriptor subface(const ReferenceCell &reference_cell, const unsigned int face_no, const unsigned int subface_no, const bool face_orientation, const bool face_flip, const bool face_rotation, const unsigned int n_quadrature_points, const internal::SubfaceCase< dim > ref_case=internal::SubfaceCase< dim >::case_isotropic)
static DataSetDescriptor face(const ReferenceCell &reference_cell, const unsigned int face_no, const bool face_orientation, const bool face_flip, const bool face_rotation, const unsigned int n_quadrature_points)
unsigned int size() const
const std::vector< unsigned int > & get_numbering() const
const std::vector< unsigned int > & get_numbering() const
const std::vector< unsigned int > & get_numbering() const
unsigned int size() const
Definition collection.h:315
std::vector< Tensor< 3, spacedim > > jacobian_pushed_forward_grads
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:498
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499
#define Assert(cond, exc)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
UpdateFlags
@ update_jacobian_pushed_forward_2nd_derivatives
@ update_jacobian_pushed_forward_grads
@ update_hessians
Second derivatives of shape functions.
@ update_values
Shape function values.
@ update_normal_vectors
Normal vectors.
@ update_3rd_derivatives
Third derivatives of shape functions.
@ update_JxW_values
Transformed quadrature weights.
@ update_covariant_transformation
Covariant transformation.
@ update_gradients
Shape function gradients.
@ update_default
No update.
@ mapping_covariant_gradient
Definition mapping.h:100
@ mapping_covariant
Definition mapping.h:89
@ mapping_covariant_hessian
Definition mapping.h:150
#define DEAL_II_NOT_IMPLEMENTED()
void reference_cell(Triangulation< dim, spacedim > &tria, const ReferenceCell &reference_cell)
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:1699