105 template <
class InputVector,
int spacedim>
108 const InputVector &solution,
109 const unsigned int component);
135 template <
class InputVector,
int spacedim>
139 const InputVector &solution,
140 const unsigned int component)
144 std::vector<typename InputVector::value_type> values(1);
150 std::vector<Vector<typename InputVector::value_type>> values(
155 return values[0](component);
166 for (
unsigned int i = 0; i < dim; ++i)
216 template <
class InputVector,
int spacedim>
219 const InputVector &solution,
220 const unsigned int component);
250 template <
class InputVector,
int spacedim>
254 const InputVector &solution,
255 const unsigned int component)
259 std::vector<Tensor<1, dim, typename InputVector::value_type>> values(
267 std::vector<Tensor<1, dim, typename InputVector::value_type>>>
283 return std::fabs(d[0][0]);
300 const double radicand =
301 ::sqr(d[0][0] - d[1][1]) + 4 * ::sqr(d[0][1]);
303 0.5 * (d[0][0] + d[1][1] +
std::sqrt(radicand)),
304 0.5 * (d[0][0] + d[1][1] -
std::sqrt(radicand))};
424 const double am =
trace(d) / 3.;
428 for (
unsigned int i = 0; i < 3; ++i)
431 const double ss01 = s[0][1] * s[0][1], ss12 = s[1][2] * s[1][2],
432 ss02 = s[0][2] * s[0][2];
434 const double J2 = (s[0][0] * s[0][0] + s[1][1] * s[1][1] +
435 s[2][2] * s[2][2] + 2 * (ss01 + ss02 + ss12)) /
441 3. * s[1][1] * (ss01 + ss12) + 3. * s[2][2] * (ss02 + ss12) +
442 6. * s[0][1] * s[0][2] * s[1][2]) /
445 const double R =
std::sqrt(4. * J2 / 3.);
447 double EE[3] = {0, 0, 0};
454 if (R <= 1e-14 * std::fabs(am))
455 EE[0] = EE[1] = EE[2] = am;
460 const double R3 = R * R * R;
461 const double XX = 4. * J3 / R3;
462 const double YY = 1. - std::fabs(XX);
469 const double a = (XX > 0 ? -1. : 1.) * R / 2;
470 EE[0] = EE[1] = am + a;
483 std::max(std::fabs(EE[1]), std::fabs(EE[2])));
516 for (
unsigned int i = 0; i < dim; ++i)
517 for (
unsigned int j = i + 1; j < dim; ++j)
519 const double s = (d[i][j] + d[j][i]) / 2;
520 d[i][j] = d[j][i] = s;
555 template <
class InputVector,
int spacedim>
558 const InputVector &solution,
559 const unsigned int component);
589 template <
class InputVector,
int spacedim>
593 const InputVector &solution,
594 const unsigned int component)
598 std::vector<Tensor<2, dim, typename InputVector::value_type>> values(
606 std::vector<Tensor<2, dim, typename InputVector::value_type>>>
622 return std::fabs(d[0][0][0]);
646 for (
unsigned int i = 0; i < dim; ++i)
647 for (
unsigned int j = i + 1; j < dim; ++j)
648 for (
unsigned int k = j + 1; k < dim; ++k)
650 const double s = (d[i][j][k] + d[i][k][j] + d[j][i][k] +
651 d[j][k][i] + d[k][i][j] + d[k][j][i]) /
653 d[i][j][k] = d[i][k][j] = d[j][i][k] = d[j][k][i] = d[k][i][j] =
658 for (
unsigned int i = 0; i < dim; ++i)
659 for (
unsigned int j = i + 1; j < dim; ++j)
663 const double s = (d[i][i][j] + d[i][j][i] + d[j][i][i]) / 3;
664 d[i][i][j] = d[i][j][i] = d[j][i][i] = s;
668 const double t = (d[i][j][j] + d[j][i][j] + d[j][j][i]) / 3;
669 d[i][j][j] = d[j][i][j] = d[j][j][i] = t;
674 template <
int order,
int dim>
739 template <
class DerivativeDescription,
747 const InputVector &solution,
748 const unsigned int component,
750 typename DerivativeDescription::Derivative &derivative)
780 std::vector<typename DoFHandler<dim, spacedim>::active_cell_iterator>
790 typename DerivativeDescription::Derivative projected_derivative;
793 x_fe_midpoint_value.
reinit(cell);
799 const typename DerivativeDescription::ProjectedDerivative
800 this_midpoint_value =
801 DerivativeDescription::get_projected_derivative(fe_midpoint_value,
824 cell, active_neighbors);
829 auto neighbor_ptr = active_neighbors.begin();
830 for (; neighbor_ptr != active_neighbors.end(); ++neighbor_ptr)
832 const auto neighbor = *neighbor_ptr;
835 x_fe_midpoint_value.
reinit(neighbor);
841 const typename DerivativeDescription::ProjectedDerivative
842 neighbor_midpoint_value =
843 DerivativeDescription::get_projected_derivative(
844 neighbor_fe_midpoint_value, solution, component);
858 const double distance = y.
norm();
872 for (
unsigned int i = 0; i < dim; ++i)
873 for (
unsigned int j = 0; j < dim; ++j)
874 Y[i][j] += y[i] * y[j];
879 typename DerivativeDescription::ProjectedDerivative
880 projected_finite_difference =
881 (neighbor_midpoint_value - this_midpoint_value);
882 projected_finite_difference /= distance;
884 projected_derivative +=
outer_product(y, projected_finite_difference);
901 derivative = Y_inverse * projected_derivative;
904 DerivativeDescription::symmetrize(derivative);
914 template <
class DerivativeDescription,
925 const InputVector &solution,
926 const unsigned int component)
929 if (std::get<0>(*cell)->is_locally_owned() ==
false)
930 *std::get<1>(*cell) = 0;
933 typename DerivativeDescription::Derivative derivative;
946 *std::get<1>(*cell) =
947 DerivativeDescription::derivative_norm(derivative);
962 template <
class DerivativeDescription,
969 const InputVector &solution,
970 const unsigned int component,
981 std::tuple<typename DoFHandler<dim, spacedim>::active_cell_iterator,
993 [&mapping, &dof_handler, &solution, component](
998 cell, mapping, dof_handler, solution, component);
1000 std::function<void(const internal::Assembler::CopyData &)>(),
1014 template <
int dim,
class InputVector,
int spacedim>
1018 const InputVector &solution,
1020 const unsigned int component)
1027 template <
int dim,
class InputVector,
int spacedim>
1030 const InputVector &solution,
1032 const unsigned int component)
1036 const auto reference_cell =
1047 template <
int dim,
class InputVector,
int spacedim>
1051 const InputVector &solution,
1053 const unsigned int component)
1060 template <
int dim,
class InputVector,
int spacedim>
1063 const InputVector &solution,
1065 const unsigned int component)
1069 const auto reference_cell =
1080 template <
int dim,
int spacedim,
class InputVector,
int order>
1085 const InputVector &solution,
1093 const unsigned int component)
1097 mapping, dof, solution, component, cell, derivative);
1102 template <
int dim,
int spacedim,
class InputVector,
int order>
1106 const InputVector &solution,
1114 const unsigned int component)
1118 cell->reference_cell()
1129 template <
int dim,
int order>
1141#include "derivative_approximation.inst"
static ProjectedDerivative get_projected_derivative(const FEValues< dim, spacedim > &fe_values, const InputVector &solution, const unsigned int component)
static const UpdateFlags update_flags
static void symmetrize(Derivative &derivative_tensor)
static double derivative_norm(const Derivative &d)
static const UpdateFlags update_flags
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)
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 double derivative_norm(const Derivative &d)
static const UpdateFlags update_flags
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.
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
const FEValuesType & get_present_fe_values() const
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)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
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.
#define DEAL_II_NOT_IMPLEMENTED()
const Mapping< dim, spacedim > & get_default_linear_mapping(const Triangulation< dim, spacedim > &triangulation)
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)
constexpr T fixed_power(const T t)
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 constexpr double PI
::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)
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 > &)