107 template <
class InputVector,
int spacedim>
110 const InputVector & solution,
111 const unsigned int component);
137 template <
class InputVector,
int spacedim>
141 const InputVector & solution,
142 const unsigned int component)
144 if (fe_values.get_fe().n_components() == 1)
146 std::vector<typename InputVector::value_type> values(1);
147 fe_values.get_function_values(solution, values);
152 std::vector<Vector<typename InputVector::value_type>> values(
154 Vector<typename InputVector::value_type>(
155 fe_values.get_fe().n_components()));
156 fe_values.get_function_values(solution, values);
157 return values[0](component);
168 for (
unsigned int i = 0; i < dim; ++i)
220 template <
class InputVector,
int spacedim>
223 const InputVector & solution,
224 const unsigned int component);
254 template <
class InputVector,
int spacedim>
258 const InputVector & solution,
259 const unsigned int component)
261 if (fe_values.get_fe().n_components() == 1)
263 std::vector<Tensor<1, dim, typename InputVector::value_type>> values(
265 fe_values.get_function_gradients(solution, values);
271 std::vector<Tensor<1, dim, typename InputVector::value_type>>>
275 fe_values.get_fe().n_components()));
276 fe_values.get_function_gradients(solution, values);
304 const double radicand =
305 ::sqr(d[0][0] - d[1][1]) + 4 * ::sqr(d[0][1]);
307 0.5 * (d[0][0] + d[1][1] + std::sqrt(radicand)),
308 0.5 * (d[0][0] + d[1][1] - std::sqrt(radicand))};
428 const double am =
trace(d) / 3.;
432 for (
unsigned int i = 0; i < 3; ++i)
435 const double ss01 = s[0][1] * s[0][1], ss12 = s[1][2] * s[1][2],
436 ss02 = s[0][2] * s[0][2];
438 const double J2 = (s[0][0] * s[0][0] + s[1][1] * s[1][1] +
439 s[2][2] * s[2][2] + 2 * (ss01 + ss02 + ss12)) /
442 (std::pow(s[0][0], 3) + std::pow(s[1][1], 3) + std::pow(s[2][2], 3) +
443 3. * s[0][0] * (ss01 + ss02) + 3. * s[1][1] * (ss01 + ss12) +
444 3. * s[2][2] * (ss02 + ss12) + 6. * s[0][1] * s[0][2] * s[1][2]) /
447 const double R = std::sqrt(4. * J2 / 3.);
449 double EE[3] = {0, 0, 0};
457 EE[0] = EE[1] = EE[2] = am;
462 const double R3 = R * R * R;
463 const double XX = 4. * J3 / R3;
471 const double a = (XX > 0 ? -1. : 1.) * R / 2;
472 EE[0] = EE[1] = am + a;
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);
518 for (
unsigned int i = 0; i < dim; ++i)
519 for (
unsigned int j = i + 1; j < dim; ++j)
521 const double s = (d[i][j] + d[j][i]) / 2;
522 d[i][j] = d[j][i] = s;
557 template <
class InputVector,
int spacedim>
560 const InputVector & solution,
561 const unsigned int component);
591 template <
class InputVector,
int spacedim>
595 const InputVector & solution,
596 const unsigned int component)
598 if (fe_values.get_fe().n_components() == 1)
600 std::vector<Tensor<2, dim, typename InputVector::value_type>> values(
602 fe_values.get_function_hessians(solution, values);
608 std::vector<Tensor<2, dim, typename InputVector::value_type>>>
612 fe_values.get_fe().n_components()));
613 fe_values.get_function_hessians(solution, values);
648 for (
unsigned int i = 0; i < dim; ++i)
649 for (
unsigned int j = i + 1; j < dim; ++j)
650 for (
unsigned int k = j + 1; k < dim; ++k)
652 const double s = (d[i][j][k] + d[i][k][j] + d[j][i][k] +
653 d[j][k][i] + d[k][i][j] + d[k][j][i]) /
655 d[i][j][k] = d[i][k][j] = d[j][i][k] = d[j][k][i] = d[k][i][j] =
660 for (
unsigned int i = 0; i < dim; ++i)
661 for (
unsigned int j = i + 1; j < dim; ++j)
665 const double s = (d[i][i][j] + d[i][j][i] + d[j][i][i]) / 3;
666 d[i][i][j] = d[i][j][i] = d[j][i][i] = s;
670 const double t = (d[i][j][j] + d[j][i][j] + d[j][j][i]) / 3;
671 d[i][j][j] = d[j][i][j] = d[j][j][i] = t;
676 template <
int order,
int dim>
741 template <
class DerivativeDescription,
743 template <
int,
int>
class DoFHandlerType,
749 const DoFHandlerType<dim, spacedim> &dof_handler,
750 const InputVector & solution,
751 const unsigned int component,
754 typename DerivativeDescription::Derivative &derivative)
767 dof_handler.get_fe_collection();
795 typename DerivativeDescription::Derivative projected_derivative;
798 x_fe_midpoint_value.
reinit(cell);
804 const typename DerivativeDescription::ProjectedDerivative
805 this_midpoint_value =
806 DerivativeDescription::get_projected_derivative(fe_midpoint_value,
824 GridTools::get_active_neighbors<DoFHandlerType<dim, spacedim>>(
825 cell, active_neighbors);
831 ::DoFCellAccessor<DoFHandlerType<dim, spacedim>,
false>>>::
832 const_iterator neighbor_ptr = active_neighbors.begin();
833 for (; neighbor_ptr != active_neighbors.end(); ++neighbor_ptr)
836 ::DoFCellAccessor<DoFHandlerType<dim, spacedim>,
false>>
837 neighbor = *neighbor_ptr;
840 x_fe_midpoint_value.
reinit(neighbor);
846 const typename DerivativeDescription::ProjectedDerivative
847 neighbor_midpoint_value =
848 DerivativeDescription::get_projected_derivative(
849 neighbor_fe_midpoint_value, solution, component);
862 const double distance = y.
norm();
876 for (
unsigned int i = 0; i < dim; ++i)
877 for (
unsigned int j = 0; j < dim; ++j)
878 Y[i][j] += y[i] * y[j];
883 typename DerivativeDescription::ProjectedDerivative
884 projected_finite_difference =
885 (neighbor_midpoint_value - this_midpoint_value);
886 projected_finite_difference /= distance;
888 projected_derivative +=
outer_product(y, projected_finite_difference);
905 derivative = Y_inverse * projected_derivative;
918 template <
class DerivativeDescription,
920 template <
int,
int>
class DoFHandlerType,
930 const DoFHandlerType<dim, spacedim> &dof_handler,
931 const InputVector & solution,
932 const unsigned int component)
935 if (std::get<0>(*cell)->is_locally_owned() ==
false)
936 *std::get<1>(*cell) = 0;
939 typename DerivativeDescription::Derivative derivative;
955 *std::get<1>(*cell) =
971 template <
class DerivativeDescription,
973 template <
int,
int>
class DoFHandlerType,
978 const DoFHandlerType<dim, spacedim> &dof_handler,
979 const InputVector & solution,
980 const unsigned int component,
983 Assert(derivative_norm.size() ==
984 dof_handler.get_triangulation().n_active_cells(),
986 derivative_norm.size(),
987 dof_handler.get_triangulation().n_active_cells()));
990 using Iterators = std::tuple<
995 Iterators(dof_handler.begin_active(), derivative_norm.begin())),
996 end(Iterators(dof_handler.end(), derivative_norm.end()));
1004 [&mapping, &dof_handler, &solution, component](
1013 cell, mapping, dof_handler, solution, component);
1015 std::function<void(internal::Assembler::CopyData const &)>(),
1030 template <
int,
int>
class DoFHandlerType,
1035 const DoFHandlerType<dim, spacedim> &dof_handler,
1036 const InputVector & solution,
1038 const unsigned int component)
1040 internal::approximate_derivative<internal::Gradient<dim>, dim>(
1046 template <
int,
int>
class DoFHandlerType,
1051 const InputVector & solution,
1053 const unsigned int component)
1055 internal::approximate_derivative<internal::Gradient<dim>, dim>(
1065 template <
int,
int>
class DoFHandlerType,
1071 const DoFHandlerType<dim, spacedim> &dof_handler,
1072 const InputVector & solution,
1074 const unsigned int component)
1076 internal::approximate_derivative<internal::SecondDerivative<dim>, dim>(
1082 template <
int,
int>
class DoFHandlerType,
1087 const DoFHandlerType<dim, spacedim> &dof_handler,
1088 const InputVector & solution,
1090 const unsigned int component)
1092 internal::approximate_derivative<internal::SecondDerivative<dim>, dim>(
1101 template <
typename DoFHandlerType,
class InputVector,
int order>
1106 const DoFHandlerType &dof,
1107 const InputVector & solution,
1109 const typename DoFHandlerType::active_cell_iterator &cell,
1115 const unsigned int component)
1119 DerivDescr>(mapping, dof, solution, component, cell, derivative);
1124 template <
typename DoFHandlerType,
class InputVector,
int order>
1127 const DoFHandlerType &dof,
1128 const InputVector & solution,
1130 const typename DoFHandlerType::active_cell_iterator &cell,
1136 const unsigned int component)
1141 DoFHandlerType::space_dimension>::mapping,
1151 template <
int dim,
int order>
1163 #include "derivative_approximation.inst"
static ::ExceptionBase & ExcInsufficientDirections()
constexpr Number determinant(const SymmetricTensor< 2, dim, Number > &)
constexpr SymmetricTensor< 2, dim, Number > symmetrize(const Tensor< 2, dim, Number > &t)
static void symmetrize(Derivative &derivative_tensor)
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
void approximate_cell(const Mapping< dim, spacedim > &mapping, const DoFHandlerType< dim, spacedim > &dof_handler, const InputVector &solution, const unsigned int component, const TriaActiveIterator< ::DoFCellAccessor< DoFHandlerType< dim, spacedim >, false >> &cell, typename DerivativeDescription::Derivative &derivative)
static const UpdateFlags update_flags
#define AssertIndexRange(index, range)
Transformed quadrature points.
#define AssertThrow(cond, exc)
static ProjectedDerivative get_projected_derivative(const FEValues< dim, spacedim > &fe_values, const InputVector &solution, const unsigned int component)
constexpr SymmetricTensor< 4, dim, Number > outer_product(const SymmetricTensor< 2, dim, Number > &t1, const SymmetricTensor< 2, dim, Number > &t2)
const Point< spacedim > & quadrature_point(const unsigned int q) const
static void symmetrize(Derivative &derivative_tensor)
static const UpdateFlags update_flags
Expression acos(const Expression &x)
#define Assert(cond, exc)
static double derivative_norm(const Derivative &d)
Abstract base class for mapping classes.
#define DEAL_II_NAMESPACE_CLOSE
static void symmetrize(Derivative &derivative_tensor)
VectorType::value_type * end(VectorType &V)
static ProjectedDerivative get_projected_derivative(const FEValues< dim, spacedim > &fe_values, const InputVector &solution, const unsigned int component)
void reinit(const TriaIterator< DoFCellAccessor< DoFHandlerType, 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)
static double derivative_norm(const Derivative &d)
Expression fabs(const Expression &x)
void approximate_derivative(const Mapping< dim, spacedim > &mapping, const DoFHandlerType< dim, spacedim > &dof_handler, const InputVector &solution, const unsigned int component, Vector< float > &derivative_norm)
Second derivatives of shape functions.
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
double derivative_norm(const Tensor< order, dim > &derivative)
void approximate_derivative_tensor(const Mapping< DoFHandlerType::dimension, DoFHandlerType::space_dimension > &mapping, const DoFHandlerType &dof, const InputVector &solution, const typename DoFHandlerType::active_cell_iterator &cell, Tensor< order, DoFHandlerType::dimension > &derivative, const unsigned int component=0)
static constexpr double PI
#define DEAL_II_NAMESPACE_OPEN
VectorType::value_type * begin(VectorType &V)
Shape function gradients.
void approximate(SynchronousIterators< std::tuple< TriaActiveIterator< ::DoFCellAccessor< DoFHandlerType< dim, spacedim >, false >>, Vector< float >::iterator >> const &cell, const Mapping< dim, spacedim > &mapping, const DoFHandlerType< dim, spacedim > &dof_handler, const InputVector &solution, const unsigned int component)
constexpr Number trace(const SymmetricTensor< 2, dim2, Number > &)
static ::ExceptionBase & ExcNotImplemented()
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)
static ProjectedDerivative get_projected_derivative(const FEValues< dim, spacedim > &fe_values, const InputVector &solution, const unsigned int component)
static double derivative_norm(const Derivative &d)
constexpr SymmetricTensor< 2, dim, Number > invert(const SymmetricTensor< 2, dim, Number > &)
static const UpdateFlags update_flags
numbers::NumberTraits< Number >::real_type norm() const
void approximate_second_derivative(const Mapping< dim, spacedim > &mapping, const DoFHandlerType< dim, spacedim > &dof, const InputVector &solution, Vector< float > &derivative_norm, const unsigned int component=0)
void approximate_gradient(const Mapping< dim, spacedim > &mapping, const DoFHandlerType< dim, spacedim > &dof, const InputVector &solution, Vector< float > &derivative_norm, const unsigned int component=0)
Eigenvalue vector is filled.
Tensor< 0, dim > ProjectedDerivative
T max(const T &t, const MPI_Comm &mpi_communicator)
static ::ExceptionBase & ExcInternalError()
const ::FEValues< dim, spacedim > & get_present_fe_values() const
static ::ExceptionBase & ExcVectorLengthVsNActiveCells(int arg1, int arg2)