16 #ifndef dealii__fe_values_h 17 #define dealii__fe_values_h 20 #include <deal.II/base/config.h> 21 #include <deal.II/base/exceptions.h> 22 #include <deal.II/base/subscriptor.h> 23 #include <deal.II/base/point.h> 24 #include <deal.II/base/derivative_form.h> 25 #include <deal.II/base/symmetric_tensor.h> 26 #include <deal.II/base/vector_slice.h> 27 #include <deal.II/base/quadrature.h> 28 #include <deal.II/base/table.h> 29 #include <deal.II/base/std_cxx11/unique_ptr.h> 30 #include <deal.II/grid/tria.h> 31 #include <deal.II/grid/tria_iterator.h> 32 #include <deal.II/dofs/dof_handler.h> 33 #include <deal.II/dofs/dof_accessor.h> 34 #include <deal.II/hp/dof_handler.h> 35 #include <deal.II/fe/fe.h> 36 #include <deal.II/fe/fe_update_flags.h> 37 #include <deal.II/fe/fe_values_extractors.h> 38 #include <deal.II/fe/mapping.h> 45 #ifdef DEAL_II_WITH_PETSC 49 DEAL_II_NAMESPACE_OPEN
136 template <
int dim,
int spacedim=dim>
228 value (
const unsigned int shape_function,
229 const unsigned int q_point)
const;
242 gradient (
const unsigned int shape_function,
243 const unsigned int q_point)
const;
256 hessian (
const unsigned int shape_function,
257 const unsigned int q_point)
const;
271 const unsigned int q_point)
const;
290 template <
class InputVector>
311 template <
class InputVector>
332 template <
class InputVector>
354 template <
class InputVector>
376 template <
class InputVector>
379 typename InputVector::value_type>::type> &third_derivatives)
const;
430 template <
int dim,
int spacedim=dim>
476 typedef typename ::internal::CurlType<spacedim>::type
curl_type;
528 unsigned int single_nonzero_component_index;
570 value (
const unsigned int shape_function,
571 const unsigned int q_point)
const;
587 gradient (
const unsigned int shape_function,
588 const unsigned int q_point)
const;
607 const unsigned int q_point)
const;
620 divergence (
const unsigned int shape_function,
621 const unsigned int q_point)
const;
644 curl (
const unsigned int shape_function,
645 const unsigned int q_point)
const;
658 hessian (
const unsigned int shape_function,
659 const unsigned int q_point)
const;
673 const unsigned int q_point)
const;
692 template <
class InputVector>
713 template <
class InputVector>
740 template <
class InputVector>
763 template <
class InputVector>
785 template <
class InputVector>
806 template <
class InputVector>
828 template <
class InputVector>
850 template <
class InputVector>
853 typename InputVector::value_type>::type> &third_derivatives)
const;
874 template <
int rank,
int dim,
int spacedim = dim>
901 template <
int dim,
int spacedim>
928 struct ShapeFunctionData
938 bool is_nonzero_shape_function_component[value_type::n_independent_components];
949 unsigned int row_index[value_type::n_independent_components];
960 unsigned int single_nonzero_component_index;
978 const unsigned int first_tensor_component);
1004 value (
const unsigned int shape_function,
1005 const unsigned int q_point)
const;
1022 divergence (
const unsigned int shape_function,
1023 const unsigned int q_point)
const;
1042 template <
class InputVector>
1043 void get_function_values (
const InputVector &fe_function,
1067 template <
class InputVector>
1068 void get_function_divergences (
const InputVector &fe_function,
1090 template <
int rank,
int dim,
int spacedim = dim>
1112 template <
int dim,
int spacedim>
1132 struct ShapeFunctionData
1142 bool is_nonzero_shape_function_component[value_type::n_independent_components];
1153 unsigned int row_index[value_type::n_independent_components];
1164 unsigned int single_nonzero_component_index;
1183 const unsigned int first_tensor_component);
1209 value (
const unsigned int shape_function,
1210 const unsigned int q_point)
const;
1226 divergence (
const unsigned int shape_function,
1227 const unsigned int q_point)
const;
1246 template <
class InputVector>
1247 void get_function_values (
const InputVector &fe_function,
1272 template <
class InputVector>
1273 void get_function_divergences (
const InputVector &fe_function,
1308 template <
int dim,
int spacedim>
1315 std::vector<::FEValuesViews::Scalar<dim,spacedim> >
scalars;
1316 std::vector<::FEValuesViews::Vector<dim,spacedim> > vectors;
1317 std::vector<::FEValuesViews::SymmetricTensor<2,dim,spacedim> >
1318 symmetric_second_order_tensors;
1319 std::vector<::FEValuesViews::Tensor<2,dim,spacedim> >
1320 second_order_tensors;
1434 template <
int dim,
int spacedim>
1505 const double &
shape_value (
const unsigned int function_no,
1506 const unsigned int point_no)
const;
1529 const unsigned int point_no,
1530 const unsigned int component)
const;
1579 const unsigned int point_no,
1580 const unsigned int component)
const;
1603 const unsigned int point_no)
const;
1623 const unsigned int point_no,
1624 const unsigned int component)
const;
1647 const unsigned int point_no)
const;
1667 const unsigned int point_no,
1668 const unsigned int component)
const;
1712 template <
class InputVector>
1714 std::vector<typename InputVector::value_type> &values)
const;
1729 template <
class InputVector>
1751 template <
class InputVector>
1753 const VectorSlice<
const std::vector<types::global_dof_index> > &indices,
1754 std::vector<typename InputVector::value_type> &values)
const;
1777 template <
class InputVector>
1779 const VectorSlice<
const std::vector<types::global_dof_index> > &indices,
1813 template <
class InputVector>
1815 const VectorSlice<
const std::vector<types::global_dof_index> > &indices,
1816 VectorSlice<std::vector<std::vector<typename InputVector::value_type> > > values,
1817 const bool quadrature_points_fastest)
const;
1861 template <
class InputVector>
1881 template <
class InputVector>
1891 template <
class InputVector>
1893 const VectorSlice<
const std::vector<types::global_dof_index> > &indices,
1902 template <
class InputVector>
1904 const VectorSlice<
const std::vector<types::global_dof_index> > &indices,
1906 bool quadrature_points_fastest =
false)
const;
1951 template <
class InputVector>
1973 template <
class InputVector>
1977 bool quadrature_points_fastest =
false)
const;
1983 template <
class InputVector>
1985 const InputVector &fe_function,
1986 const VectorSlice<
const std::vector<types::global_dof_index> > &indices,
1995 template <
class InputVector>
1997 const InputVector &fe_function,
1998 const VectorSlice<
const std::vector<types::global_dof_index> > &indices,
2000 bool quadrature_points_fastest =
false)
const;
2044 template <
class InputVector>
2047 std::vector<typename InputVector::value_type> &laplacians)
const;
2068 template <
class InputVector>
2079 template <
class InputVector>
2081 const InputVector &fe_function,
2082 const VectorSlice<
const std::vector<types::global_dof_index> > &indices,
2083 std::vector<typename InputVector::value_type> &laplacians)
const;
2091 template <
class InputVector>
2093 const InputVector &fe_function,
2094 const VectorSlice<
const std::vector<types::global_dof_index> > &indices,
2103 template <
class InputVector>
2105 const InputVector &fe_function,
2106 const VectorSlice<
const std::vector<types::global_dof_index> > &indices,
2107 std::vector<std::vector<typename InputVector::value_type> > &laplacians,
2108 bool quadrature_points_fastest =
false)
const;
2154 template <
class InputVector>
2177 template <
class InputVector>
2181 bool quadrature_points_fastest =
false)
const;
2187 template <
class InputVector>
2189 const InputVector &fe_function,
2190 const VectorSlice<
const std::vector<types::global_dof_index> > &indices,
2199 template <
class InputVector>
2201 const InputVector &fe_function,
2202 const VectorSlice<
const std::vector<types::global_dof_index> > &indices,
2204 bool quadrature_points_fastest =
false)
const;
2261 const std::vector<DerivativeForm<1,dim,spacedim> > &
get_jacobians ()
const;
2438 const
std::vector<
Tensor<1,dim> > &original,
2541 << "
object for which this kind of information has not been computed. What "
2542 << "information these objects compute is determined by the update_* flags you "
2543 << "pass to the constructor. Here, the operation you are attempting requires "
2544 << "the <" << arg1 << "> flag to be set, but it was apparently not specified "
2545 << "upon construction.");
2571 << "The shape function with index " << arg1
2572 << " is not primitive, i.e. it is vector-valued and "
2573 << "has more than one non-zero vector component. This "
2574 << "function cannot be called for these shape functions. "
2575 << "Maybe you want to use the same function with the "
2576 << "_component suffix?");
2613 class CellIteratorBase;
2775 template <
int dim,
int spacedim=dim>
2783 static const unsigned int integral_dimension = dim;
2809 template <
template <
int,
int>
class DoFHandlerType,
bool level_dof_access>
2886 template <
int dim,
int spacedim=dim>
2894 static const unsigned int integral_dimension = dim-1;
2929 const std::vector<Tensor<1,spacedim> > &get_boundary_forms ()
const;
2935 unsigned int get_face_index()
const;
2941 const Quadrature<dim-1> & get_quadrature ()
const;
2979 template <
int dim,
int spacedim=dim>
2995 static const unsigned int integral_dimension = dim-1;
3019 template <
template <
int,
int>
class DoFHandlerType,
bool level_dof_access>
3021 const unsigned int face_no);
3037 const unsigned int face_no);
3067 void do_reinit (
const unsigned int face_no);
3088 template <
int dim,
int spacedim=dim>
3106 static const unsigned int integral_dimension = dim-1;
3132 template <
template <
int,
int>
class DoFHandlerType,
bool level_dof_access>
3134 const unsigned int face_no,
3135 const unsigned int subface_no);
3151 const unsigned int face_no,
3152 const unsigned int subface_no);
3197 void do_reinit (
const unsigned int face_no,
3198 const unsigned int subface_no);
3209 template <
int dim,
int spacedim>
3211 typename Scalar<dim,spacedim>::value_type
3212 Scalar<dim,spacedim>::value (
const unsigned int shape_function,
3213 const unsigned int q_point)
const 3215 Assert (shape_function < fe_values->
fe->dofs_per_cell,
3216 ExcIndexRange (shape_function, 0, fe_values->fe->dofs_per_cell));
3223 if (shape_function_data[shape_function].is_nonzero_shape_function_component)
3224 return fe_values->finite_element_output.shape_values(shape_function_data[shape_function]
3234 template <
int dim,
int spacedim>
3236 typename Scalar<dim,spacedim>::gradient_type
3237 Scalar<dim,spacedim>::gradient (
const unsigned int shape_function,
3238 const unsigned int q_point)
const 3240 Assert (shape_function < fe_values->
fe->dofs_per_cell,
3241 ExcIndexRange (shape_function, 0, fe_values->fe->dofs_per_cell));
3251 if (shape_function_data[shape_function].is_nonzero_shape_function_component)
3252 return fe_values->finite_element_output.shape_gradients[shape_function_data[shape_function]
3253 .row_index][q_point];
3255 return gradient_type();
3260 template <
int dim,
int spacedim>
3262 typename Scalar<dim,spacedim>::hessian_type
3263 Scalar<dim,spacedim>::hessian (
const unsigned int shape_function,
3264 const unsigned int q_point)
const 3266 Assert (shape_function < fe_values->
fe->dofs_per_cell,
3267 ExcIndexRange (shape_function, 0, fe_values->fe->dofs_per_cell));
3277 if (shape_function_data[shape_function].is_nonzero_shape_function_component)
3278 return fe_values->finite_element_output.shape_hessians[shape_function_data[shape_function].row_index][q_point];
3280 return hessian_type();
3285 template <
int dim,
int spacedim>
3287 typename Scalar<dim,spacedim>::third_derivative_type
3288 Scalar<dim,spacedim>::third_derivative (
const unsigned int shape_function,
3289 const unsigned int q_point)
const 3291 Assert (shape_function < fe_values->
fe->dofs_per_cell,
3292 ExcIndexRange (shape_function, 0, fe_values->fe->dofs_per_cell));
3302 if (shape_function_data[shape_function].is_nonzero_shape_function_component)
3303 return fe_values->finite_element_output.shape_3rd_derivatives[shape_function_data[shape_function].row_index][q_point];
3305 return third_derivative_type();
3310 template <
int dim,
int spacedim>
3314 const unsigned int q_point)
const 3316 Assert (shape_function < fe_values->
fe->dofs_per_cell,
3317 ExcIndexRange (shape_function, 0, fe_values->fe->dofs_per_cell));
3323 const int snc = shape_function_data[shape_function].single_nonzero_component;
3325 return value_type();
3328 value_type return_value;
3329 return_value[shape_function_data[shape_function].single_nonzero_component_index]
3330 = fe_values->finite_element_output.shape_values(snc,q_point);
3331 return return_value;
3335 value_type return_value;
3336 for (
unsigned int d=0;
d<dim; ++
d)
3337 if (shape_function_data[shape_function].is_nonzero_shape_function_component[d])
3339 = fe_values->finite_element_output.shape_values(shape_function_data[shape_function].row_index[d],q_point);
3341 return return_value;
3347 template <
int dim,
int spacedim>
3351 const unsigned int q_point)
const 3353 Assert (shape_function < fe_values->
fe->dofs_per_cell,
3354 ExcIndexRange (shape_function, 0, fe_values->fe->dofs_per_cell));
3360 const int snc = shape_function_data[shape_function].single_nonzero_component;
3362 return gradient_type();
3365 gradient_type return_value;
3366 return_value[shape_function_data[shape_function].single_nonzero_component_index]
3367 = fe_values->finite_element_output.shape_gradients[snc][q_point];
3368 return return_value;
3372 gradient_type return_value;
3373 for (
unsigned int d=0;
d<dim; ++
d)
3374 if (shape_function_data[shape_function].is_nonzero_shape_function_component[d])
3376 = fe_values->finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[
d]][q_point];
3378 return return_value;
3384 template <
int dim,
int spacedim>
3388 const unsigned int q_point)
const 3392 Assert (shape_function < fe_values->
fe->dofs_per_cell,
3393 ExcIndexRange (shape_function, 0, fe_values->fe->dofs_per_cell));
3399 const int snc = shape_function_data[shape_function].single_nonzero_component;
3401 return divergence_type();
3404 fe_values->finite_element_output.shape_gradients[snc][q_point][shape_function_data[shape_function].single_nonzero_component_index];
3407 divergence_type return_value = 0;
3408 for (
unsigned int d=0;
d<dim; ++
d)
3409 if (shape_function_data[shape_function].is_nonzero_shape_function_component[d])
3411 += fe_values->finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[
d]][q_point][
d];
3413 return return_value;
3419 template <
int dim,
int spacedim>
3426 Assert (shape_function < fe_values->
fe->dofs_per_cell,
3427 ExcIndexRange (shape_function, 0, fe_values->fe->dofs_per_cell));
3431 const int snc = shape_function_data[shape_function].single_nonzero_component;
3434 return curl_type ();
3441 Assert (
false,
ExcMessage(
"Computing the curl in 1d is not a useful operation"));
3442 return curl_type ();
3449 curl_type return_value;
3455 if (shape_function_data[shape_function].single_nonzero_component_index == 0)
3456 return_value[0] = -1.0 * fe_values->finite_element_output.shape_gradients[snc][q_point][1];
3458 return_value[0] = fe_values->finite_element_output.shape_gradients[snc][q_point][0];
3460 return return_value;
3465 curl_type return_value;
3467 return_value[0] = 0.0;
3469 if (shape_function_data[shape_function].is_nonzero_shape_function_component[0])
3471 -= fe_values->finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[0]][q_point][1];
3473 if (shape_function_data[shape_function].is_nonzero_shape_function_component[1])
3475 += fe_values->finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[1]][q_point][0];
3477 return return_value;
3485 curl_type return_value;
3487 switch (shape_function_data[shape_function].single_nonzero_component_index)
3491 return_value[0] = 0;
3492 return_value[1] = fe_values->finite_element_output.shape_gradients[snc][q_point][2];
3493 return_value[2] = -1.0 * fe_values->finite_element_output.shape_gradients[snc][q_point][1];
3494 return return_value;
3499 return_value[0] = -1.0 * fe_values->finite_element_output.shape_gradients[snc][q_point][2];
3500 return_value[1] = 0;
3501 return_value[2] = fe_values->finite_element_output.shape_gradients[snc][q_point][0];
3502 return return_value;
3507 return_value[0] = fe_values->finite_element_output.shape_gradients[snc][q_point][1];
3508 return_value[1] = -1.0 * fe_values->finite_element_output.shape_gradients[snc][q_point][0];
3509 return_value[2] = 0;
3510 return return_value;
3517 curl_type return_value;
3519 for (
unsigned int i = 0; i < dim; ++i)
3520 return_value[i] = 0.0;
3522 if (shape_function_data[shape_function].is_nonzero_shape_function_component[0])
3525 += fe_values->finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[0]][q_point][2];
3527 -= fe_values->finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[0]][q_point][1];
3530 if (shape_function_data[shape_function].is_nonzero_shape_function_component[1])
3533 -= fe_values->finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[1]][q_point][2];
3535 += fe_values->finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[1]][q_point][0];
3538 if (shape_function_data[shape_function].is_nonzero_shape_function_component[2])
3541 += fe_values->finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[2]][q_point][1];
3543 -= fe_values->finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[2]][q_point][0];
3546 return return_value;
3555 template <
int dim,
int spacedim>
3559 const unsigned int q_point)
const 3563 Assert (shape_function < fe_values->
fe->dofs_per_cell,
3564 ExcIndexRange (shape_function, 0, fe_values->fe->dofs_per_cell));
3570 const int snc = shape_function_data[shape_function].single_nonzero_component;
3572 return hessian_type();
3575 hessian_type return_value;
3576 return_value[shape_function_data[shape_function].single_nonzero_component_index]
3577 = fe_values->finite_element_output.shape_hessians[snc][q_point];
3578 return return_value;
3582 hessian_type return_value;
3583 for (
unsigned int d=0;
d<dim; ++
d)
3584 if (shape_function_data[shape_function].is_nonzero_shape_function_component[d])
3586 = fe_values->finite_element_output.shape_hessians[shape_function_data[shape_function].row_index[
d]][q_point];
3588 return return_value;
3592 template <
int dim,
int spacedim>
3596 const unsigned int q_point)
const 3600 Assert (shape_function < fe_values->
fe->dofs_per_cell,
3601 ExcIndexRange (shape_function, 0, fe_values->fe->dofs_per_cell));
3607 const int snc = shape_function_data[shape_function].single_nonzero_component;
3609 return third_derivative_type();
3612 third_derivative_type return_value;
3613 return_value[shape_function_data[shape_function].single_nonzero_component_index]
3614 = fe_values->finite_element_output.shape_3rd_derivatives[snc][q_point];
3615 return return_value;
3619 third_derivative_type return_value;
3620 for (
unsigned int d=0;
d<dim; ++
d)
3621 if (shape_function_data[shape_function].is_nonzero_shape_function_component[d])
3623 = fe_values->finite_element_output.shape_3rd_derivatives[shape_function_data[shape_function].row_index[
d]][q_point];
3625 return return_value;
3637 ::SymmetricTensor<2,1>
3638 symmetrize_single_row (
const unsigned int n,
3639 const ::Tensor<1,1> &t)
3644 const double array[1] = { t[0] };
3645 return ::SymmetricTensor<2,1>(array);
3650 ::SymmetricTensor<2,2>
3651 symmetrize_single_row (
const unsigned int n,
3652 const ::Tensor<1,2> &t)
3658 const double array[3] = { t[0], 0, t[1]/2 };
3659 return ::SymmetricTensor<2,2>(array);
3663 const double array[3] = { 0, t[1], t[0]/2 };
3664 return ::SymmetricTensor<2,2>(array);
3669 return ::SymmetricTensor<2,2>();
3676 ::SymmetricTensor<2,3>
3677 symmetrize_single_row (
const unsigned int n,
3678 const ::Tensor<1,3> &t)
3684 const double array[6] = { t[0], 0, 0, t[1]/2, t[2]/2, 0 };
3685 return ::SymmetricTensor<2,3>(array);
3689 const double array[6] = { 0, t[1], 0, t[0]/2, 0, t[2]/2 };
3690 return ::SymmetricTensor<2,3>(array);
3694 const double array[6] = { 0, 0, t[2], 0, t[0]/2, t[1]/2 };
3695 return ::SymmetricTensor<2,3>(array);
3700 return ::SymmetricTensor<2,3>();
3707 template <
int dim,
int spacedim>
3711 const unsigned int q_point)
const 3713 Assert (shape_function < fe_values->
fe->dofs_per_cell,
3714 ExcIndexRange (shape_function, 0, fe_values->fe->dofs_per_cell));
3720 const int snc = shape_function_data[shape_function].single_nonzero_component;
3722 return symmetric_gradient_type();
3724 return symmetrize_single_row (shape_function_data[shape_function].single_nonzero_component_index,
3725 fe_values->finite_element_output.shape_gradients[snc][q_point]);
3728 gradient_type return_value;
3729 for (
unsigned int d=0;
d<dim; ++
d)
3730 if (shape_function_data[shape_function].is_nonzero_shape_function_component[d])
3732 = fe_values->finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[
d]][q_point];
3734 return symmetrize(return_value);
3740 template <
int dim,
int spacedim>
3744 const unsigned int q_point)
const 3746 Assert (shape_function < fe_values->
fe->dofs_per_cell,
3747 ExcIndexRange (shape_function, 0, fe_values->fe->dofs_per_cell));
3756 = shape_function_data[shape_function].single_nonzero_component;
3762 return value_type();
3767 value_type return_value;
3768 const unsigned int comp =
3769 shape_function_data[shape_function].single_nonzero_component_index;
3770 return_value[value_type::unrolled_to_component_indices(comp)]
3771 = fe_values->finite_element_output.shape_values(snc,q_point);
3772 return return_value;
3776 value_type return_value;
3777 for (
unsigned int d = 0;
d < value_type::n_independent_components; ++
d)
3778 if (shape_function_data[shape_function].is_nonzero_shape_function_component[d])
3779 return_value[value_type::unrolled_to_component_indices(d)]
3780 = fe_values->finite_element_output.shape_values(shape_function_data[shape_function].row_index[d],q_point);
3781 return return_value;
3786 template <
int dim,
int spacedim>
3790 const unsigned int q_point)
const 3792 Assert (shape_function < fe_values->
fe->dofs_per_cell,
3793 ExcIndexRange (shape_function, 0, fe_values->fe->dofs_per_cell));
3797 const int snc = shape_function_data[shape_function].single_nonzero_component;
3803 return divergence_type();
3836 const unsigned int comp =
3837 shape_function_data[shape_function].single_nonzero_component_index;
3838 const unsigned int ii = value_type::unrolled_to_component_indices(comp)[0];
3839 const unsigned int jj = value_type::unrolled_to_component_indices(comp)[1];
3856 const ::Tensor<1, spacedim> phi_grad = fe_values->finite_element_output.shape_gradients[snc][q_point];
3858 divergence_type return_value;
3859 return_value[ii] = phi_grad[jj];
3862 return_value[jj] = phi_grad[ii];
3864 return return_value;
3870 divergence_type return_value;
3871 return return_value;
3875 template <
int dim,
int spacedim>
3879 const unsigned int q_point)
const 3881 Assert (shape_function < fe_values->
fe->dofs_per_cell,
3882 ExcIndexRange (shape_function, 0, fe_values->fe->dofs_per_cell));
3891 = shape_function_data[shape_function].single_nonzero_component;
3897 return value_type();
3902 value_type return_value;
3903 const unsigned int comp =
3904 shape_function_data[shape_function].single_nonzero_component_index;
3906 return_value[indices] = fe_values->finite_element_output.shape_values(snc,q_point);
3907 return return_value;
3911 value_type return_value;
3912 for (
unsigned int d = 0;
d < dim*dim; ++
d)
3913 if (shape_function_data[shape_function].is_nonzero_shape_function_component[d])
3916 return_value[indices]
3917 = fe_values->finite_element_output.shape_values(shape_function_data[shape_function].row_index[d],q_point);
3919 return return_value;
3924 template <
int dim,
int spacedim>
3928 const unsigned int q_point)
const 3930 Assert (shape_function < fe_values->
fe->dofs_per_cell,
3931 ExcIndexRange (shape_function, 0, fe_values->fe->dofs_per_cell));
3935 const int snc = shape_function_data[shape_function].single_nonzero_component;
3941 return divergence_type();
3961 const unsigned int comp =
3962 shape_function_data[shape_function].single_nonzero_component_index;
3964 const unsigned int ii = indices[0];
3965 const unsigned int jj = indices[1];
3967 const ::Tensor<1, spacedim> phi_grad = fe_values->finite_element_output.shape_gradients[snc][q_point];
3969 divergence_type return_value;
3970 return_value[jj] = phi_grad[ii];
3972 return return_value;
3978 divergence_type return_value;
3979 return return_value;
3990 template <
int dim,
int spacedim>
4005 template <
int dim,
int spacedim>
4019 template <
int dim,
int spacedim>
4033 template <
int dim,
int spacedim>
4050 template <
int dim,
int spacedim>
4054 const unsigned int j)
const 4065 if (
fe->is_primitive())
4078 row = this->
finite_element_output.shape_function_to_row_table[i *
fe->n_components() +
fe->system_to_component_index(i).first];
4085 template <
int dim,
int spacedim>
4089 const unsigned int j,
4090 const unsigned int component)
const 4096 Assert (component < fe->n_components(),
4102 if (
fe->get_nonzero_components(i)[component] ==
false)
4115 template <
int dim,
int spacedim>
4119 const unsigned int j)
const 4130 if (
fe->is_primitive())
4143 row = this->finite_element_output.shape_function_to_row_table[i *
fe->n_components() +
fe->system_to_component_index(i).first];
4150 template <
int dim,
int spacedim>
4154 const unsigned int j,
4155 const unsigned int component)
const 4157 Assert (i < fe->dofs_per_cell,
4167 if (fe->get_nonzero_components(i)[component] ==
false)
4174 row = this->finite_element_output.shape_function_to_row_table[i * fe->n_components() + component];
4175 return this->finite_element_output.shape_gradients[row][j];
4180 template <
int dim,
int spacedim>
4184 const unsigned int j)
const 4186 Assert (i < fe->dofs_per_cell,
4190 Assert (fe->is_primitive (i),
4191 ExcShapeFunctionNotPrimitive(i));
4195 if (fe->is_primitive())
4196 return this->finite_element_output.shape_hessians[i][j];
4208 row = this->finite_element_output.shape_function_to_row_table[i * fe->n_components() + fe->system_to_component_index(i).first];
4209 return this->finite_element_output.shape_hessians[row][j];
4215 template <
int dim,
int spacedim>
4219 const unsigned int j,
4220 const unsigned int component)
const 4222 Assert (i < fe->dofs_per_cell,
4232 if (fe->get_nonzero_components(i)[component] ==
false)
4239 row = this->finite_element_output.shape_function_to_row_table[i * fe->n_components() + component];
4240 return this->finite_element_output.shape_hessians[row][j];
4245 template <
int dim,
int spacedim>
4249 const unsigned int j)
const 4251 Assert (i < fe->dofs_per_cell,
4255 Assert (fe->is_primitive (i),
4256 ExcShapeFunctionNotPrimitive(i));
4260 if (fe->is_primitive())
4261 return this->finite_element_output.shape_3rd_derivatives[i][j];
4273 row = this->finite_element_output.shape_function_to_row_table[i * fe->n_components() + fe->system_to_component_index(i).first];
4274 return this->finite_element_output.shape_3rd_derivatives[row][j];
4280 template <
int dim,
int spacedim>
4284 const unsigned int j,
4285 const unsigned int component)
const 4287 Assert (i < fe->dofs_per_cell,
4297 if (fe->get_nonzero_components(i)[component] ==
false)
4304 row = this->finite_element_output.shape_function_to_row_table[i * fe->n_components() + component];
4305 return this->finite_element_output.shape_3rd_derivatives[row][j];
4310 template <
int dim,
int spacedim>
4319 template <
int dim,
int spacedim>
4329 template <
int dim,
int spacedim>
4334 return this->update_flags;
4339 template <
int dim,
int spacedim>
4341 const std::vector<Point<spacedim> > &
4346 return this->mapping_output.quadrature_points;
4351 template <
int dim,
int spacedim>
4353 const std::vector<double> &
4358 return this->mapping_output.JxW_values;
4363 template <
int dim,
int spacedim>
4365 const std::vector<DerivativeForm<1,dim,spacedim> > &
4370 return this->mapping_output.jacobians;
4375 template <
int dim,
int spacedim>
4377 const std::vector<DerivativeForm<2,dim,spacedim> > &
4382 return this->mapping_output.jacobian_grads;
4387 template <
int dim,
int spacedim>
4394 return this->mapping_output.jacobian_pushed_forward_grads[i];
4399 template <
int dim,
int spacedim>
4401 const std::vector<Tensor<3,spacedim> > &
4406 return this->mapping_output.jacobian_pushed_forward_grads;
4411 template <
int dim,
int spacedim>
4418 return this->mapping_output.jacobian_2nd_derivatives[i];
4423 template <
int dim,
int spacedim>
4425 const std::vector<DerivativeForm<3,dim,spacedim> > &
4430 return this->mapping_output.jacobian_2nd_derivatives;
4435 template <
int dim,
int spacedim>
4442 return this->mapping_output.jacobian_pushed_forward_2nd_derivatives[i];
4447 template <
int dim,
int spacedim>
4449 const std::vector<Tensor<4,spacedim> > &
4454 return this->mapping_output.jacobian_pushed_forward_2nd_derivatives;
4459 template <
int dim,
int spacedim>
4466 return this->mapping_output.jacobian_3rd_derivatives[i];
4471 template <
int dim,
int spacedim>
4473 const std::vector<DerivativeForm<4,dim,spacedim> > &
4478 return this->mapping_output.jacobian_3rd_derivatives;
4483 template <
int dim,
int spacedim>
4490 return this->mapping_output.jacobian_pushed_forward_3rd_derivatives[i];
4495 template <
int dim,
int spacedim>
4497 const std::vector<Tensor<5,spacedim> > &
4502 return this->mapping_output.jacobian_pushed_forward_3rd_derivatives;
4506 template <
int dim,
int spacedim>
4508 const std::vector<DerivativeForm<1,spacedim,dim> > &
4513 return this->mapping_output.inverse_jacobians;
4518 template <
int dim,
int spacedim>
4525 Assert (i<this->mapping_output.quadrature_points.size(),
4526 ExcIndexRange(i, 0, this->mapping_output.quadrature_points.size()));
4528 return this->mapping_output.quadrature_points[i];
4534 template <
int dim,
int spacedim>
4541 Assert (i<this->mapping_output.JxW_values.size(),
4542 ExcIndexRange(i, 0, this->mapping_output.JxW_values.size()));
4544 return this->mapping_output.JxW_values[i];
4549 template <
int dim,
int spacedim>
4556 Assert (i<this->mapping_output.jacobians.size(),
4557 ExcIndexRange(i, 0, this->mapping_output.jacobians.size()));
4559 return this->mapping_output.jacobians[i];
4564 template <
int dim,
int spacedim>
4571 Assert (i<this->mapping_output.jacobian_grads.size(),
4572 ExcIndexRange(i, 0, this->mapping_output.jacobian_grads.size()));
4574 return this->mapping_output.jacobian_grads[i];
4579 template <
int dim,
int spacedim>
4586 Assert (i<this->mapping_output.inverse_jacobians.size(),
4587 ExcIndexRange(i, 0, this->mapping_output.inverse_jacobians.size()));
4589 return this->mapping_output.inverse_jacobians[i];
4593 template <
int dim,
int spacedim>
4600 Assert (i<this->mapping_output.normal_vectors.size(),
4601 ExcIndexRange(i, 0, this->mapping_output.normal_vectors.size()));
4603 return this->mapping_output.normal_vectors[i];
4611 template <
int dim,
int spacedim>
4621 template <
int dim,
int spacedim>
4633 template <
int dim,
int spacedim>
4638 return present_face_index;
4644 template <
int dim,
int spacedim>
4654 template <
int dim,
int spacedim>
4664 template <
int dim,
int spacedim>
4674 template <
int dim,
int spacedim>
4679 Assert (i<this->mapping_output.boundary_forms.size(),
4680 ExcIndexRange(i, 0, this->mapping_output.boundary_forms.size()));
4684 return this->mapping_output.boundary_forms[i];
4689 DEAL_II_NAMESPACE_CLOSE
void get_function_third_derivatives(const InputVector &fe_function, std::vector< typename ProductType< third_derivative_type, typename InputVector::value_type >::type > &third_derivatives) const
Transformed quadrature weights.
void get_function_gradients(const InputVector &fe_function, std::vector< Tensor< 1, spacedim, typename InputVector::value_type > > &gradients) const
std_cxx11::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > fe_data
value_type value(const unsigned int shape_function, const unsigned int q_point) const
static ::ExceptionBase & ExcCannotInitializeField()
const DerivativeForm< 1, dim, spacedim > & jacobian(const unsigned int quadrature_point) const
const Tensor< 3, spacedim > & jacobian_pushed_forward_grad(const unsigned int quadrature_point) const
CellSimilarity::Similarity cell_similarity
::Tensor< 1, spacedim > divergence_type
Cache(const FEValuesBase< dim, spacedim > &fe_values)
const Tensor< 5, spacedim > & jacobian_pushed_forward_3rd_derivative(const unsigned int quadrature_point) const
unsigned int present_face_index
const FEValues< dim, spacedim > & get_present_fe_values() const
std::vector< ShapeFunctionData > shape_function_data
static ::ExceptionBase & ExcAccessToUninitializedField()
value_type value(const unsigned int shape_function, const unsigned int q_point) const
std::vector<::FEValuesViews::Scalar< dim, spacedim > > scalars
const unsigned int dofs_per_cell
const unsigned int component
void get_function_values(const InputVector &fe_function, std::vector< typename ProductType< value_type, typename InputVector::value_type >::type > &values) const
gradient_type gradient(const unsigned int shape_function, const unsigned int q_point) const
const Tensor< 3, spacedim > & shape_3rd_derivative(const unsigned int function_no, const unsigned int point_no) const
FEValuesBase(const unsigned int n_q_points, const unsigned int dofs_per_cell, const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe)
::Tensor< 3, spacedim > hessian_type
const Quadrature< dim-1 > quadrature
void get_function_third_derivatives(const InputVector &fe_function, std::vector< typename ProductType< third_derivative_type, typename InputVector::value_type >::type > &third_derivatives) const
int single_nonzero_component
void get_function_values(const InputVector &fe_function, std::vector< typename InputVector::value_type > &values) const
::internal::CurlType< spacedim >::type curl_type
Outer normal vector, not normalized.
const std::vector< DerivativeForm< 4, dim, spacedim > > & get_jacobian_3rd_derivatives() const
const FiniteElement< dim, spacedim > & get_fe() const
static ::ExceptionBase & ExcFEDontMatch()
Scalar & operator=(const Scalar< dim, spacedim > &)
void get_function_divergences(const InputVector &fe_function, std::vector< typename ProductType< divergence_type, typename InputVector::value_type >::type > &divergences) const
const unsigned int first_tensor_component
void get_function_laplacians(const InputVector &fe_function, std::vector< typename ProductType< value_type, typename InputVector::value_type >::type > &laplacians) const
UpdateFlags get_update_flags() const
third_derivative_type third_derivative(const unsigned int shape_function, const unsigned int q_point) const
Transformed quadrature points.
::Tensor< 1, spacedim > gradient_type
::SymmetricTensor< 2, spacedim > value_type
void transform(std::vector< Tensor< 1, spacedim > > &transformed, const std::vector< Tensor< 1, dim > > &original, MappingType mapping) const 1
unsigned int row_index[spacedim]
const DerivativeForm< 3, dim, spacedim > & jacobian_2nd_derivative(const unsigned int quadrature_point) const
void check_cell_similarity(const typename Triangulation< dim, spacedim >::cell_iterator &cell)
static const unsigned int space_dimension
::internal::FEValuesViews::Cache< dim, spacedim > fe_values_views_cache
const FEValuesViews::Scalar< dim, spacedim > & operator[](const FEValuesExtractors::Scalar &scalar) const
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
const DerivativeForm< 2, dim, spacedim > & jacobian_grad(const unsigned int quadrature_point) const
const std::vector< DerivativeForm< 1, dim, spacedim > > & get_jacobians() const
void get_function_hessians(const InputVector &fe_function, std::vector< Tensor< 2, spacedim, typename InputVector::value_type > > &hessians) const
const Mapping< dim, spacedim > & get_mapping() const
static ::ExceptionBase & ExcShapeFunctionNotPrimitive(int arg1)
const DerivativeForm< 4, dim, spacedim > & jacobian_3rd_derivative(const unsigned int quadrature_point) const
double shape_value_component(const unsigned int function_no, const unsigned int point_no, const unsigned int component) const
static ::ExceptionBase & ExcInvalidUpdateFlag()
const std::vector< Tensor< 5, spacedim > > & get_jacobian_pushed_forward_3rd_derivatives() const
void get_function_laplacians(const InputVector &fe_function, std::vector< typename ProductType< value_type, typename InputVector::value_type >::type > &laplacians) const
const std::vector< Point< spacedim > > & get_quadrature_points() const
const Point< spacedim > & quadrature_point(const unsigned int q) const
void get_function_gradients(const InputVector &fe_function, std::vector< typename ProductType< gradient_type, typename InputVector::value_type >::type > &gradients) const
const std::vector< DerivativeForm< 2, dim, spacedim > > & get_jacobian_grads() const
void get_function_hessians(const InputVector &fe_function, std::vector< typename ProductType< hessian_type, typename InputVector::value_type >::type > &hessians) const
void get_function_laplacians(const InputVector &fe_function, std::vector< typename InputVector::value_type > &laplacians) const
const SmartPointer< const FiniteElement< dim, spacedim >, FEValuesBase< dim, spacedim > > fe
static ::ExceptionBase & ExcMessage(std::string arg1)
const unsigned int first_tensor_component
::Tensor< 1, spacedim > value_type
#define DeclException1(Exception1, type1, outsequence)
std_cxx11::unique_ptr< const CellIteratorBase > present_cell
const double & shape_value(const unsigned int function_no, const unsigned int point_no) const
CellSimilarity::Similarity get_cell_similarity() const
Third derivatives of shape functions.
void get_function_third_derivatives(const InputVector &fe_function, std::vector< Tensor< 3, spacedim, typename InputVector::value_type > > &third_derivatives) const
UpdateFlags compute_update_flags(const UpdateFlags update_flags) const
std::vector< ShapeFunctionData > shape_function_data
unsigned int get_face_index() const
#define Assert(cond, exc)
::Tensor< 2, spacedim > value_type
const std::vector< Tensor< 3, spacedim > > & get_jacobian_pushed_forward_grads() const
const DerivativeForm< 1, spacedim, dim > & inverse_jacobian(const unsigned int quadrature_point) const
Abstract base class for mapping classes.
const Quadrature< dim-1 > & get_quadrature() const
std::vector< ShapeFunctionData > shape_function_data
const Quadrature< dim > quadrature
const unsigned int first_vector_component
#define DeclException0(Exception0)
::internal::FEValues::FiniteElementRelatedData< dim, spacedim > finite_element_output
const std::vector< DerivativeForm< 1, spacedim, dim > > & get_inverse_jacobians() const
void invalidate_present_cell()
const SmartPointer< const Mapping< dim, spacedim >, FEValuesBase< dim, spacedim > > mapping
Vector & operator=(const Vector< dim, spacedim > &)
::Tensor< 1, spacedim > divergence_type
const FEFaceValues< dim, spacedim > & get_present_fe_values() const
symmetric_gradient_type symmetric_gradient(const unsigned int shape_function, const unsigned int q_point) const
bool is_nonzero_shape_function_component[spacedim]
void get_function_symmetric_gradients(const InputVector &fe_function, std::vector< typename ProductType< symmetric_gradient_type, typename InputVector::value_type >::type > &symmetric_gradients) const
Second derivatives of shape functions.
static ::ExceptionBase & ExcFENotPrimitive()
static ::ExceptionBase & ExcAccessToUninitializedField(char *arg1)
Gradient of volume element.
const std::vector< Tensor< 1, spacedim > > & get_all_normal_vectors() const
const Tensor< 2, spacedim > & shape_hessian(const unsigned int function_no, const unsigned int point_no) const
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
const Tensor< 1, spacedim > & boundary_form(const unsigned int i) const
static TableIndices< rank_ > unrolled_to_component_indices(const unsigned int i)
divergence_type divergence(const unsigned int shape_function, const unsigned int q_point) const
void get_function_hessians(const InputVector &fe_function, std::vector< typename ProductType< hessian_type, typename InputVector::value_type >::type > &hessians) const
const Quadrature< dim > & get_quadrature() const
const unsigned int n_quadrature_points
const std::vector< DerivativeForm< 3, dim, spacedim > > & get_jacobian_2nd_derivatives() const
const SmartPointer< const FEValuesBase< dim, spacedim > > fe_values
void get_function_curls(const InputVector &fe_function, std::vector< typename ProductType< curl_type, typename InputVector::value_type >::type > &curls) const
static const unsigned int dimension
const std::vector< Tensor< 4, spacedim > > & get_jacobian_pushed_forward_2nd_derivatives() const
::Tensor< 2, spacedim > hessian_type
boost::signals2::connection tria_listener_mesh_transform
int single_nonzero_component
const SmartPointer< const FEValuesBase< dim, spacedim > > fe_values
const std::vector< double > & get_JxW_values() const
::SymmetricTensor< 2, spacedim > symmetric_gradient_type
double JxW(const unsigned int quadrature_point) const
Shape function gradients.
void maybe_invalidate_previous_present_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell)
::Tensor< 2, spacedim > gradient_type
Tensor< 3, spacedim > shape_3rd_derivative_component(const unsigned int function_no, const unsigned int point_no, const unsigned int component) const
std_cxx11::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > mapping_data
const FESubfaceValues< dim, spacedim > & get_present_fe_values() const
hessian_type hessian(const unsigned int shape_function, const unsigned int q_point) const
const SmartPointer< const FEValuesBase< dim, spacedim > > fe_values
const SmartPointer< const FEValuesBase< dim, spacedim > > fe_values
curl_type curl(const unsigned int shape_function, const unsigned int q_point) const
static ::ExceptionBase & ExcNotImplemented()
bool is_nonzero_shape_function_component
const Triangulation< dim, spacedim >::cell_iterator get_cell() const
void get_function_values(const InputVector &fe_function, std::vector< typename ProductType< value_type, typename InputVector::value_type >::type > &values) const
boost::signals2::connection tria_listener_refinement
const Tensor< 4, spacedim > & jacobian_pushed_forward_2nd_derivative(const unsigned int quadrature_point) const
::internal::FEValues::MappingRelatedData< dim, spacedim > mapping_output
gradient_type gradient(const unsigned int shape_function, const unsigned int q_point) const
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
third_derivative_type third_derivative(const unsigned int shape_function, const unsigned int q_point) const
void get_function_gradients(const InputVector &fe_function, std::vector< typename ProductType< gradient_type, typename InputVector::value_type >::type > &gradients) const
hessian_type hessian(const unsigned int shape_function, const unsigned int q_point) const
std::vector< Point< spacedim > > get_normal_vectors() const 1
Tensor< 2, spacedim > shape_hessian_component(const unsigned int function_no, const unsigned int point_no, const unsigned int component) const
int single_nonzero_component
::Tensor< 3, spacedim > third_derivative_type
std::size_t memory_consumption() const
const Tensor< 1, spacedim > & shape_grad(const unsigned int function_no, const unsigned int quadrature_point) const
std::vector< ShapeFunctionData > shape_function_data
static ::ExceptionBase & ExcInternalError()
::Tensor< 4, spacedim > third_derivative_type
const Tensor< 1, spacedim > & normal_vector(const unsigned int i) const
Tensor< 1, spacedim > shape_grad_component(const unsigned int function_no, const unsigned int point_no, const unsigned int component) const