16 #include <deal.II/base/quadrature_lib.h>
17 #include <deal.II/base/work_stream.h>
19 #include <deal.II/dofs/dof_accessor.h>
20 #include <deal.II/dofs/dof_handler.h>
22 #include <deal.II/fe/fe.h>
23 #include <deal.II/fe/fe_values.h>
24 #include <deal.II/fe/mapping_q1.h>
26 #include <deal.II/grid/filtered_iterator.h>
27 #include <deal.II/grid/grid_tools.h>
28 #include <deal.II/grid/tria_iterator.h>
30 #include <deal.II/hp/fe_collection.h>
31 #include <deal.II/hp/fe_values.h>
32 #include <deal.II/hp/mapping_collection.h>
33 #include <deal.II/hp/q_collection.h>
35 #include <deal.II/lac/block_vector.h>
36 #include <deal.II/lac/la_parallel_block_vector.h>
37 #include <deal.II/lac/la_parallel_vector.h>
38 #include <deal.II/lac/la_vector.h>
39 #include <deal.II/lac/petsc_block_vector.h>
40 #include <deal.II/lac/petsc_vector.h>
41 #include <deal.II/lac/trilinos_parallel_block_vector.h>
42 #include <deal.II/lac/trilinos_tpetra_vector.h>
43 #include <deal.II/lac/trilinos_vector.h>
44 #include <deal.II/lac/vector.h>
46 #include <deal.II/numerics/derivative_approximation.h>
50 DEAL_II_NAMESPACE_OPEN
106 template <
class InputVector,
int spacedim>
109 const InputVector & solution,
110 const unsigned int component);
136 template <
class InputVector,
int spacedim>
140 const InputVector & solution,
141 const unsigned int component)
143 if (fe_values.
get_fe().n_components() == 1)
145 std::vector<typename InputVector::value_type> values(1);
151 std::vector<Vector<typename InputVector::value_type>> values(
153 Vector<typename InputVector::value_type>(
154 fe_values.
get_fe().n_components()));
156 return values[0](component);
167 for (
unsigned int i = 0; i < dim; ++i)
219 template <
class InputVector,
int spacedim>
222 const InputVector & solution,
223 const unsigned int component);
253 template <
class InputVector,
int spacedim>
257 const InputVector & solution,
258 const unsigned int component)
260 if (fe_values.
get_fe().n_components() == 1)
262 std::vector<Tensor<1, dim, typename InputVector::value_type>> values(
270 std::vector<Tensor<1, dim, typename InputVector::value_type>>>
274 fe_values.
get_fe().n_components()));
286 return std::fabs(d[0][0]);
303 const double radicand =
304 ::sqr(d[0][0] - d[1][1]) + 4 * ::sqr(d[0][1]);
306 0.5 * (d[0][0] + d[1][1] + std::sqrt(radicand)),
307 0.5 * (d[0][0] + d[1][1] - std::sqrt(radicand))};
427 const double am =
trace(d) / 3.;
431 for (
unsigned int i = 0; i < 3; ++i)
434 const double ss01 = s[0][1] * s[0][1], ss12 = s[1][2] * s[1][2],
435 ss02 = s[0][2] * s[0][2];
437 const double J2 = (s[0][0] * s[0][0] + s[1][1] * s[1][1] +
438 s[2][2] * s[2][2] + 2 * (ss01 + ss02 + ss12)) /
441 (std::pow(s[0][0], 3) + std::pow(s[1][1], 3) + std::pow(s[2][2], 3) +
442 3. * s[0][0] * (ss01 + ss02) + 3. * s[1][1] * (ss01 + ss12) +
443 3. * s[2][2] * (ss02 + ss12) + 6. * s[0][1] * s[0][2] * s[1][2]) /
446 const double R = std::sqrt(4. * J2 / 3.);
448 double EE[3] = {0, 0, 0};
455 if (R <= 1e-14 * std::fabs(am))
456 EE[0] = EE[1] = EE[2] = am;
461 const double R3 = R * R * R;
462 const double XX = 4. * J3 / R3;
463 const double YY = 1. - std::fabs(XX);
470 const double a = (XX > 0 ? -1. : 1.) * R / 2;
471 EE[0] = EE[1] = am + a;
476 const double theta = std::acos(XX) / 3.;
477 EE[0] = am + R * std::cos(theta);
478 EE[1] = am + R * std::cos(theta + 2. / 3. *
numbers::PI);
479 EE[2] = am + R * std::cos(theta + 4. / 3. *
numbers::PI);
483 return std::max(std::fabs(EE[0]),
484 std::max(std::fabs(EE[1]), std::fabs(EE[2])));
517 for (
unsigned int i = 0; i < dim; ++i)
518 for (
unsigned int j = i + 1; j < dim; ++j)
520 const double s = (d[i][j] + d[j][i]) / 2;
521 d[i][j] = d[j][i] = s;
528 class ThirdDerivative
535 static const UpdateFlags update_flags;
556 template <
class InputVector,
int spacedim>
557 static ProjectedDerivative
559 const InputVector & solution,
560 const unsigned int component);
586 const UpdateFlags ThirdDerivative<dim>::update_flags = update_hessians;
590 template <
class InputVector,
int spacedim>
591 inline typename ThirdDerivative<dim>::ProjectedDerivative
592 ThirdDerivative<dim>::get_projected_derivative(
594 const InputVector & solution,
595 const unsigned int component)
597 if (fe_values.
get_fe().n_components() == 1)
599 std::vector<Tensor<2, dim, typename InputVector::value_type>> values(
602 return ProjectedDerivative(values[0]);
607 std::vector<Tensor<2, dim, typename InputVector::value_type>>>
611 fe_values.
get_fe().n_components()));
613 return ProjectedDerivative(values[0][component]);
621 ThirdDerivative<1>::derivative_norm(
const Derivative &d)
623 return std::fabs(d[0][0][0]);
630 ThirdDerivative<dim>::derivative_norm(
const Derivative &d)
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)
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]) /
654 d[i][j][k] =
d[i][k][j] =
d[j][i][k] =
d[j][k][i] =
d[k][i][j] =
659 for (
unsigned int i = 0; i < dim; ++i)
660 for (
unsigned int j = i + 1; j < dim; ++j)
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;
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;
675 template <
int order,
int dim>
676 class DerivativeSelector
684 using DerivDescr = void;
688 class DerivativeSelector<1, dim>
691 using DerivDescr = Gradient<dim>;
695 class DerivativeSelector<2, dim>
698 using DerivDescr = SecondDerivative<dim>;
702 class DerivativeSelector<3, dim>
705 using DerivDescr = ThirdDerivative<dim>;
724 CopyData() =
default;
740 template <
class DerivativeDescription,
742 template <
int,
int>
class DoFHandlerType,
748 const DoFHandlerType<dim, spacedim> &dof_handler,
749 const InputVector & solution,
750 const unsigned int component,
753 typename DerivativeDescription::Derivative &derivative)
766 dof_handler.get_fe_collection();
773 DerivativeDescription::update_flags | update_quadrature_points);
794 typename DerivativeDescription::Derivative projected_derivative;
797 x_fe_midpoint_value.reinit(cell);
803 const typename DerivativeDescription::ProjectedDerivative
804 this_midpoint_value =
805 DerivativeDescription::get_projected_derivative(fe_midpoint_value,
823 GridTools::get_active_neighbors<DoFHandlerType<dim, spacedim>>(
824 cell, active_neighbors);
831 const_iterator neighbor_ptr = active_neighbors.begin();
832 for (; neighbor_ptr != active_neighbors.end(); ++neighbor_ptr)
836 neighbor = *neighbor_ptr;
839 x_fe_midpoint_value.reinit(neighbor);
845 const typename DerivativeDescription::ProjectedDerivative
846 neighbor_midpoint_value =
847 DerivativeDescription::get_projected_derivative(
848 neighbor_fe_midpoint_value, solution, component);
861 const double distance = y.
norm();
875 for (
unsigned int i = 0; i < dim; ++i)
876 for (
unsigned int j = 0; j < dim; ++j)
877 Y[i][j] += y[i] * y[j];
882 typename DerivativeDescription::ProjectedDerivative
883 projected_finite_difference =
884 (neighbor_midpoint_value - this_midpoint_value);
885 projected_finite_difference /= distance;
887 projected_derivative +=
outer_product(y, projected_finite_difference);
904 derivative = Y_inverse * projected_derivative;
917 template <
class DerivativeDescription,
919 template <
int,
int>
class DoFHandlerType,
927 Vector<float>::iterator>>
const & cell,
929 const DoFHandlerType<dim, spacedim> &dof_handler,
930 const InputVector & solution,
931 const unsigned int component)
934 if (std::get<0>(*cell)->is_locally_owned() ==
false)
935 *std::get<1>(*cell) = 0;
938 typename DerivativeDescription::Derivative derivative;
941 approximate_cell<DerivativeDescription,
954 *std::get<1>(*cell) =
955 DerivativeDescription::derivative_norm(derivative);
970 template <
class DerivativeDescription,
972 template <
int,
int>
class DoFHandlerType,
977 const DoFHandlerType<dim, spacedim> &dof_handler,
978 const InputVector & solution,
979 const unsigned int component,
983 dof_handler.get_triangulation().n_active_cells(),
986 dof_handler.get_triangulation().n_active_cells()));
987 Assert(component < dof_handler.get_fe(0).n_components(),
988 ExcIndexRange(component, 0, dof_handler.get_fe(0).n_components()));
990 using Iterators = std::tuple<
993 Vector<float>::iterator>;
1005 Assembler::Scratch
const &,
1006 Assembler::CopyData &)
>>(
1007 std::bind(&approximate<DerivativeDescription,
1012 std::placeholders::_1,
1014 std::cref(dof_handler),
1015 std::cref(solution),
1017 std::function<
void(internal::Assembler::CopyData
const &)>(),
1018 internal::Assembler::Scratch(),
1019 internal::Assembler::CopyData());
1032 template <
int,
int>
class DoFHandlerType,
1037 const DoFHandlerType<dim, spacedim> &dof_handler,
1038 const InputVector & solution,
1040 const unsigned int component)
1042 internal::approximate_derivative<internal::Gradient<dim>, dim>(
1048 template <
int,
int>
class DoFHandlerType,
1053 const InputVector & solution,
1055 const unsigned int component)
1057 internal::approximate_derivative<internal::Gradient<dim>, dim>(
1067 template <
int,
int>
class DoFHandlerType,
1073 const DoFHandlerType<dim, spacedim> &dof_handler,
1074 const InputVector & solution,
1076 const unsigned int component)
1078 internal::approximate_derivative<internal::SecondDerivative<dim>, dim>(
1084 template <
int,
int>
class DoFHandlerType,
1089 const DoFHandlerType<dim, spacedim> &dof_handler,
1090 const InputVector & solution,
1092 const unsigned int component)
1094 internal::approximate_derivative<internal::SecondDerivative<dim>, dim>(
1103 template <
typename DoFHandlerType,
class InputVector,
int order>
1108 const DoFHandlerType &dof,
1109 const InputVector & solution,
1111 const typename DoFHandlerType::active_cell_iterator &cell,
1117 const unsigned int component)
1119 internal::approximate_cell<
1120 typename internal::DerivativeSelector<order, DoFHandlerType::dimension>::
1121 DerivDescr>(mapping, dof, solution, component, cell, derivative);
1126 template <
typename DoFHandlerType,
class InputVector,
int order>
1129 const DoFHandlerType &dof,
1130 const InputVector & solution,
1132 const typename DoFHandlerType::active_cell_iterator &cell,
1138 const unsigned int component)
1143 DoFHandlerType::space_dimension>::mapping,
1153 template <
int dim,
int order>
1157 return internal::DerivativeSelector<order, dim>::DerivDescr::
1158 derivative_norm(derivative);
1165 #include "derivative_approximation.inst"
1168 DEAL_II_NAMESPACE_CLOSE