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/parallel_vector.h> 21 #include <deal.II/lac/parallel_block_vector.h> 22 #include <deal.II/lac/petsc_vector.h> 23 #include <deal.II/lac/petsc_block_vector.h> 24 #include <deal.II/lac/trilinos_vector.h> 25 #include <deal.II/lac/trilinos_block_vector.h> 26 #include <deal.II/grid/tria_iterator.h> 27 #include <deal.II/grid/grid_tools.h> 28 #include <deal.II/grid/filtered_iterator.h> 29 #include <deal.II/dofs/dof_accessor.h> 30 #include <deal.II/dofs/dof_handler.h> 31 #include <deal.II/fe/fe.h> 32 #include <deal.II/fe/fe_values.h> 33 #include <deal.II/hp/fe_values.h> 34 #include <deal.II/fe/mapping_q1.h> 35 #include <deal.II/hp/q_collection.h> 36 #include <deal.II/hp/fe_collection.h> 37 #include <deal.II/hp/mapping_collection.h> 38 #include <deal.II/numerics/derivative_approximation.h> 42 DEAL_II_NAMESPACE_OPEN
49 inline T sqr (
const T t)
98 template <
class InputVector,
int spacedim>
99 static ProjectedDerivative
101 const InputVector &solution,
102 const unsigned int component);
117 static void symmetrize (Derivative &derivative_tensor);
126 template <
class InputVector,
int spacedim>
131 const InputVector &solution,
132 const unsigned int component)
134 if (fe_values.
get_fe().n_components() == 1)
136 std::vector<typename InputVector::value_type> values (1);
142 std::vector<Vector<typename InputVector::value_type> > values
145 return values[0](component);
157 for (
unsigned int i=0; i<dim; ++i)
210 template <
class InputVector,
int spacedim>
211 static ProjectedDerivative
213 const InputVector &solution,
214 const unsigned int component);
234 static void symmetrize (Derivative &derivative_tensor);
242 template <
class InputVector,
int spacedim>
247 const InputVector &solution,
248 const unsigned int component)
250 if (fe_values.
get_fe().n_components() == 1)
252 std::vector<Tensor<1,dim,typename InputVector::value_type> > values (1);
258 std::vector<std::vector<Tensor<1,dim,typename InputVector::value_type> > > values
273 return std::fabs (d[0][0]);
292 const double radicand = ::sqr(d[0][0] - d[1][1]) +
294 const double eigenvalues[2]
295 = { 0.5*(d[0][0] + d[1][1] + std::sqrt(radicand)),
296 0.5*(d[0][0] + d[1][1] - std::sqrt(radicand))
299 return std::max (std::fabs (eigenvalues[0]),
300 std::fabs (eigenvalues[1]));
420 const double am = trace(d) / 3.;
424 for (
unsigned int i=0; i<3; ++i)
427 const double ss01 = s[0][1] * s[0][1],
428 ss12 = s[1][2] * s[1][2],
429 ss02 = s[0][2] * s[0][2];
431 const double J2 = (s[0][0]*s[0][0] + s[1][1]*s[1][1] + s[2][2]*s[2][2]
432 + 2 * (ss01 + ss02 + ss12)) / 2.;
433 const double J3 = (std::pow(s[0][0],3) + std::pow(s[1][1],3) + std::pow(s[2][2],3)
434 + 3. * s[0][0] * (ss01 + ss02)
435 + 3. * s[1][1] * (ss01 + ss12)
436 + 3. * s[2][2] * (ss02 + ss12)
437 + 6. * s[0][1] * s[0][2] * s[1][2]) / 3.;
439 const double R = std::sqrt (4. * J2 / 3.);
441 double EE[3] = { 0, 0, 0 };
448 if (R <= 1e-14*std::fabs(am))
449 EE[0] = EE[1] = EE[2] = am;
454 const double R3 = R*R*R;
455 const double XX = 4. * J3 / R3;
456 const double YY = 1. - std::fabs(XX);
458 Assert (YY > -1e-14, ExcInternalError());
463 const double a = (XX>0 ? -1. : 1.) * R / 2;
464 EE[0] = EE[1] = am + a;
469 const double theta = std::acos(XX) / 3.;
470 EE[0] = am + R*std::cos(theta);
471 EE[1] = am + R*std::cos(theta + 2./3.*
numbers::PI);
472 EE[2] = am + R*std::cos(theta + 4./3.*
numbers::PI);
476 return std::max (std::fabs (EE[0]),
477 std::max (std::fabs (EE[1]),
502 Assert (
false, ExcNotImplemented());
514 for (
unsigned int i=0; i<dim; ++i)
515 for (
unsigned int j=i+1; j<dim; ++j)
517 const double s = (d[i][j] + d[j][i]) / 2;
518 d[i][j] = d[j][i] = s;
525 class ThirdDerivative
553 template <
class InputVector,
int spacedim>
554 static ProjectedDerivative
556 const InputVector &solution,
557 const unsigned int component);
577 static void symmetrize (Derivative &derivative_tensor);
585 template <
class InputVector,
int spacedim>
588 ThirdDerivative<dim>::
590 const InputVector &solution,
591 const unsigned int component)
593 if (fe_values.
get_fe().n_components() == 1)
595 std::vector<Tensor<2,dim,typename InputVector::value_type> > values (1);
601 std::vector<std::vector<Tensor<2,dim,typename InputVector::value_type> > > values
616 return std::fabs (d[0][0][0]);
624 ThirdDerivative<dim>::
636 ThirdDerivative<dim>::symmetrize (
Derivative &d)
643 for (
unsigned int i=0; i<dim; ++i)
644 for (
unsigned int j=i+1; j<dim; ++j)
645 for (
unsigned int k=j+1; k<dim; ++k)
647 const double s = (d[i][j][k] +
663 for (
unsigned int i=0; i<dim; ++i)
664 for (
unsigned int j=i+1; j<dim; ++j)
668 const double s = (d[i][i][j] +
678 const double t = (d[i][j][j] +
689 template <
int order,
int dim>
690 class DerivativeSelector
698 typedef void DerivDescr;
703 class DerivativeSelector<1,dim>
711 class DerivativeSelector<2,dim>
719 class DerivativeSelector<3,dim>
723 typedef ThirdDerivative<dim> DerivDescr;
759 template <
class DerivativeDescription,
int dim,
template <
int,
int>
class DoFHandlerType,
760 class InputVector,
int spacedim>
763 const DoFHandlerType<dim,spacedim> &dof_handler,
764 const InputVector &solution,
765 const unsigned int component,
768 typename DerivativeDescription::Derivative &derivative)
785 DerivativeDescription::update_flags |
796 std::vector<TriaActiveIterator<::DoFCellAccessor<DoFHandlerType<dim, spacedim>,
797 false> > > active_neighbors;
806 typename DerivativeDescription::Derivative projected_derivative;
809 x_fe_midpoint_value.
reinit (cell);
815 const typename DerivativeDescription::ProjectedDerivative
817 = DerivativeDescription::get_projected_derivative (fe_midpoint_value,
835 GridTools::get_active_neighbors<DoFHandlerType<dim,spacedim> >(cell, active_neighbors);
840 typename std::vector<TriaActiveIterator<::DoFCellAccessor<DoFHandlerType<dim, spacedim>,
841 false> > >::const_iterator
842 neighbor_ptr = active_neighbors.begin();
843 for (; neighbor_ptr!=active_neighbors.end(); ++neighbor_ptr)
846 neighbor = *neighbor_ptr;
849 x_fe_midpoint_value.
reinit (neighbor);
855 const typename DerivativeDescription::ProjectedDerivative
856 neighbor_midpoint_value
857 = DerivativeDescription::get_projected_derivative (neighbor_fe_midpoint_value,
858 solution, component);
871 const double distance = y.
norm();
885 for (
unsigned int i=0; i<dim; ++i)
886 for (
unsigned int j=0; j<dim; ++j)
887 Y[i][j] += y[i] * y[j];
892 typename DerivativeDescription::ProjectedDerivative
893 projected_finite_difference
894 = (neighbor_midpoint_value -
895 this_midpoint_value);
896 projected_finite_difference /= distance;
898 projected_derivative += outer_product(y, projected_finite_difference);
911 ExcInsufficientDirections());
916 derivative = Y_inverse * projected_derivative;
919 DerivativeDescription::symmetrize (derivative);
929 template <
class DerivativeDescription,
int dim,
930 template <
int,
int>
class DoFHandlerType,
class InputVector,
int spacedim>
935 const DoFHandlerType<dim,spacedim> &dof_handler,
936 const InputVector &solution,
937 const unsigned int component)
940 if (std_cxx11::get<0>(cell.iterators)->is_locally_owned() ==
false)
941 *std_cxx11::get<1>(cell.iterators) = 0;
944 typename DerivativeDescription::Derivative derivative;
947 approximate_cell<DerivativeDescription,dim,DoFHandlerType,InputVector, spacedim>
948 (mapping,dof_handler,solution,component,std_cxx11::get<0>(cell.iterators),derivative);
952 *std_cxx11::get<1>(cell.iterators) = DerivativeDescription::derivative_norm (derivative);
967 template <
class DerivativeDescription,
int dim,
968 template <
int,
int>
class DoFHandlerType,
class InputVector,
int spacedim>
971 const DoFHandlerType<dim,spacedim> &dof_handler,
972 const InputVector &solution,
973 const unsigned int component,
976 Assert (derivative_norm.
size() == dof_handler.get_triangulation().n_active_cells(),
977 ExcVectorLengthVsNActiveCells (derivative_norm.
size(),
978 dof_handler.get_triangulation().n_active_cells()));
979 Assert (component < dof_handler.get_fe().n_components(),
980 ExcIndexRange (component, 0, dof_handler.get_fe().n_components()));
983 <DoFHandlerType<dim, spacedim>,
false> >,
984 Vector<float>::iterator> Iterators;
986 derivative_norm.
begin())),
987 end(Iterators(dof_handler.end(),
988 derivative_norm.
end()));
994 static_cast<std_cxx11::function<void (SynchronousIterators<Iterators>
const &,
995 Assembler::Scratch
const &, Assembler::CopyData &)
> >
996 (std_cxx11::bind(&approximate<DerivativeDescription,dim,DoFHandlerType,
997 InputVector,spacedim>,
999 std_cxx11::cref(mapping),
1000 std_cxx11::cref(dof_handler),
1001 std_cxx11::cref(solution),component)),
1002 std_cxx11::function<
void (internal::Assembler::CopyData
const &)> (),
1003 internal::Assembler::Scratch (),internal::Assembler::CopyData ());
1015 template <
int dim,
template <
int,
int>
class DoFHandlerType,
class InputVector,
int spacedim>
1018 const DoFHandlerType<dim,spacedim> &dof_handler,
1019 const InputVector &solution,
1021 const unsigned int component)
1023 internal::approximate_derivative<internal::Gradient<dim>,dim> (mapping,
1031 template <
int dim,
template <
int,
int>
class DoFHandlerType,
class InputVector,
int spacedim>
1034 const InputVector &solution,
1036 const unsigned int component)
1046 template <
int dim,
template <
int,
int>
class DoFHandlerType,
class InputVector,
int spacedim>
1049 const DoFHandlerType<dim,spacedim> &dof_handler,
1050 const InputVector &solution,
1052 const unsigned int component)
1054 internal::approximate_derivative<internal::SecondDerivative<dim>,dim> (mapping,
1062 template <
int dim,
template <
int,
int>
class DoFHandlerType,
class InputVector,
int spacedim>
1065 const InputVector &solution,
1067 const unsigned int component)
1077 template <
typename DoFHandlerType,
int dim,
int spacedim,
class InputVector,
int order>
1081 const DoFHandlerType &dof,
1082 const InputVector &solution,
1084 const typename DoFHandlerType::active_cell_iterator &cell,
1089 const unsigned int component)
1091 internal::approximate_cell<typename internal::DerivativeSelector<order,DoFHandlerType::dimension>::DerivDescr>
1102 template <
typename DoFHandlerType,
int dim,
int spacedim,
class InputVector,
int order>
1105 (
const DoFHandlerType &dof,
1106 const InputVector &solution,
1108 const typename DoFHandlerType::active_cell_iterator &cell,
1113 const unsigned int component)
1116 approximate_derivative_tensor<DoFHandlerType, dim, spacedim, InputVector, order>
1129 template <
int dim,
int order>
1133 return internal::DerivativeSelector<order,dim>::DerivDescr::derivative_norm(derivative);
1140 #include "derivative_approximation.inst" 1143 DEAL_II_NAMESPACE_CLOSE
void get_function_gradients(const InputVector &fe_function, std::vector< Tensor< 1, spacedim, typename InputVector::value_type > > &gradients) const
void approximate_derivative_tensor(const Mapping< dim, spacedim > &mapping, const DoFHandlerType &dof, const InputVector &solution, const typename DoFHandlerType::active_cell_iterator &cell, Tensor< order, dim > &derivative, const unsigned int component=0)
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 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 ::FEValues< dim, spacedim > & get_present_fe_values() const
const Point< spacedim > & quadrature_point(const unsigned int q) const
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 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.
Tensor< 0, dim > ProjectedDerivative
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 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)