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
derivative_approximation.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2000 - 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
18
21
22#include <deal.II/fe/fe.h>
24#include <deal.II/fe/mapping.h>
25
30
35
46#include <deal.II/lac/vector.h>
47
49
50#include <algorithm>
51#include <cmath>
52
54
55
56
57namespace
58{
59 template <typename T>
60 inline T
61 sqr(const T t)
62 {
63 return t * t;
64 }
65} // namespace
66
67// --- First the classes and functions that describe individual derivatives ---
68
70{
71 namespace internal
72 {
79 template <int dim>
81 {
82 public:
87 static const UpdateFlags update_flags;
88
94
100
107 template <class InputVector, int spacedim>
110 const InputVector &solution,
111 const unsigned int component);
112
117 static double
118 derivative_norm(const Derivative &d);
119
127 static void
128 symmetrize(Derivative &derivative_tensor);
129 };
130
131 // static variables
132 template <int dim>
134
135
136 template <int dim>
137 template <class InputVector, int spacedim>
140 const FEValues<dim, spacedim> &fe_values,
141 const InputVector &solution,
142 const unsigned int component)
143 {
144 if (fe_values.get_fe().n_components() == 1)
145 {
146 std::vector<typename InputVector::value_type> values(1);
147 fe_values.get_function_values(solution, values);
148 return values[0];
149 }
150 else
151 {
152 std::vector<Vector<typename InputVector::value_type>> values(
153 1,
155 fe_values.get_fe().n_components()));
156 fe_values.get_function_values(solution, values);
157 return values[0](component);
158 }
159 }
160
161
162
163 template <int dim>
164 inline double
166 {
167 double s = 0;
168 for (unsigned int i = 0; i < dim; ++i)
169 s += d[i] * d[i];
170 return std::sqrt(s);
171 }
172
173
174
175 template <int dim>
176 inline void
178 {
179 // nothing to do here
180 }
181
182
183
190 template <int dim>
192 {
193 public:
198 static const UpdateFlags update_flags;
199
205
211
218 template <class InputVector, int spacedim>
221 const InputVector &solution,
222 const unsigned int component);
223
231 static double
232 derivative_norm(const Derivative &d);
233
243 static void
244 symmetrize(Derivative &derivative_tensor);
245 };
246
247 template <int dim>
249
250
251 template <int dim>
252 template <class InputVector, int spacedim>
255 const FEValues<dim, spacedim> &fe_values,
256 const InputVector &solution,
257 const unsigned int component)
258 {
259 if (fe_values.get_fe().n_components() == 1)
260 {
261 std::vector<Tensor<1, dim, typename InputVector::value_type>> values(
262 1);
263 fe_values.get_function_gradients(solution, values);
264 return ProjectedDerivative(values[0]);
265 }
266 else
267 {
268 std::vector<
269 std::vector<Tensor<1, dim, typename InputVector::value_type>>>
270 values(
271 1,
273 fe_values.get_fe().n_components()));
274 fe_values.get_function_gradients(solution, values);
275 return ProjectedDerivative(values[0][component]);
276 }
277 }
278
279
280
281 template <>
282 inline double
284 {
285 return std::fabs(d[0][0]);
286 }
287
288
289
290 template <>
291 inline double
293 {
294 // note that d should be a
295 // symmetric 2x2 tensor, so the
296 // eigenvalues are:
297 //
298 // 1/2(a+b\pm\sqrt((a-b)^2+4c^2))
299 //
300 // if the d_11=a, d_22=b,
301 // d_12=d_21=c
302 const double radicand =
303 ::sqr(d[0][0] - d[1][1]) + 4 * ::sqr(d[0][1]);
304 const double eigenvalues[2] = {
305 0.5 * (d[0][0] + d[1][1] + std::sqrt(radicand)),
306 0.5 * (d[0][0] + d[1][1] - std::sqrt(radicand))};
307
308 return std::max(std::fabs(eigenvalues[0]), std::fabs(eigenvalues[1]));
309 }
310
311
312
313 template <>
314 inline double
316 {
317 /*
318 compute the three eigenvalues of the tensor @p{d} and take the
319 largest. one could use the following maple script to generate C
320 code:
321
322 with(linalg);
323 readlib(C);
324 A:=matrix(3,3,[[a00,a01,a02],[a01,a11,a12],[a02,a12,a22]]);
325 E:=eigenvals(A);
326 EE:=vector(3,[E[1],E[2],E[3]]);
327 C(EE);
328
329 Unfortunately, with both optimized and non-optimized output, at some
330 places the code `sqrt(-1.0)' is emitted, and I don't know what
331 Maple intends to do with it. This happens both with Maple4 and
332 Maple5.
333
334 Fortunately, Roger Young provided the following Fortran code, which
335 is transcribed below to C. The code uses an algorithm that uses the
336 invariants of a symmetric matrix. (The translated algorithm is
337 augmented by a test for R>0, since R==0 indicates that all three
338 eigenvalues are equal.)
339
340
341 PROGRAM MAIN
342
343 C FIND EIGENVALUES OF REAL SYMMETRIC MATRIX
344 C (ROGER YOUNG, 2001)
345
346 IMPLICIT NONE
347
348 REAL*8 A11, A12, A13, A22, A23, A33
349 REAL*8 I1, J2, J3, AM
350 REAL*8 S11, S12, S13, S22, S23, S33
351 REAL*8 SS12, SS23, SS13
352 REAL*8 R,R3, XX,YY, THETA
353 REAL*8 A1,A2,A3
354 REAL*8 PI
355 PARAMETER (PI=3.141592653587932384D0)
356 REAL*8 A,B,C, TOL
357 PARAMETER (TOL=1.D-14)
358
359 C DEFINE A TEST MATRIX
360
361 A11 = -1.D0
362 A12 = 5.D0
363 A13 = 3.D0
364 A22 = -2.D0
365 A23 = 0.5D0
366 A33 = 4.D0
367
368
369 I1 = A11 + A22 + A33
370 AM = I1/3.D0
371
372 S11 = A11 - AM
373 S22 = A22 - AM
374 S33 = A33 - AM
375 S12 = A12
376 S13 = A13
377 S23 = A23
378
379 SS12 = S12*S12
380 SS23 = S23*S23
381 SS13 = S13*S13
382
383 J2 = S11*S11 + S22*S22 + S33*S33
384 J2 = J2 + 2.D0*(SS12 + SS23 + SS13)
385 J2 = J2/2.D0
386
387 J3 = S11**3 + S22**3 + S33**3
388 J3 = J3 + 3.D0*S11*(SS12 + SS13)
389 J3 = J3 + 3.D0*S22*(SS12 + SS23)
390 J3 = J3 + 3.D0*S33*(SS13 + SS23)
391 J3 = J3 + 6.D0*S12*S23*S13
392 J3 = J3/3.D0
393
394 R = SQRT(4.D0*J2/3.D0)
395 R3 = R*R*R
396 XX = 4.D0*J3/R3
397
398 YY = 1.D0 - DABS(XX)
399 IF(YY.LE.0.D0)THEN
400 IF(YY.GT.(-TOL))THEN
401 WRITE(6,*)'Equal roots: XX= ',XX
402 A = -(XX/DABS(XX))*SQRT(J2/3.D0)
403 B = AM + A
404 C = AM - 2.D0*A
405 WRITE(6,*)B,' (twice) ',C
406 STOP
407 ELSE
408 WRITE(6,*)'Error: XX= ',XX
409 STOP
410 ENDIF
411 ENDIF
412
413 THETA = (ACOS(XX))/3.D0
414
415 A1 = AM + R*COS(THETA)
416 A2 = AM + R*COS(THETA + 2.D0*PI/3.D0)
417 A3 = AM + R*COS(THETA + 4.D0*PI/3.D0)
418
419 WRITE(6,*)A1,A2,A3
420
421 STOP
422 END
423
424 */
425
426 const double am = trace(d) / 3.;
427
428 // s := d - trace(d) I
429 Tensor<2, 3> s = d;
430 for (unsigned int i = 0; i < 3; ++i)
431 s[i][i] -= am;
432
433 const double ss01 = s[0][1] * s[0][1], ss12 = s[1][2] * s[1][2],
434 ss02 = s[0][2] * s[0][2];
435
436 const double J2 = (s[0][0] * s[0][0] + s[1][1] * s[1][1] +
437 s[2][2] * s[2][2] + 2 * (ss01 + ss02 + ss12)) /
438 2.;
439 const double J3 =
440 (Utilities::fixed_power<3>(s[0][0]) +
442 Utilities::fixed_power<3>(s[2][2]) + 3. * s[0][0] * (ss01 + ss02) +
443 3. * s[1][1] * (ss01 + ss12) + 3. * s[2][2] * (ss02 + ss12) +
444 6. * s[0][1] * s[0][2] * s[1][2]) /
445 3.;
446
447 const double R = std::sqrt(4. * J2 / 3.);
448
449 double EE[3] = {0, 0, 0};
450 // the eigenvalues are away from
451 // @p{am} in the order of R. thus,
452 // if R<<AM, then we have the
453 // degenerate case with three
454 // identical eigenvalues. check
455 // this first
456 if (R <= 1e-14 * std::fabs(am))
457 EE[0] = EE[1] = EE[2] = am;
458 else
459 {
460 // at least two eigenvalues are
461 // distinct
462 const double R3 = R * R * R;
463 const double XX = 4. * J3 / R3;
464 const double YY = 1. - std::fabs(XX);
465
466 Assert(YY > -1e-14, ExcInternalError());
467
468 if (YY < 0)
469 {
470 // two roots are equal
471 const double a = (XX > 0 ? -1. : 1.) * R / 2;
472 EE[0] = EE[1] = am + a;
473 EE[2] = am - 2. * a;
474 }
475 else
476 {
477 const double theta = std::acos(XX) / 3.;
478 EE[0] = am + R * std::cos(theta);
479 EE[1] = am + R * std::cos(theta + 2. / 3. * numbers::PI);
480 EE[2] = am + R * std::cos(theta + 4. / 3. * numbers::PI);
481 }
482 }
483
484 return std::max({std::fabs(EE[0]), std::fabs(EE[1]), std::fabs(EE[2])});
485 }
486
487
488
489 template <int dim>
490 inline double
492 {
493 // computing the spectral norm is
494 // not so simple in general. it is
495 // feasible for dim==3 as shown
496 // above, since then there are
497 // still closed form expressions of
498 // the roots of the characteristic
499 // polynomial, and they can easily
500 // be computed using
501 // maple. however, for higher
502 // dimensions, some other method
503 // needs to be employed. maybe some
504 // steps of the power method would
505 // suffice?
507 return 0;
508 }
509
510
511
512 template <int dim>
513 inline void
515 {
516 // symmetrize non-diagonal entries
517 for (unsigned int i = 0; i < dim; ++i)
518 for (unsigned int j = i + 1; j < dim; ++j)
519 {
520 const double s = (d[i][j] + d[j][i]) / 2;
521 d[i][j] = d[j][i] = s;
522 }
523 }
524
525
526
527 template <int dim>
529 {
530 public:
535 static const UpdateFlags update_flags;
536
543
549
556 template <class InputVector, int spacedim>
559 const InputVector &solution,
560 const unsigned int component);
561
569 static double
570 derivative_norm(const Derivative &d);
571
581 static void
582 symmetrize(Derivative &derivative_tensor);
583 };
584
585 template <int dim>
587
588
589 template <int dim>
590 template <class InputVector, int spacedim>
593 const FEValues<dim, spacedim> &fe_values,
594 const InputVector &solution,
595 const unsigned int component)
596 {
597 if (fe_values.get_fe().n_components() == 1)
598 {
599 std::vector<Tensor<2, dim, typename InputVector::value_type>> values(
600 1);
601 fe_values.get_function_hessians(solution, values);
602 return ProjectedDerivative(values[0]);
603 }
604 else
605 {
606 std::vector<
607 std::vector<Tensor<2, dim, typename InputVector::value_type>>>
608 values(
609 1,
611 fe_values.get_fe().n_components()));
612 fe_values.get_function_hessians(solution, values);
613 return ProjectedDerivative(values[0][component]);
614 }
615 }
616
617
618
619 template <>
620 inline double
622 {
623 return std::fabs(d[0][0][0]);
624 }
625
626
627
628 template <int dim>
629 inline double
631 {
632 // return the Frobenius-norm. this is a
633 // member function of Tensor<rank_,dim>
634 return d.norm();
635 }
636
637
638 template <int dim>
639 inline void
641 {
642 // symmetrize non-diagonal entries
643
644 // first do it in the case, that i,j,k are
645 // pairwise different (which can onlky happen
646 // in dim >= 3)
647 for (unsigned int i = 0; i < dim; ++i)
648 for (unsigned int j = i + 1; j < dim; ++j)
649 for (unsigned int k = j + 1; k < dim; ++k)
650 {
651 const double s = (d[i][j][k] + d[i][k][j] + d[j][i][k] +
652 d[j][k][i] + d[k][i][j] + d[k][j][i]) /
653 6;
654 d[i][j][k] = d[i][k][j] = d[j][i][k] = d[j][k][i] = d[k][i][j] =
655 d[k][j][i] = s;
656 }
657 // now do the case, where two indices are
658 // equal
659 for (unsigned int i = 0; i < dim; ++i)
660 for (unsigned int j = i + 1; j < dim; ++j)
661 {
662 // case 1: index i (lower one) is
663 // double
664 const double s = (d[i][i][j] + d[i][j][i] + d[j][i][i]) / 3;
665 d[i][i][j] = d[i][j][i] = d[j][i][i] = s;
666
667 // case 2: index j (higher one) is
668 // double
669 const double t = (d[i][j][j] + d[j][i][j] + d[j][j][i]) / 3;
670 d[i][j][j] = d[j][i][j] = d[j][j][i] = t;
671 }
672 }
673
674
675 template <int order, int dim>
677 {
678 public:
684 using DerivDescr = void;
685 };
686
687 template <int dim>
688 class DerivativeSelector<1, dim>
689 {
690 public:
692 };
693
694 template <int dim>
695 class DerivativeSelector<2, dim>
696 {
697 public:
699 };
700
701 template <int dim>
702 class DerivativeSelector<3, dim>
703 {
704 public:
706 };
707 } // namespace internal
708} // namespace DerivativeApproximation
709
710// Dummy structures and dummy function used for WorkStream
712{
713 namespace internal
714 {
715 namespace Assembler
716 {
717 struct Scratch
718 {
719 Scratch() = default;
720 };
721
722 struct CopyData
723 {
724 CopyData() = default;
725 };
726 } // namespace Assembler
727 } // namespace internal
728} // namespace DerivativeApproximation
729
730// --------------- now for the functions that do the actual work --------------
731
733{
734 namespace internal
735 {
740 template <class DerivativeDescription,
741 int dim,
742 class InputVector,
743 int spacedim>
744 void
747 const DoFHandler<dim, spacedim> &dof_handler,
748 const InputVector &solution,
749 const unsigned int component,
751 typename DerivativeDescription::Derivative &derivative)
752 {
753 const QMidpoint<dim> midpoint_rule;
754
755 // create collection objects from
756 // single quadratures, mappings,
757 // and finite elements. if we have
758 // an hp-DoFHandler,
759 // dof_handler.get_fe() returns a
760 // collection of which we do a
761 // shallow copy instead
762 const hp::QCollection<dim> q_collection(midpoint_rule);
763 const hp::FECollection<dim> &fe_collection =
764 dof_handler.get_fe_collection();
765 const hp::MappingCollection<dim> mapping_collection(mapping);
766
767 hp::FEValues<dim> x_fe_midpoint_value(
768 mapping_collection,
769 fe_collection,
770 q_collection,
771 DerivativeDescription::update_flags | update_quadrature_points);
772
773 // matrix Y=sum_i y_i y_i^T
775
776
777 // vector to hold iterators to all
778 // active neighbors of a cell
779 // reserve the maximal number of
780 // active neighbors
781 std::vector<typename DoFHandler<dim, spacedim>::active_cell_iterator>
782 active_neighbors;
783
784 active_neighbors.reserve(GeometryInfo<dim>::faces_per_cell *
786
787 // vector
788 // g=sum_i y_i (f(x+y_i)-f(x))/|y_i|
789 // or related type for higher
790 // derivatives
791 typename DerivativeDescription::Derivative projected_derivative;
792
793 // reinit FE values object...
794 x_fe_midpoint_value.reinit(cell);
795 const FEValues<dim> &fe_midpoint_value =
796 x_fe_midpoint_value.get_present_fe_values();
797
798 // ...and get the value of the
799 // projected derivative...
800 const typename DerivativeDescription::ProjectedDerivative
801 this_midpoint_value =
802 DerivativeDescription::get_projected_derivative(fe_midpoint_value,
803 solution,
804 component);
805 // ...and the place where it lives
806 // This needs to be a copy. If it was a reference, it would be changed
807 // after the next `reinit` call of the FEValues object. clang-tidy
808 // complains about this not being a reference, so we suppress the warning.
809 const Point<dim> this_center =
810 fe_midpoint_value.quadrature_point(0); // NOLINT
811
812 // loop over all neighbors and
813 // accumulate the difference
814 // quotients from them. note
815 // that things get a bit more
816 // complicated if the neighbor
817 // is more refined than the
818 // present one
819 //
820 // to make processing simpler,
821 // first collect all neighbor
822 // cells in a vector, and then
823 // collect the data from them
825 cell, active_neighbors);
826
827 // now loop over all active
828 // neighbors and collect the
829 // data we need
830 auto neighbor_ptr = active_neighbors.begin();
831 for (; neighbor_ptr != active_neighbors.end(); ++neighbor_ptr)
832 {
833 const auto &neighbor = *neighbor_ptr;
834
835 // reinit FE values object...
836 x_fe_midpoint_value.reinit(neighbor);
837 const FEValues<dim> &neighbor_fe_midpoint_value =
838 x_fe_midpoint_value.get_present_fe_values();
839
840 // ...and get the value of the
841 // solution...
842 const typename DerivativeDescription::ProjectedDerivative
843 neighbor_midpoint_value =
844 DerivativeDescription::get_projected_derivative(
845 neighbor_fe_midpoint_value, solution, component);
846
847 // ...and the place where it lives
848 const Point<dim> &neighbor_center =
849 neighbor_fe_midpoint_value.quadrature_point(0);
850
851
852 // vector for the
853 // normalized
854 // direction between
855 // the centers of two
856 // cells
857
858 Tensor<1, dim> y = neighbor_center - this_center;
859 const double distance = y.norm();
860 // normalize y
861 y /= distance;
862 // *** note that unlike in
863 // the docs, y denotes the
864 // normalized vector
865 // connecting the centers
866 // of the two cells, rather
867 // than the normal
868 // difference! ***
869
870 // add up the
871 // contribution of
872 // this cell to Y
873 for (unsigned int i = 0; i < dim; ++i)
874 for (unsigned int j = 0; j < dim; ++j)
875 Y[i][j] += y[i] * y[j];
876
877 // then update the sum
878 // of difference
879 // quotients
880 typename DerivativeDescription::ProjectedDerivative
881 projected_finite_difference =
882 (neighbor_midpoint_value - this_midpoint_value);
883 projected_finite_difference /= distance;
884
885 projected_derivative += outer_product(y, projected_finite_difference);
886 }
887
888 // can we determine an
889 // approximation of the
890 // gradient for the present
891 // cell? if so, then we need to
892 // have passed over vectors y_i
893 // which span the whole space,
894 // otherwise we would not have
895 // all components of the
896 // gradient
898
899 // compute Y^-1 g
900 const Tensor<2, dim> Y_inverse = invert(Y);
901
902 derivative = Y_inverse * projected_derivative;
903
904 // finally symmetrize the derivative
905 DerivativeDescription::symmetrize(derivative);
906 }
907
908
909
915 template <class DerivativeDescription,
916 int dim,
917 class InputVector,
918 int spacedim>
919 void
925 const DoFHandler<dim, spacedim> &dof_handler,
926 const InputVector &solution,
927 const unsigned int component)
928 {
929 // if the cell is not locally owned, then there is nothing to do
930 if (std::get<0>(*cell)->is_locally_owned() == false)
931 *std::get<1>(*cell) = 0;
932 else
933 {
934 typename DerivativeDescription::Derivative derivative;
935 // call the function doing the actual
936 // work on this cell
938 mapping,
939 dof_handler,
940 solution,
941 component,
942 std::get<0>(*cell),
943 derivative);
944
945 // evaluate the norm and fill the vector
946 //*derivative_norm_on_this_cell
947 *std::get<1>(*cell) =
948 DerivativeDescription::derivative_norm(derivative);
949 }
950 }
951
952
963 template <class DerivativeDescription,
964 int dim,
965 class InputVector,
966 int spacedim>
967 void
969 const DoFHandler<dim, spacedim> &dof_handler,
970 const InputVector &solution,
971 const unsigned int component,
973 {
974 Assert(derivative_norm.size() ==
975 dof_handler.get_triangulation().n_active_cells(),
977 derivative_norm.size(),
978 dof_handler.get_triangulation().n_active_cells()));
979 AssertIndexRange(component, dof_handler.get_fe(0).n_components());
980
981 using Iterators =
982 std::tuple<typename DoFHandler<dim, spacedim>::active_cell_iterator,
985 Iterators(dof_handler.begin_active(), derivative_norm.begin())),
986 end(Iterators(dof_handler.end(), derivative_norm.end()));
987
988 // There is no need for a copier because there is no conflict between
989 // threads to write in derivative_norm. Scratch and CopyData are also
990 // useless.
992 begin,
993 end,
994 [&mapping, &dof_handler, &solution, component](
996 const Assembler::Scratch &,
999 cell, mapping, dof_handler, solution, component);
1000 },
1001 std::function<void(const internal::Assembler::CopyData &)>(),
1004 }
1005
1006 } // namespace internal
1007
1008} // namespace DerivativeApproximation
1009
1010
1011// ------------------------ finally for the public interface of this namespace
1012
1014{
1015 template <int dim, class InputVector, int spacedim>
1016 void
1018 const DoFHandler<dim, spacedim> &dof_handler,
1019 const InputVector &solution,
1021 const unsigned int component)
1022 {
1024 mapping, dof_handler, solution, component, derivative_norm);
1025 }
1026
1027
1028 template <int dim, class InputVector, int spacedim>
1029 void
1031 const InputVector &solution,
1033 const unsigned int component)
1034 {
1035 Assert(!dof_handler.get_triangulation().is_mixed_mesh(),
1037 const auto reference_cell =
1038 dof_handler.get_triangulation().get_reference_cells()[0];
1040 reference_cell.template get_default_linear_mapping<dim, spacedim>(),
1041 dof_handler,
1042 solution,
1043 component,
1045 }
1046
1047
1048 template <int dim, class InputVector, int spacedim>
1049 void
1051 const DoFHandler<dim, spacedim> &dof_handler,
1052 const InputVector &solution,
1054 const unsigned int component)
1055 {
1057 mapping, dof_handler, solution, component, derivative_norm);
1058 }
1059
1060
1061 template <int dim, class InputVector, int spacedim>
1062 void
1064 const InputVector &solution,
1066 const unsigned int component)
1067 {
1068 Assert(!dof_handler.get_triangulation().is_mixed_mesh(),
1070 const auto reference_cell =
1071 dof_handler.get_triangulation().get_reference_cells()[0];
1073 reference_cell.template get_default_linear_mapping<dim, spacedim>(),
1074 dof_handler,
1075 solution,
1076 component,
1078 }
1079
1080
1081 template <int dim, int spacedim, class InputVector, int order>
1082 void
1085 const DoFHandler<dim, spacedim> &dof,
1086 const InputVector &solution,
1087#ifndef _MSC_VER
1089#else
1091 &cell,
1092#endif
1093 Tensor<order, dim> &derivative,
1094 const unsigned int component)
1095 {
1098 mapping, dof, solution, component, cell, derivative);
1099 }
1100
1101
1102
1103 template <int dim, int spacedim, class InputVector, int order>
1104 void
1106 const DoFHandler<dim, spacedim> &dof,
1107 const InputVector &solution,
1108#ifndef _MSC_VER
1110#else
1112 &cell,
1113#endif
1114 Tensor<order, dim> &derivative,
1115 const unsigned int component)
1116 {
1117 // just call the respective function with Q1 mapping
1119 cell->reference_cell()
1121 dof,
1122 solution,
1123 cell,
1124 derivative,
1125 component);
1126 }
1127
1128
1129
1130 template <int dim, int order>
1131 double
1137
1138} // namespace DerivativeApproximation
1139
1140
1141// --------------------------- explicit instantiations ---------------------
1142#include "numerics/derivative_approximation.inst"
1143
1144
static ProjectedDerivative get_projected_derivative(const FEValues< dim, spacedim > &fe_values, const InputVector &solution, const unsigned int component)
static void symmetrize(Derivative &derivative_tensor)
static double derivative_norm(const Derivative &d)
static ProjectedDerivative get_projected_derivative(const FEValues< dim, spacedim > &fe_values, const InputVector &solution, const unsigned int component)
static void symmetrize(Derivative &derivative_tensor)
static ProjectedDerivative get_projected_derivative(const FEValues< dim, spacedim > &fe_values, const InputVector &solution, const unsigned int component)
static void symmetrize(Derivative &derivative_tensor)
cell_iterator end() const
const hp::FECollection< dim, spacedim > & get_fe_collection() const
const FiniteElement< dim, spacedim > & get_fe(const types::fe_index index=0) const
const Triangulation< dim, spacedim > & get_triangulation() const
active_cell_iterator begin_active(const unsigned int level=0) const
void get_function_values(const ReadVector< Number > &fe_function, std::vector< Number > &values) const
void get_function_hessians(const ReadVector< Number > &fe_function, std::vector< Tensor< 2, spacedim, Number > > &hessians) const
const Point< spacedim > & quadrature_point(const unsigned int q_point) const
void get_function_gradients(const ReadVector< Number > &fe_function, std::vector< Tensor< 1, spacedim, Number > > &gradients) const
const FiniteElement< dim, spacedim > & get_fe() const
unsigned int n_components() const
Abstract base class for mapping classes.
Definition mapping.h:320
Definition point.h:113
numbers::NumberTraits< Number >::real_type norm() const
unsigned int n_active_cells() const
const std::vector< ReferenceCell > & get_reference_cells() const
bool is_mixed_mesh() const
value_type * iterator
Definition vector.h:140
const FEValuesType & get_present_fe_values() const
Definition fe_values.h:695
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, lda > > &cell, const unsigned int q_index=numbers::invalid_unsigned_int, const unsigned int mapping_index=numbers::invalid_unsigned_int, const unsigned int fe_index=numbers::invalid_unsigned_int)
Definition fe_values.cc:296
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:40
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:41
#define DEAL_II_NOT_IMPLEMENTED()
static ::ExceptionBase & ExcInsufficientDirections()
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcVectorLengthVsNActiveCells(int arg1, int arg2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
#define AssertThrow(cond, exc)
typename ActiveSelector::active_cell_iterator active_cell_iterator
@ update_hessians
Second derivatives of shape functions.
@ update_values
Shape function values.
@ update_gradients
Shape function gradients.
@ update_quadrature_points
Transformed quadrature points.
const Mapping< dim, spacedim > & get_default_linear_mapping(const Triangulation< dim, spacedim > &triangulation)
Definition mapping.cc:316
MappingQ< dim, spacedim > StaticMappingQ1< dim, spacedim >::mapping
Definition mapping_q1.h:104
void approximate_derivative(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof_handler, const InputVector &solution, const unsigned int component, Vector< float > &derivative_norm)
void approximate(const SynchronousIterators< std::tuple< typename DoFHandler< dim, spacedim >::active_cell_iterator, Vector< float >::iterator > > &cell, const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof_handler, const InputVector &solution, const unsigned int component)
void approximate_cell(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof_handler, const InputVector &solution, const unsigned int component, const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell, typename DerivativeDescription::Derivative &derivative)
void approximate_gradient(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const InputVector &solution, Vector< float > &derivative_norm, const unsigned int component=0)
double derivative_norm(const Tensor< order, dim > &derivative)
void approximate_derivative_tensor(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const InputVector &solution, const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell, Tensor< order, dim > &derivative, const unsigned int component=0)
void approximate_second_derivative(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const InputVector &solution, Vector< float > &derivative_norm, const unsigned int component=0)
void get_active_neighbors(const typename MeshType::active_cell_iterator &cell, std::vector< typename MeshType::active_cell_iterator > &active_neighbors)
constexpr char T
constexpr T fixed_power(const T t)
Definition utilities.h:943
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)
constexpr double PI
Definition numbers.h:239
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > cos(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
inline ::VectorizedArray< Number, width > acos(const ::VectorizedArray< Number, width > &x)
static constexpr unsigned int max_children_per_face
static constexpr unsigned int faces_per_cell
DEAL_II_HOST constexpr Number determinant(const SymmetricTensor< 2, dim, Number > &)
DEAL_II_HOST constexpr SymmetricTensor< 2, dim, Number > invert(const SymmetricTensor< 2, dim, Number > &)
std::array< Number, 1 > eigenvalues(const SymmetricTensor< 2, 1, Number > &T)
DEAL_II_HOST constexpr SymmetricTensor< 4, dim, Number > outer_product(const SymmetricTensor< 2, dim, Number > &t1, const SymmetricTensor< 2, dim, Number > &t2)
DEAL_II_HOST constexpr Number trace(const SymmetricTensor< 2, dim2, Number > &)