16 #include <deal.II/base/quadrature_lib.h> 17 #include <deal.II/base/work_stream.h> 18 #include <deal.II/lac/vector.h> 19 #include <deal.II/lac/block_vector.h> 20 #include <deal.II/lac/la_vector.h> 21 #include <deal.II/lac/la_parallel_vector.h> 22 #include <deal.II/lac/la_parallel_block_vector.h> 23 #include <deal.II/lac/petsc_vector.h> 24 #include <deal.II/lac/petsc_block_vector.h> 25 #include <deal.II/lac/trilinos_vector.h> 26 #include <deal.II/lac/trilinos_block_vector.h> 27 #include <deal.II/grid/tria_iterator.h> 28 #include <deal.II/grid/grid_tools.h> 29 #include <deal.II/grid/filtered_iterator.h> 30 #include <deal.II/dofs/dof_accessor.h> 31 #include <deal.II/dofs/dof_handler.h> 32 #include <deal.II/fe/fe.h> 33 #include <deal.II/fe/fe_values.h> 34 #include <deal.II/hp/fe_values.h> 35 #include <deal.II/fe/mapping_q1.h> 36 #include <deal.II/hp/q_collection.h> 37 #include <deal.II/hp/fe_collection.h> 38 #include <deal.II/hp/mapping_collection.h> 39 #include <deal.II/numerics/derivative_approximation.h> 43 DEAL_II_NAMESPACE_OPEN
50 inline T sqr (
const T t)
99 template <
class InputVector,
int spacedim>
102 const InputVector &solution,
103 const unsigned int component);
127 template <
class InputVector,
int spacedim>
132 const InputVector &solution,
133 const unsigned int component)
135 if (fe_values.
get_fe().n_components() == 1)
137 std::vector<typename InputVector::value_type> values (1);
143 std::vector<Vector<typename InputVector::value_type> > values
146 return values[0](component);
158 for (
unsigned int i=0; i<dim; ++i)
211 template <
class InputVector,
int spacedim>
214 const InputVector &solution,
215 const unsigned int component);
243 template <
class InputVector,
int spacedim>
248 const InputVector &solution,
249 const unsigned int component)
251 if (fe_values.
get_fe().n_components() == 1)
253 std::vector<Tensor<1,dim,typename InputVector::value_type> > values (1);
259 std::vector<std::vector<Tensor<1,dim,typename InputVector::value_type> > > values
274 return std::fabs (d[0][0]);
293 const double radicand = ::sqr(d[0][0] - d[1][1]) +
295 const double eigenvalues[2]
296 = { 0.5*(d[0][0] + d[1][1] + std::sqrt(radicand)),
297 0.5*(d[0][0] + d[1][1] - std::sqrt(radicand))
300 return std::max (std::fabs (eigenvalues[0]),
301 std::fabs (eigenvalues[1]));
421 const double am = trace(d) / 3.;
425 for (
unsigned int i=0; i<3; ++i)
428 const double ss01 = s[0][1] * s[0][1],
429 ss12 = s[1][2] * s[1][2],
430 ss02 = s[0][2] * s[0][2];
432 const double J2 = (s[0][0]*s[0][0] + s[1][1]*s[1][1] + s[2][2]*s[2][2]
433 + 2 * (ss01 + ss02 + ss12)) / 2.;
434 const double J3 = (std::pow(s[0][0],3) + std::pow(s[1][1],3) + std::pow(s[2][2],3)
435 + 3. * s[0][0] * (ss01 + ss02)
436 + 3. * s[1][1] * (ss01 + ss12)
437 + 3. * s[2][2] * (ss02 + ss12)
438 + 6. * s[0][1] * s[0][2] * s[1][2]) / 3.;
440 const double R = std::sqrt (4. * J2 / 3.);
442 double EE[3] = { 0, 0, 0 };
449 if (R <= 1e-14*std::fabs(am))
450 EE[0] = EE[1] = EE[2] = am;
455 const double R3 = R*R*R;
456 const double XX = 4. * J3 / R3;
457 const double YY = 1. - std::fabs(XX);
464 const double a = (XX>0 ? -1. : 1.) * R / 2;
465 EE[0] = EE[1] = am + a;
470 const double theta = std::acos(XX) / 3.;
471 EE[0] = am + R*std::cos(theta);
472 EE[1] = am + R*std::cos(theta + 2./3.*
numbers::PI);
473 EE[2] = am + R*std::cos(theta + 4./3.*
numbers::PI);
477 return std::max (std::fabs (EE[0]),
478 std::max (std::fabs (EE[1]),
515 for (
unsigned int i=0; i<dim; ++i)
516 for (
unsigned int j=i+1; j<dim; ++j)
518 const double s = (d[i][j] + d[j][i]) / 2;
519 d[i][j] = d[j][i] = s;
526 class ThirdDerivative
554 template <
class InputVector,
int spacedim>
555 static ProjectedDerivative
557 const InputVector &solution,
558 const unsigned int component);
578 static void symmetrize (Derivative &derivative_tensor);
586 template <
class InputVector,
int spacedim>
588 typename ThirdDerivative<dim>::ProjectedDerivative
589 ThirdDerivative<dim>::
591 const InputVector &solution,
592 const unsigned int component)
594 if (fe_values.
get_fe().n_components() == 1)
596 std::vector<Tensor<2,dim,typename InputVector::value_type> > values (1);
598 return ProjectedDerivative(values[0]);
602 std::vector<std::vector<Tensor<2,dim,typename InputVector::value_type> > > values
605 return ProjectedDerivative(values[0][component]);
615 derivative_norm (
const Derivative &d)
617 return std::fabs (d[0][0][0]);
625 ThirdDerivative<dim>::
626 derivative_norm (
const Derivative &d)
637 ThirdDerivative<dim>::symmetrize (Derivative &d)
644 for (
unsigned int i=0; i<dim; ++i)
645 for (
unsigned int j=i+1; j<dim; ++j)
646 for (
unsigned int k=j+1; k<dim; ++k)
648 const double s = (
d[i][j][k] +
664 for (
unsigned int i=0; i<dim; ++i)
665 for (
unsigned int j=i+1; j<dim; ++j)
669 const double s = (
d[i][i][j] +
679 const double t = (
d[i][j][j] +
690 template <
int order,
int dim>
691 class DerivativeSelector
699 typedef void DerivDescr;
704 class DerivativeSelector<1,dim>
708 typedef Gradient<dim> DerivDescr;
712 class DerivativeSelector<2,dim>
716 typedef SecondDerivative<dim> DerivDescr;
720 class DerivativeSelector<3,dim>
724 typedef ThirdDerivative<dim> DerivDescr;
760 template <
class DerivativeDescription,
int dim,
template <
int,
int>
class DoFHandlerType,
761 class InputVector,
int spacedim>
764 const DoFHandlerType<dim,spacedim> &dof_handler,
765 const InputVector &solution,
766 const unsigned int component,
769 typename DerivativeDescription::Derivative &derivative)
786 DerivativeDescription::update_flags |
797 std::vector<TriaActiveIterator<::DoFCellAccessor<DoFHandlerType<dim, spacedim>,
798 false> > > active_neighbors;
807 typename DerivativeDescription::Derivative projected_derivative;
810 x_fe_midpoint_value.reinit (cell);
816 const typename DerivativeDescription::ProjectedDerivative
818 = DerivativeDescription::get_projected_derivative (fe_midpoint_value,
836 GridTools::get_active_neighbors<DoFHandlerType<dim,spacedim> >(cell, active_neighbors);
841 typename std::vector<TriaActiveIterator<::DoFCellAccessor<DoFHandlerType<dim, spacedim>,
842 false> > >::const_iterator
843 neighbor_ptr = active_neighbors.begin();
844 for (; neighbor_ptr!=active_neighbors.end(); ++neighbor_ptr)
847 neighbor = *neighbor_ptr;
850 x_fe_midpoint_value.reinit (neighbor);
856 const typename DerivativeDescription::ProjectedDerivative
857 neighbor_midpoint_value
858 = DerivativeDescription::get_projected_derivative (neighbor_fe_midpoint_value,
859 solution, component);
872 const double distance = y.
norm();
886 for (
unsigned int i=0; i<dim; ++i)
887 for (
unsigned int j=0; j<dim; ++j)
888 Y[i][j] += y[i] * y[j];
893 typename DerivativeDescription::ProjectedDerivative
894 projected_finite_difference
895 = (neighbor_midpoint_value -
896 this_midpoint_value);
897 projected_finite_difference /= distance;
899 projected_derivative += outer_product(y, projected_finite_difference);
917 derivative = Y_inverse * projected_derivative;
920 DerivativeDescription::symmetrize (derivative);
930 template <
class DerivativeDescription,
int dim,
931 template <
int,
int>
class DoFHandlerType,
class InputVector,
int spacedim>
936 const DoFHandlerType<dim,spacedim> &dof_handler,
937 const InputVector &solution,
938 const unsigned int component)
941 if (std_cxx11::get<0>(*cell)->is_locally_owned() ==
false)
942 *std_cxx11::get<1>(*cell) = 0;
945 typename DerivativeDescription::Derivative derivative;
948 approximate_cell<DerivativeDescription,dim,DoFHandlerType,InputVector, spacedim>
949 (mapping,dof_handler,solution,component,std_cxx11::get<0>(*cell),derivative);
953 *std_cxx11::get<1>(*cell) = DerivativeDescription::derivative_norm (derivative);
968 template <
class DerivativeDescription,
int dim,
969 template <
int,
int>
class DoFHandlerType,
class InputVector,
int spacedim>
972 const DoFHandlerType<dim,spacedim> &dof_handler,
973 const InputVector &solution,
974 const unsigned int component,
979 dof_handler.get_triangulation().n_active_cells()));
980 Assert (component < dof_handler.get_fe().n_components(),
981 ExcIndexRange (component, 0, dof_handler.get_fe().n_components()));
984 <DoFHandlerType<dim, spacedim>,
false> >,
985 Vector<float>::iterator> Iterators;
988 end(Iterators(dof_handler.end(),
995 static_cast<std_cxx11::function<void (SynchronousIterators<Iterators>
const &,
996 Assembler::Scratch
const &, Assembler::CopyData &)
> >
997 (std_cxx11::bind(&approximate<DerivativeDescription,dim,DoFHandlerType,
998 InputVector,spacedim>,
1000 std_cxx11::cref(mapping),
1001 std_cxx11::cref(dof_handler),
1002 std_cxx11::cref(solution),component)),
1003 std_cxx11::function<
void (internal::Assembler::CopyData
const &)> (),
1004 internal::Assembler::Scratch (),internal::Assembler::CopyData ());
1016 template <
int dim,
template <
int,
int>
class DoFHandlerType,
class InputVector,
int spacedim>
1019 const DoFHandlerType<dim,spacedim> &dof_handler,
1020 const InputVector &solution,
1022 const unsigned int component)
1024 internal::approximate_derivative<internal::Gradient<dim>,dim> (mapping,
1032 template <
int dim,
template <
int,
int>
class DoFHandlerType,
class InputVector,
int spacedim>
1035 const InputVector &solution,
1037 const unsigned int component)
1047 template <
int dim,
template <
int,
int>
class DoFHandlerType,
class InputVector,
int spacedim>
1050 const DoFHandlerType<dim,spacedim> &dof_handler,
1051 const InputVector &solution,
1053 const unsigned int component)
1055 internal::approximate_derivative<internal::SecondDerivative<dim>,dim> (mapping,
1063 template <
int dim,
template <
int,
int>
class DoFHandlerType,
class InputVector,
int spacedim>
1066 const InputVector &solution,
1068 const unsigned int component)
1078 template <
typename DoFHandlerType,
class InputVector,
int order>
1082 const DoFHandlerType &dof,
1083 const InputVector &solution,
1085 const typename DoFHandlerType::active_cell_iterator &cell,
1090 const unsigned int component)
1092 internal::approximate_cell<typename internal::DerivativeSelector<order,DoFHandlerType::dimension>::DerivDescr>
1103 template <
typename DoFHandlerType,
class InputVector,
int order>
1106 (
const DoFHandlerType &dof,
1107 const InputVector &solution,
1109 const typename DoFHandlerType::active_cell_iterator &cell,
1114 const unsigned int component)
1130 template <
int dim,
int order>
1134 return internal::DerivativeSelector<order,dim>::DerivDescr::derivative_norm(derivative);
1141 #include "derivative_approximation.inst" 1144 DEAL_II_NAMESPACE_CLOSE
void get_function_gradients(const InputVector &fe_function, std::vector< Tensor< 1, spacedim, typename InputVector::value_type > > &gradients) const
static ::ExceptionBase & ExcInsufficientDirections()
const FEValues< dim, spacedim > & get_present_fe_values() const
static void symmetrize(Derivative &derivative_tensor)
static const UpdateFlags update_flags
void get_function_values(const InputVector &fe_function, std::vector< typename InputVector::value_type > &values) const
const FiniteElement< dim, spacedim > & get_fe() const
Transformed quadrature points.
#define AssertThrow(cond, exc)
Tensor< 1, dim > ProjectedDerivative
Tensor< 1, dim > Derivative
numbers::NumberTraits< Number >::real_type norm() const
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
static ProjectedDerivative get_projected_derivative(const FEValues< dim, spacedim > &fe_values, const InputVector &solution, const unsigned int component)
void get_function_hessians(const InputVector &fe_function, std::vector< Tensor< 2, spacedim, typename InputVector::value_type > > &hessians) const
const Point< spacedim > & quadrature_point(const unsigned int q) const
static const UpdateFlags update_flags
Tensor< 2, dim > Derivative
#define Assert(cond, exc)
static double derivative_norm(const Derivative &d)
Abstract base class for mapping classes.
static void symmetrize(Derivative &derivative_tensor)
static ProjectedDerivative get_projected_derivative(const FEValues< dim, spacedim > &fe_values, const InputVector &solution, const unsigned int component)
Second derivatives of shape functions.
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
Tensor< 0, dim > ProjectedDerivative
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)
Shape function gradients.
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 ::ExceptionBase & ExcNotImplemented()
static double derivative_norm(const Derivative &d)
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)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcVectorLengthVsNActiveCells(int arg1, int arg2)