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
54 template <
typename Number>
class Vector;
141 template <
int dim,
int spacedim=dim>
211 const unsigned int component);
233 value (
const unsigned int shape_function,
234 const unsigned int q_point)
const;
247 gradient (
const unsigned int shape_function,
248 const unsigned int q_point)
const;
261 hessian (
const unsigned int shape_function,
262 const unsigned int q_point)
const;
274 third_derivative_type
275 third_derivative (
const unsigned int shape_function,
276 const unsigned int q_point)
const;
295 template <
class InputVector>
296 void get_function_values (
const InputVector &fe_function,
316 template <
class InputVector>
317 void get_function_gradients (
const InputVector &fe_function,
337 template <
class InputVector>
338 void get_function_hessians (
const InputVector &fe_function,
359 template <
class InputVector>
360 void get_function_laplacians (
const InputVector &fe_function,
381 template <
class InputVector>
382 void get_function_third_derivatives (
const InputVector &fe_function,
383 std::vector<
typename ProductType<third_derivative_type,
384 typename InputVector::value_type>::type> &third_derivatives)
const;
435 template <
int dim,
int spacedim=dim>
481 typedef typename ::internal::CurlType<spacedim>::type
curl_type;
511 bool is_nonzero_shape_function_component[spacedim];
522 unsigned int row_index[spacedim];
533 unsigned int single_nonzero_component_index;
550 const unsigned int first_vector_component);
575 value (
const unsigned int shape_function,
576 const unsigned int q_point)
const;
592 gradient (
const unsigned int shape_function,
593 const unsigned int q_point)
const;
610 symmetric_gradient_type
611 symmetric_gradient (
const unsigned int shape_function,
612 const unsigned int q_point)
const;
625 divergence (
const unsigned int shape_function,
626 const unsigned int q_point)
const;
649 curl (
const unsigned int shape_function,
650 const unsigned int q_point)
const;
663 hessian (
const unsigned int shape_function,
664 const unsigned int q_point)
const;
676 third_derivative_type
677 third_derivative (
const unsigned int shape_function,
678 const unsigned int q_point)
const;
697 template <
class InputVector>
698 void get_function_values (
const InputVector &fe_function,
718 template <
class InputVector>
719 void get_function_gradients (
const InputVector &fe_function,
745 template <
class InputVector>
747 get_function_symmetric_gradients (
const InputVector &fe_function,
768 template <
class InputVector>
769 void get_function_divergences (
const InputVector &fe_function,
790 template <
class InputVector>
791 void get_function_curls (
const InputVector &fe_function,
811 template <
class InputVector>
812 void get_function_hessians (
const InputVector &fe_function,
833 template <
class InputVector>
834 void get_function_laplacians (
const InputVector &fe_function,
855 template <
class InputVector>
856 void get_function_third_derivatives (
const InputVector &fe_function,
857 std::vector<
typename ProductType<third_derivative_type,
858 typename InputVector::value_type>::type> &third_derivatives)
const;
879 template <
int rank,
int dim,
int spacedim = dim>
906 template <
int dim,
int spacedim>
933 struct ShapeFunctionData
943 bool is_nonzero_shape_function_component[value_type::n_independent_components];
954 unsigned int row_index[value_type::n_independent_components];
965 unsigned int single_nonzero_component_index;
983 const unsigned int first_tensor_component);
1009 value (
const unsigned int shape_function,
1010 const unsigned int q_point)
const;
1027 divergence (
const unsigned int shape_function,
1028 const unsigned int q_point)
const;
1047 template <
class InputVector>
1048 void get_function_values (
const InputVector &fe_function,
1072 template <
class InputVector>
1073 void get_function_divergences (
const InputVector &fe_function,
1095 template <
int rank,
int dim,
int spacedim = dim>
1117 template <
int dim,
int spacedim>
1137 struct ShapeFunctionData
1147 bool is_nonzero_shape_function_component[value_type::n_independent_components];
1158 unsigned int row_index[value_type::n_independent_components];
1169 unsigned int single_nonzero_component_index;
1188 const unsigned int first_tensor_component);
1214 value (
const unsigned int shape_function,
1215 const unsigned int q_point)
const;
1231 divergence (
const unsigned int shape_function,
1232 const unsigned int q_point)
const;
1251 template <
class InputVector>
1252 void get_function_values (
const InputVector &fe_function,
1277 template <
class InputVector>
1278 void get_function_divergences (
const InputVector &fe_function,
1313 template <
int dim,
int spacedim>
1320 std::vector<::FEValuesViews::Scalar<dim,spacedim> >
scalars;
1321 std::vector<::FEValuesViews::Vector<dim,spacedim> > vectors;
1322 std::vector<::FEValuesViews::SymmetricTensor<2,dim,spacedim> >
1323 symmetric_second_order_tensors;
1324 std::vector<::FEValuesViews::Tensor<2,dim,spacedim> >
1325 second_order_tensors;
1439 template <
int dim,
int spacedim>
1447 static const unsigned int dimension = dim;
1452 static const unsigned int space_dimension = spacedim;
1475 const unsigned int dofs_per_cell,
1510 const double &shape_value (
const unsigned int function_no,
1511 const unsigned int point_no)
const;
1533 double shape_value_component (
const unsigned int function_no,
1534 const unsigned int point_no,
1535 const unsigned int component)
const;
1563 shape_grad (
const unsigned int function_no,
1564 const unsigned int quadrature_point)
const;
1583 shape_grad_component (
const unsigned int function_no,
1584 const unsigned int point_no,
1585 const unsigned int component)
const;
1607 shape_hessian (
const unsigned int function_no,
1608 const unsigned int point_no)
const;
1627 shape_hessian_component (
const unsigned int function_no,
1628 const unsigned int point_no,
1629 const unsigned int component)
const;
1651 shape_3rd_derivative (
const unsigned int function_no,
1652 const unsigned int point_no)
const;
1671 shape_3rd_derivative_component (
const unsigned int function_no,
1672 const unsigned int point_no,
1673 const unsigned int component)
const;
1717 template <
class InputVector>
1718 void get_function_values (
const InputVector &fe_function,
1719 std::vector<typename InputVector::value_type> &values)
const;
1734 template <
class InputVector>
1735 void get_function_values (
const InputVector &fe_function,
1756 template <
class InputVector>
1757 void get_function_values (
const InputVector &fe_function,
1758 const VectorSlice<
const std::vector<types::global_dof_index> > &indices,
1759 std::vector<typename InputVector::value_type> &values)
const;
1782 template <
class InputVector>
1783 void get_function_values (
const InputVector &fe_function,
1784 const VectorSlice<
const std::vector<types::global_dof_index> > &indices,
1818 template <
class InputVector>
1819 void get_function_values (
const InputVector &fe_function,
1820 const VectorSlice<
const std::vector<types::global_dof_index> > &indices,
1821 VectorSlice<std::vector<std::vector<typename InputVector::value_type> > > values,
1822 const bool quadrature_points_fastest)
const;
1866 template <
class InputVector>
1867 void get_function_gradients (
const InputVector &fe_function,
1886 template <
class InputVector>
1887 void get_function_gradients (
const InputVector &fe_function,
1896 template <
class InputVector>
1897 void get_function_gradients (
const InputVector &fe_function,
1898 const VectorSlice<
const std::vector<types::global_dof_index> > &indices,
1907 template <
class InputVector>
1908 void get_function_gradients (
const InputVector &fe_function,
1909 const VectorSlice<
const std::vector<types::global_dof_index> > &indices,
1911 bool quadrature_points_fastest =
false)
const;
1956 template <
class InputVector>
1958 get_function_hessians (
const InputVector &fe_function,
1978 template <
class InputVector>
1980 get_function_hessians (
const InputVector &fe_function,
1982 bool quadrature_points_fastest =
false)
const;
1988 template <
class InputVector>
1989 void get_function_hessians (
1990 const InputVector &fe_function,
1991 const VectorSlice<
const std::vector<types::global_dof_index> > &indices,
2000 template <
class InputVector>
2001 void get_function_hessians (
2002 const InputVector &fe_function,
2003 const VectorSlice<
const std::vector<types::global_dof_index> > &indices,
2005 bool quadrature_points_fastest =
false)
const;
2051 template <
class InputVector>
2053 get_function_laplacians (
const InputVector &fe_function,
2054 std::vector<typename InputVector::value_type> &laplacians)
const;
2075 template <
class InputVector>
2077 get_function_laplacians (
const InputVector &fe_function,
2086 template <
class InputVector>
2087 void get_function_laplacians (
2088 const InputVector &fe_function,
2089 const VectorSlice<
const std::vector<types::global_dof_index> > &indices,
2090 std::vector<typename InputVector::value_type> &laplacians)
const;
2098 template <
class InputVector>
2099 void get_function_laplacians (
2100 const InputVector &fe_function,
2101 const VectorSlice<
const std::vector<types::global_dof_index> > &indices,
2110 template <
class InputVector>
2111 void get_function_laplacians (
2112 const InputVector &fe_function,
2113 const VectorSlice<
const std::vector<types::global_dof_index> > &indices,
2114 std::vector<std::vector<typename InputVector::value_type> > &laplacians,
2115 bool quadrature_points_fastest =
false)
const;
2161 template <
class InputVector>
2163 get_function_third_derivatives (
const InputVector &fe_function,
2184 template <
class InputVector>
2186 get_function_third_derivatives (
const InputVector &fe_function,
2188 bool quadrature_points_fastest =
false)
const;
2194 template <
class InputVector>
2195 void get_function_third_derivatives (
2196 const InputVector &fe_function,
2197 const VectorSlice<
const std::vector<types::global_dof_index> > &indices,
2206 template <
class InputVector>
2207 void get_function_third_derivatives (
2208 const InputVector &fe_function,
2209 const VectorSlice<
const std::vector<types::global_dof_index> > &indices,
2211 bool quadrature_points_fastest =
false)
const;
2223 quadrature_point (
const unsigned int q)
const;
2230 const std::vector<Point<spacedim> > &get_quadrature_points ()
const;
2247 double JxW (
const unsigned int quadrature_point)
const;
2252 const std::vector<double> &get_JxW_values ()
const;
2268 const std::vector<DerivativeForm<1,dim,spacedim> > &get_jacobians ()
const;
2285 const std::vector<DerivativeForm<2,dim,spacedim> > &get_jacobian_grads ()
const;
2295 const Tensor<3,spacedim> &jacobian_pushed_forward_grad (
const unsigned int quadrature_point)
const;
2303 const std::vector<Tensor<3,spacedim> > &get_jacobian_pushed_forward_grads ()
const;
2320 const std::vector<DerivativeForm<3,dim,spacedim> > &get_jacobian_2nd_derivatives ()
const;
2331 const Tensor<4,spacedim> &jacobian_pushed_forward_2nd_derivative (
const unsigned int quadrature_point)
const;
2339 const std::vector<Tensor<4,spacedim> > &get_jacobian_pushed_forward_2nd_derivatives ()
const;
2357 const std::vector<DerivativeForm<4,dim,spacedim> > &get_jacobian_3rd_derivatives ()
const;
2368 const Tensor<5,spacedim> &jacobian_pushed_forward_3rd_derivative (
const unsigned int quadrature_point)
const;
2376 const std::vector<Tensor<5,spacedim> > &get_jacobian_pushed_forward_3rd_derivatives ()
const;
2392 const std::vector<DerivativeForm<1,spacedim,dim> > &get_inverse_jacobians ()
const;
2424 const std::vector<Tensor<1,spacedim> > &get_all_normal_vectors ()
const;
2444 void transform (
std::vector<
Tensor<1,spacedim> > &transformed,
2445 const
std::vector<Tensor<1,dim> > &original,
2507 const
Mapping<dim,spacedim> &get_mapping () const;
2522 const typename
Triangulation<dim,spacedim>::cell_iterator get_cell () const;
2535 std::
size_t memory_consumption () const;
2548 "
object for which this kind of information has not been computed. What "
2549 "information these objects compute is determined by the update_* flags you "
2550 "pass to the constructor. Here, the operation you are attempting requires "
2553 << "> flag to be set, but it was apparently not specified upon construction.");
2559 DeclException0 (ExcCannotInitializeField);
2565 DeclException0 (ExcInvalidUpdateFlag);
2571 DeclException0 (ExcFEDontMatch);
2579 << "The shape function with index " << arg1
2580 << " is not primitive, i.e. it is vector-valued and "
2581 << "has more than one non-zero vector component. This "
2582 << "function cannot be called for these shape functions. "
2583 << "Maybe you want to use the same function with the "
2584 << "_component suffix?");
2590 DeclException0 (ExcFENotPrimitive);
2621 class CellIteratorBase;
2635 std_cxx11::unique_ptr<const CellIteratorBase> present_cell;
2644 boost::signals2::connection tria_listener;
2651 void invalidate_present_cell ();
2663 maybe_invalidate_previous_present_cell (const typename
Triangulation<dim,spacedim>::cell_iterator &cell);
2732 check_cell_similarity (const typename
Triangulation<dim,spacedim>::cell_iterator &cell);
2745 FEValuesBase &operator= (const FEValuesBase &);
2758 template <
int,
int,
int> friend class
FEValuesViews::SymmetricTensor;
2759 template <
int,
int,
int> friend class
FEValuesViews::Tensor;
2774 template <
int dim,
int spacedim=dim>
2775 class
FEValues : public FEValuesBase<dim,spacedim>
2782 static const unsigned int integral_dimension = dim;
2808 template <
template <
int,
int>
class DoFHandlerType,
bool level_dof_access>
2836 std::size_t memory_consumption ()
const;
2884 template <
int dim,
int spacedim=dim>
2892 static const unsigned int integral_dimension = dim-1;
2906 const unsigned int dofs_per_cell,
2927 const std::vector<Tensor<1,spacedim> > &get_boundary_forms ()
const;
2933 unsigned int get_face_index()
const;
2939 const Quadrature<dim-1> & get_quadrature ()
const;
2945 std::size_t memory_consumption ()
const;
2977 template <
int dim,
int spacedim=dim>
2985 static const unsigned int dimension = dim;
2987 static const unsigned int space_dimension = spacedim;
2993 static const unsigned int integral_dimension = dim-1;
3017 template <
template <
int,
int>
class DoFHandlerType,
bool level_dof_access>
3019 const unsigned int face_no);
3035 const unsigned int face_no);
3064 void do_reinit (
const unsigned int face_no);
3085 template <
int dim,
int spacedim=dim>
3092 static const unsigned int dimension = dim;
3097 static const unsigned int space_dimension = spacedim;
3103 static const unsigned int integral_dimension = dim-1;
3129 template <
template <
int,
int>
class DoFHandlerType,
bool level_dof_access>
3131 const unsigned int face_no,
3132 const unsigned int subface_no);
3148 const unsigned int face_no,
3149 const unsigned int subface_no);
3193 void do_reinit (
const unsigned int face_no,
3194 const unsigned int subface_no);
3205 template <
int dim,
int spacedim>
3207 typename Scalar<dim,spacedim>::value_type
3208 Scalar<dim,spacedim>::value (
const unsigned int shape_function,
3209 const unsigned int q_point)
const 3212 Assert (shape_function < fe_values.fe->dofs_per_cell,
3213 ExcIndexRange (shape_function, 0, fe_values.fe->dofs_per_cell));
3215 typename FVB::ExcAccessToUninitializedField(
"update_values"));
3220 if (shape_function_data[shape_function].is_nonzero_shape_function_component)
3221 return fe_values.finite_element_output.shape_values(shape_function_data[shape_function]
3231 template <
int dim,
int spacedim>
3233 typename Scalar<dim,spacedim>::gradient_type
3234 Scalar<dim,spacedim>::gradient (
const unsigned int shape_function,
3235 const unsigned int q_point)
const 3238 Assert (shape_function < fe_values.fe->dofs_per_cell,
3239 ExcIndexRange (shape_function, 0, fe_values.fe->dofs_per_cell));
3241 typename FVB::ExcAccessToUninitializedField(
"update_gradients"));
3249 if (shape_function_data[shape_function].is_nonzero_shape_function_component)
3250 return fe_values.finite_element_output.shape_gradients[shape_function_data[shape_function]
3251 .row_index][q_point];
3253 return gradient_type();
3258 template <
int dim,
int spacedim>
3260 typename Scalar<dim,spacedim>::hessian_type
3261 Scalar<dim,spacedim>::hessian (
const unsigned int shape_function,
3262 const unsigned int q_point)
const 3265 Assert (shape_function < fe_values.fe->dofs_per_cell,
3266 ExcIndexRange (shape_function, 0, fe_values.fe->dofs_per_cell));
3268 typename FVB::ExcAccessToUninitializedField(
"update_hessians"));
3276 if (shape_function_data[shape_function].is_nonzero_shape_function_component)
3277 return fe_values.finite_element_output.shape_hessians[shape_function_data[shape_function].row_index][q_point];
3279 return hessian_type();
3284 template <
int dim,
int spacedim>
3286 typename Scalar<dim,spacedim>::third_derivative_type
3287 Scalar<dim,spacedim>::third_derivative (
const unsigned int shape_function,
3288 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));
3294 typename FVB::ExcAccessToUninitializedField(
"update_3rd_derivatives"));
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 3317 Assert (shape_function < fe_values.fe->dofs_per_cell,
3318 ExcIndexRange (shape_function, 0, fe_values.fe->dofs_per_cell));
3320 typename FVB::ExcAccessToUninitializedField(
"update_values"));
3324 const int snc = shape_function_data[shape_function].single_nonzero_component;
3326 return value_type();
3329 value_type return_value;
3330 return_value[shape_function_data[shape_function].single_nonzero_component_index]
3331 = fe_values.finite_element_output.shape_values(snc,q_point);
3332 return return_value;
3336 value_type return_value;
3337 for (
unsigned int d=0; d<dim; ++d)
3338 if (shape_function_data[shape_function].is_nonzero_shape_function_component[d])
3340 = fe_values.finite_element_output.shape_values(shape_function_data[shape_function].row_index[d],q_point);
3342 return return_value;
3348 template <
int dim,
int spacedim>
3352 const unsigned int q_point)
const 3355 Assert (shape_function < fe_values.fe->dofs_per_cell,
3356 ExcIndexRange (shape_function, 0, fe_values.fe->dofs_per_cell));
3358 typename FVB::ExcAccessToUninitializedField(
"update_gradients"));
3362 const int snc = shape_function_data[shape_function].single_nonzero_component;
3364 return gradient_type();
3367 gradient_type return_value;
3368 return_value[shape_function_data[shape_function].single_nonzero_component_index]
3369 = fe_values.finite_element_output.shape_gradients[snc][q_point];
3370 return return_value;
3374 gradient_type return_value;
3375 for (
unsigned int d=0; d<dim; ++d)
3376 if (shape_function_data[shape_function].is_nonzero_shape_function_component[d])
3378 = fe_values.finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[d]][q_point];
3380 return return_value;
3386 template <
int dim,
int spacedim>
3390 const unsigned int q_point)
const 3395 Assert (shape_function < fe_values.fe->dofs_per_cell,
3396 ExcIndexRange (shape_function, 0, fe_values.fe->dofs_per_cell));
3398 typename FVB::ExcAccessToUninitializedField(
"update_gradients"));
3402 const int snc = shape_function_data[shape_function].single_nonzero_component;
3404 return divergence_type();
3407 fe_values.finite_element_output.shape_gradients[snc][q_point][shape_function_data[shape_function].single_nonzero_component_index];
3410 divergence_type return_value = 0;
3411 for (
unsigned int d=0; d<dim; ++d)
3412 if (shape_function_data[shape_function].is_nonzero_shape_function_component[d])
3414 += fe_values.finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[d]][q_point][d];
3416 return return_value;
3422 template <
int dim,
int spacedim>
3430 Assert (shape_function < fe_values.fe->dofs_per_cell,
3431 ExcIndexRange (shape_function, 0, fe_values.fe->dofs_per_cell));
3433 typename FVB::ExcAccessToUninitializedField(
"update_gradients"));
3435 const int snc = shape_function_data[shape_function].single_nonzero_component;
3438 return curl_type ();
3445 Assert (
false,
ExcMessage(
"Computing the curl in 1d is not a useful operation"));
3446 return curl_type ();
3453 curl_type return_value;
3459 if (shape_function_data[shape_function].single_nonzero_component_index == 0)
3460 return_value[0] = -1.0 * fe_values.finite_element_output.shape_gradients[snc][q_point][1];
3462 return_value[0] = fe_values.finite_element_output.shape_gradients[snc][q_point][0];
3464 return return_value;
3469 curl_type return_value;
3471 return_value[0] = 0.0;
3473 if (shape_function_data[shape_function].is_nonzero_shape_function_component[0])
3475 -= fe_values.finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[0]][q_point][1];
3477 if (shape_function_data[shape_function].is_nonzero_shape_function_component[1])
3479 += fe_values.finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[1]][q_point][0];
3481 return return_value;
3489 curl_type return_value;
3491 switch (shape_function_data[shape_function].single_nonzero_component_index)
3495 return_value[0] = 0;
3496 return_value[1] = fe_values.finite_element_output.shape_gradients[snc][q_point][2];
3497 return_value[2] = -1.0 * fe_values.finite_element_output.shape_gradients[snc][q_point][1];
3498 return return_value;
3503 return_value[0] = -1.0 * fe_values.finite_element_output.shape_gradients[snc][q_point][2];
3504 return_value[1] = 0;
3505 return_value[2] = fe_values.finite_element_output.shape_gradients[snc][q_point][0];
3506 return return_value;
3511 return_value[0] = fe_values.finite_element_output.shape_gradients[snc][q_point][1];
3512 return_value[1] = -1.0 * fe_values.finite_element_output.shape_gradients[snc][q_point][0];
3513 return_value[2] = 0;
3514 return return_value;
3521 curl_type return_value;
3523 for (
unsigned int i = 0; i < dim; ++i)
3524 return_value[i] = 0.0;
3526 if (shape_function_data[shape_function].is_nonzero_shape_function_component[0])
3529 += fe_values.finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[0]][q_point][2];
3531 -= fe_values.finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[0]][q_point][1];
3534 if (shape_function_data[shape_function].is_nonzero_shape_function_component[1])
3537 -= fe_values.finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[1]][q_point][2];
3539 += fe_values.finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[1]][q_point][0];
3542 if (shape_function_data[shape_function].is_nonzero_shape_function_component[2])
3545 += fe_values.finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[2]][q_point][1];
3547 -= fe_values.finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[2]][q_point][0];
3550 return return_value;
3555 Assert (
false, ExcInternalError());
3559 template <
int dim,
int spacedim>
3563 const unsigned int q_point)
const 3568 Assert (shape_function < fe_values.fe->dofs_per_cell,
3569 ExcIndexRange (shape_function, 0, fe_values.fe->dofs_per_cell));
3571 typename FVB::ExcAccessToUninitializedField(
"update_hessians"));
3575 const int snc = shape_function_data[shape_function].single_nonzero_component;
3577 return hessian_type();
3580 hessian_type return_value;
3581 return_value[shape_function_data[shape_function].single_nonzero_component_index]
3582 = fe_values.finite_element_output.shape_hessians[snc][q_point];
3583 return return_value;
3587 hessian_type return_value;
3588 for (
unsigned int d=0; d<dim; ++d)
3589 if (shape_function_data[shape_function].is_nonzero_shape_function_component[d])
3591 = fe_values.finite_element_output.shape_hessians[shape_function_data[shape_function].row_index[d]][q_point];
3593 return return_value;
3597 template <
int dim,
int spacedim>
3601 const unsigned int q_point)
const 3606 Assert (shape_function < fe_values.fe->dofs_per_cell,
3607 ExcIndexRange (shape_function, 0, fe_values.fe->dofs_per_cell));
3609 typename FVB::ExcAccessToUninitializedField(
"update_3rd_derivatives"));
3613 const int snc = shape_function_data[shape_function].single_nonzero_component;
3615 return third_derivative_type();
3618 third_derivative_type return_value;
3619 return_value[shape_function_data[shape_function].single_nonzero_component_index]
3620 = fe_values.finite_element_output.shape_3rd_derivatives[snc][q_point];
3621 return return_value;
3625 third_derivative_type return_value;
3626 for (
unsigned int d=0; d<dim; ++d)
3627 if (shape_function_data[shape_function].is_nonzero_shape_function_component[d])
3629 = fe_values.finite_element_output.shape_3rd_derivatives[shape_function_data[shape_function].row_index[d]][q_point];
3631 return return_value;
3643 ::SymmetricTensor<2,1>
3644 symmetrize_single_row (
const unsigned int n,
3645 const ::Tensor<1,1> &t)
3647 Assert (n < 1, ExcIndexRange (n, 0, 1));
3650 const double array[1] = { t[0] };
3651 return ::SymmetricTensor<2,1>(array);
3656 ::SymmetricTensor<2,2>
3657 symmetrize_single_row (
const unsigned int n,
3658 const ::Tensor<1,2> &t)
3664 const double array[3] = { t[0], 0, t[1]/2 };
3665 return ::SymmetricTensor<2,2>(array);
3669 const double array[3] = { 0, t[1], t[0]/2 };
3670 return ::SymmetricTensor<2,2>(array);
3674 Assert (
false, ExcIndexRange (n, 0, 2));
3675 return ::SymmetricTensor<2,2>();
3682 ::SymmetricTensor<2,3>
3683 symmetrize_single_row (
const unsigned int n,
3684 const ::Tensor<1,3> &t)
3690 const double array[6] = { t[0], 0, 0, t[1]/2, t[2]/2, 0 };
3691 return ::SymmetricTensor<2,3>(array);
3695 const double array[6] = { 0, t[1], 0, t[0]/2, 0, t[2]/2 };
3696 return ::SymmetricTensor<2,3>(array);
3700 const double array[6] = { 0, 0, t[2], 0, t[0]/2, t[1]/2 };
3701 return ::SymmetricTensor<2,3>(array);
3705 Assert (
false, ExcIndexRange (n, 0, 3));
3706 return ::SymmetricTensor<2,3>();
3713 template <
int dim,
int spacedim>
3717 const unsigned int q_point)
const 3720 Assert (shape_function < fe_values.fe->dofs_per_cell,
3721 ExcIndexRange (shape_function, 0, fe_values.fe->dofs_per_cell));
3723 typename FVB::ExcAccessToUninitializedField(
"update_gradients"));
3727 const int snc = shape_function_data[shape_function].single_nonzero_component;
3729 return symmetric_gradient_type();
3731 return symmetrize_single_row (shape_function_data[shape_function].single_nonzero_component_index,
3732 fe_values.finite_element_output.shape_gradients[snc][q_point]);
3735 gradient_type return_value;
3736 for (
unsigned int d=0; d<dim; ++d)
3737 if (shape_function_data[shape_function].is_nonzero_shape_function_component[d])
3739 = fe_values.finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[d]][q_point];
3741 return symmetrize(return_value);
3747 template <
int dim,
int spacedim>
3751 const unsigned int q_point)
const 3754 Assert (shape_function < fe_values.fe->dofs_per_cell,
3755 ExcIndexRange (shape_function, 0, fe_values.fe->dofs_per_cell));
3757 typename FVB::ExcAccessToUninitializedField(
"update_values"));
3764 = shape_function_data[shape_function].single_nonzero_component;
3770 return value_type();
3775 value_type return_value;
3776 const unsigned int comp =
3777 shape_function_data[shape_function].single_nonzero_component_index;
3778 return_value[value_type::unrolled_to_component_indices(comp)]
3779 = fe_values.finite_element_output.shape_values(snc,q_point);
3780 return return_value;
3784 value_type return_value;
3785 for (
unsigned int d = 0; d < value_type::n_independent_components; ++d)
3786 if (shape_function_data[shape_function].is_nonzero_shape_function_component[d])
3787 return_value[value_type::unrolled_to_component_indices(d)]
3788 = fe_values.finite_element_output.shape_values(shape_function_data[shape_function].row_index[d],q_point);
3789 return return_value;
3794 template <
int dim,
int spacedim>
3798 const unsigned int q_point)
const 3801 Assert (shape_function < fe_values.fe->dofs_per_cell,
3802 ExcIndexRange (shape_function, 0, fe_values.fe->dofs_per_cell));
3804 typename FVB::ExcAccessToUninitializedField(
"update_gradients"));
3806 const int snc = shape_function_data[shape_function].single_nonzero_component;
3812 return divergence_type();
3845 const unsigned int comp =
3846 shape_function_data[shape_function].single_nonzero_component_index;
3847 const unsigned int ii = value_type::unrolled_to_component_indices(comp)[0];
3848 const unsigned int jj = value_type::unrolled_to_component_indices(comp)[1];
3865 const ::Tensor<1, spacedim> phi_grad = fe_values.finite_element_output.shape_gradients[snc][q_point];
3867 divergence_type return_value;
3868 return_value[ii] = phi_grad[jj];
3871 return_value[jj] = phi_grad[ii];
3873 return return_value;
3878 Assert (
false, ExcNotImplemented());
3879 divergence_type return_value;
3880 return return_value;
3884 template <
int dim,
int spacedim>
3888 const unsigned int q_point)
const 3891 Assert (shape_function < fe_values.fe->dofs_per_cell,
3892 ExcIndexRange (shape_function, 0, fe_values.fe->dofs_per_cell));
3894 typename FVB::ExcAccessToUninitializedField(
"update_values"));
3901 = shape_function_data[shape_function].single_nonzero_component;
3907 return value_type();
3912 value_type return_value;
3913 const unsigned int comp =
3914 shape_function_data[shape_function].single_nonzero_component_index;
3916 return_value[indices] = fe_values.finite_element_output.shape_values(snc,q_point);
3917 return return_value;
3921 value_type return_value;
3922 for (
unsigned int d = 0; d < dim*dim; ++d)
3923 if (shape_function_data[shape_function].is_nonzero_shape_function_component[d])
3926 return_value[indices]
3927 = fe_values.finite_element_output.shape_values(shape_function_data[shape_function].row_index[d],q_point);
3929 return return_value;
3934 template <
int dim,
int spacedim>
3938 const unsigned int q_point)
const 3941 Assert (shape_function < fe_values.fe->dofs_per_cell,
3942 ExcIndexRange (shape_function, 0, fe_values.fe->dofs_per_cell));
3944 typename FVB::ExcAccessToUninitializedField(
"update_gradients"));
3946 const int snc = shape_function_data[shape_function].single_nonzero_component;
3952 return divergence_type();
3972 const unsigned int comp =
3973 shape_function_data[shape_function].single_nonzero_component_index;
3975 const unsigned int ii = indices[0];
3976 const unsigned int jj = indices[1];
3978 const ::Tensor<1, spacedim> phi_grad = fe_values.finite_element_output.shape_gradients[snc][q_point];
3980 divergence_type return_value;
3981 return_value[jj] = phi_grad[ii];
3983 return return_value;
3988 Assert (
false, ExcNotImplemented());
3989 divergence_type return_value;
3990 return return_value;
4001 template <
int dim,
int spacedim>
4009 0, fe_values_views_cache.scalars.size()));
4011 return fe_values_views_cache.scalars[scalar.
component];
4016 template <
int dim,
int spacedim>
4023 fe_values_views_cache.vectors.size(),
4025 0, fe_values_views_cache.vectors.size()));
4030 template <
int dim,
int spacedim>
4037 fe_values_views_cache.symmetric_second_order_tensors.size(),
4039 0, fe_values_views_cache.symmetric_second_order_tensors.size()));
4044 template <
int dim,
int spacedim>
4051 fe_values_views_cache.second_order_tensors.size(),
4053 0, fe_values_views_cache.second_order_tensors.size()));
4061 template <
int dim,
int spacedim>
4065 const unsigned int j)
const 4067 Assert (i < fe->dofs_per_cell,
4068 ExcIndexRange (i, 0, fe->dofs_per_cell));
4070 ExcAccessToUninitializedField(
"update_values"));
4071 Assert (fe->is_primitive (i),
4072 ExcShapeFunctionNotPrimitive(i));
4076 if (fe->is_primitive())
4077 return this->finite_element_output.shape_values(i,j);
4089 row = this->finite_element_output.shape_function_to_row_table[i * fe->n_components() + fe->system_to_component_index(i).first];
4090 return this->finite_element_output.shape_values(row, j);
4096 template <
int dim,
int spacedim>
4100 const unsigned int j,
4101 const unsigned int component)
const 4103 Assert (i < fe->dofs_per_cell,
4104 ExcIndexRange (i, 0, fe->dofs_per_cell));
4105 Assert (this->update_flags & update_values,
4106 ExcAccessToUninitializedField(
"update_values"));
4107 Assert (component < fe->n_components(),
4108 ExcIndexRange(component, 0, fe->n_components()));
4113 if (fe->get_nonzero_components(i)[component] ==
false)
4120 row = this->finite_element_output.shape_function_to_row_table[i * fe->n_components() + component];
4121 return this->finite_element_output.shape_values(row, j);
4126 template <
int dim,
int spacedim>
4130 const unsigned int j)
const 4132 Assert (i < fe->dofs_per_cell,
4133 ExcIndexRange (i, 0, fe->dofs_per_cell));
4135 ExcAccessToUninitializedField(
"update_gradients"));
4136 Assert (fe->is_primitive (i),
4137 ExcShapeFunctionNotPrimitive(i));
4141 if (fe->is_primitive())
4142 return this->finite_element_output.shape_gradients[i][j];
4154 row = this->finite_element_output.shape_function_to_row_table[i * fe->n_components() + fe->system_to_component_index(i).first];
4155 return this->finite_element_output.shape_gradients[row][j];
4161 template <
int dim,
int spacedim>
4165 const unsigned int j,
4166 const unsigned int component)
const 4168 Assert (i < fe->dofs_per_cell,
4169 ExcIndexRange (i, 0, fe->dofs_per_cell));
4170 Assert (this->update_flags & update_gradients,
4171 ExcAccessToUninitializedField(
"update_gradients"));
4172 Assert (component < fe->n_components(),
4173 ExcIndexRange(component, 0, fe->n_components()));
4178 if (fe->get_nonzero_components(i)[component] ==
false)
4185 row = this->finite_element_output.shape_function_to_row_table[i * fe->n_components() + component];
4186 return this->finite_element_output.shape_gradients[row][j];
4191 template <
int dim,
int spacedim>
4195 const unsigned int j)
const 4197 Assert (i < fe->dofs_per_cell,
4198 ExcIndexRange (i, 0, fe->dofs_per_cell));
4200 ExcAccessToUninitializedField(
"update_hessians"));
4201 Assert (fe->is_primitive (i),
4202 ExcShapeFunctionNotPrimitive(i));
4206 if (fe->is_primitive())
4207 return this->finite_element_output.shape_hessians[i][j];
4219 row = this->finite_element_output.shape_function_to_row_table[i * fe->n_components() + fe->system_to_component_index(i).first];
4220 return this->finite_element_output.shape_hessians[row][j];
4226 template <
int dim,
int spacedim>
4230 const unsigned int j,
4231 const unsigned int component)
const 4233 Assert (i < fe->dofs_per_cell,
4234 ExcIndexRange (i, 0, fe->dofs_per_cell));
4235 Assert (this->update_flags & update_hessians,
4236 ExcAccessToUninitializedField(
"update_hessians"));
4237 Assert (component < fe->n_components(),
4238 ExcIndexRange(component, 0, fe->n_components()));
4243 if (fe->get_nonzero_components(i)[component] ==
false)
4250 row = this->finite_element_output.shape_function_to_row_table[i * fe->n_components() + component];
4251 return this->finite_element_output.shape_hessians[row][j];
4256 template <
int dim,
int spacedim>
4260 const unsigned int j)
const 4262 Assert (i < fe->dofs_per_cell,
4263 ExcIndexRange (i, 0, fe->dofs_per_cell));
4264 Assert (this->update_flags & update_hessians,
4265 ExcAccessToUninitializedField(
"update_3rd_derivatives"));
4266 Assert (fe->is_primitive (i),
4267 ExcShapeFunctionNotPrimitive(i));
4271 if (fe->is_primitive())
4272 return this->finite_element_output.shape_3rd_derivatives[i][j];
4284 row = this->finite_element_output.shape_function_to_row_table[i * fe->n_components() + fe->system_to_component_index(i).first];
4285 return this->finite_element_output.shape_3rd_derivatives[row][j];
4291 template <
int dim,
int spacedim>
4295 const unsigned int j,
4296 const unsigned int component)
const 4298 Assert (i < fe->dofs_per_cell,
4299 ExcIndexRange (i, 0, fe->dofs_per_cell));
4300 Assert (this->update_flags & update_hessians,
4301 ExcAccessToUninitializedField(
"update_3rd_derivatives"));
4302 Assert (component < fe->n_components(),
4303 ExcIndexRange(component, 0, fe->n_components()));
4308 if (fe->get_nonzero_components(i)[component] ==
false)
4315 row = this->finite_element_output.shape_function_to_row_table[i * fe->n_components() + component];
4316 return this->finite_element_output.shape_3rd_derivatives[row][j];
4321 template <
int dim,
int spacedim>
4330 template <
int dim,
int spacedim>
4340 template <
int dim,
int spacedim>
4345 return this->update_flags;
4350 template <
int dim,
int spacedim>
4352 const std::vector<Point<spacedim> > &
4356 ExcAccessToUninitializedField(
"update_quadrature_points"));
4357 return this->mapping_output.quadrature_points;
4362 template <
int dim,
int spacedim>
4364 const std::vector<double> &
4368 ExcAccessToUninitializedField(
"update_JxW_values"));
4369 return this->mapping_output.JxW_values;
4374 template <
int dim,
int spacedim>
4376 const std::vector<DerivativeForm<1,dim,spacedim> > &
4380 ExcAccessToUninitializedField(
"update_jacobians"));
4381 return this->mapping_output.jacobians;
4386 template <
int dim,
int spacedim>
4388 const std::vector<DerivativeForm<2,dim,spacedim> > &
4392 ExcAccessToUninitializedField(
"update_jacobians_grads"));
4393 return this->mapping_output.jacobian_grads;
4398 template <
int dim,
int spacedim>
4404 ExcAccessToUninitializedField(
"update_jacobian_pushed_forward_grads"));
4405 return this->mapping_output.jacobian_pushed_forward_grads[i];
4410 template <
int dim,
int spacedim>
4412 const std::vector<Tensor<3,spacedim> > &
4416 ExcAccessToUninitializedField(
"update_jacobian_pushed_forward_grads"));
4417 return this->mapping_output.jacobian_pushed_forward_grads;
4422 template <
int dim,
int spacedim>
4428 ExcAccessToUninitializedField(
"update_jacobian_2nd_derivatives"));
4429 return this->mapping_output.jacobian_2nd_derivatives[i];
4434 template <
int dim,
int spacedim>
4436 const std::vector<DerivativeForm<3,dim,spacedim> > &
4440 ExcAccessToUninitializedField(
"update_jacobian_2nd_derivatives"));
4441 return this->mapping_output.jacobian_2nd_derivatives;
4446 template <
int dim,
int spacedim>
4452 ExcAccessToUninitializedField(
"update_jacobian_pushed_forward_2nd_derivatives"));
4453 return this->mapping_output.jacobian_pushed_forward_2nd_derivatives[i];
4458 template <
int dim,
int spacedim>
4460 const std::vector<Tensor<4,spacedim> > &
4464 ExcAccessToUninitializedField(
"update_jacobian_pushed_forward_2nd_derivatives"));
4465 return this->mapping_output.jacobian_pushed_forward_2nd_derivatives;
4470 template <
int dim,
int spacedim>
4476 ExcAccessToUninitializedField(
"update_jacobian_3rd_derivatives"));
4477 return this->mapping_output.jacobian_3rd_derivatives[i];
4482 template <
int dim,
int spacedim>
4484 const std::vector<DerivativeForm<4,dim,spacedim> > &
4488 ExcAccessToUninitializedField(
"update_jacobian_3rd_derivatives"));
4489 return this->mapping_output.jacobian_3rd_derivatives;
4494 template <
int dim,
int spacedim>
4500 ExcAccessToUninitializedField(
"update_jacobian_pushed_forward_3rd_derivatives"));
4501 return this->mapping_output.jacobian_pushed_forward_3rd_derivatives[i];
4506 template <
int dim,
int spacedim>
4508 const std::vector<Tensor<5,spacedim> > &
4512 ExcAccessToUninitializedField(
"update_jacobian_pushed_forward_3rd_derivatives"));
4513 return this->mapping_output.jacobian_pushed_forward_3rd_derivatives;
4517 template <
int dim,
int spacedim>
4519 const std::vector<DerivativeForm<1,spacedim,dim> > &
4523 ExcAccessToUninitializedField(
"update_inverse_jacobians"));
4524 return this->mapping_output.inverse_jacobians;
4529 template <
int dim,
int spacedim>
4535 ExcAccessToUninitializedField(
"update_quadrature_points"));
4536 Assert (i<this->mapping_output.quadrature_points.size(),
4537 ExcIndexRange(i, 0, this->mapping_output.quadrature_points.size()));
4539 return this->mapping_output.quadrature_points[i];
4545 template <
int dim,
int spacedim>
4551 ExcAccessToUninitializedField(
"update_JxW_values"));
4552 Assert (i<this->mapping_output.JxW_values.size(),
4553 ExcIndexRange(i, 0, this->mapping_output.JxW_values.size()));
4555 return this->mapping_output.JxW_values[i];
4560 template <
int dim,
int spacedim>
4566 ExcAccessToUninitializedField(
"update_jacobians"));
4567 Assert (i<this->mapping_output.jacobians.size(),
4568 ExcIndexRange(i, 0, this->mapping_output.jacobians.size()));
4570 return this->mapping_output.jacobians[i];
4575 template <
int dim,
int spacedim>
4581 ExcAccessToUninitializedField(
"update_jacobians_grads"));
4582 Assert (i<this->mapping_output.jacobian_grads.size(),
4583 ExcIndexRange(i, 0, this->mapping_output.jacobian_grads.size()));
4585 return this->mapping_output.jacobian_grads[i];
4590 template <
int dim,
int spacedim>
4596 ExcAccessToUninitializedField(
"update_inverse_jacobians"));
4597 Assert (i<this->mapping_output.inverse_jacobians.size(),
4598 ExcIndexRange(i, 0, this->mapping_output.inverse_jacobians.size()));
4600 return this->mapping_output.inverse_jacobians[i];
4604 template <
int dim,
int spacedim>
4611 typename FVB::ExcAccessToUninitializedField(
"update_normal_vectors"));
4612 Assert (i<this->mapping_output.normal_vectors.size(),
4613 ExcIndexRange(i, 0, this->mapping_output.normal_vectors.size()));
4615 return this->mapping_output.normal_vectors[i];
4623 template <
int dim,
int spacedim>
4633 template <
int dim,
int spacedim>
4645 template <
int dim,
int spacedim>
4650 return present_face_index;
4656 template <
int dim,
int spacedim>
4666 template <
int dim,
int spacedim>
4676 template <
int dim,
int spacedim>
4686 template <
int dim,
int spacedim>
4692 Assert (i<this->mapping_output.boundary_forms.size(),
4693 ExcIndexRange(i, 0, this->mapping_output.boundary_forms.size()));
4695 typename FVB::ExcAccessToUninitializedField(
"update_boundary_forms"));
4697 return this->mapping_output.boundary_forms[i];
4702 DEAL_II_NAMESPACE_CLOSE
Transformed quadrature weights.
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
::Tensor< 1, spacedim > divergence_type
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
std::vector<::FEValuesViews::Scalar< dim, spacedim > > scalars
const unsigned int dofs_per_cell
const unsigned int component
const Tensor< 3, spacedim > & shape_3rd_derivative(const unsigned int function_no, const unsigned int point_no) const
::Tensor< 3, spacedim > hessian_type
const Quadrature< dim-1 > quadrature
::ExceptionBase & ExcMessage(std::string arg1)
int single_nonzero_component
::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
const unsigned int first_tensor_component
UpdateFlags get_update_flags() const
Transformed quadrature points.
::Tensor< 1, spacedim > gradient_type
::SymmetricTensor< 2, spacedim > value_type
const DerivativeForm< 3, dim, spacedim > & jacobian_2nd_derivative(const unsigned int quadrature_point) const
const FEValuesViews::Scalar< dim, spacedim > & operator[](const FEValuesExtractors::Scalar &scalar) const
const DerivativeForm< 2, dim, spacedim > & jacobian_grad(const unsigned int quadrature_point) const
const std::vector< DerivativeForm< 1, dim, spacedim > > & get_jacobians() const
const Mapping< dim, spacedim > & get_mapping() const
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
const std::vector< Tensor< 5, spacedim > > & get_jacobian_pushed_forward_3rd_derivatives() const
const std::vector< Point< spacedim > > & get_quadrature_points() const
BlockDynamicSparsityPattern BlockCompressedSparsityPattern DEAL_II_DEPRECATED
const Point< spacedim > & quadrature_point(const unsigned int q) const
const std::vector< DerivativeForm< 2, dim, spacedim > > & get_jacobian_grads() const
const unsigned int first_tensor_component
::Tensor< 1, spacedim > value_type
#define DeclException1(Exception1, type1, outsequence)
const double & shape_value(const unsigned int function_no, const unsigned int point_no) const
Third derivatives of shape functions.
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)
const std::vector< DerivativeForm< 1, spacedim, dim > > & get_inverse_jacobians() const
::Tensor< 1, spacedim > divergence_type
const FEFaceValues< dim, spacedim > & get_present_fe_values() const
Second derivatives of shape functions.
Gradient of volume element.
const Tensor< 2, spacedim > & shape_hessian(const unsigned int function_no, const unsigned int point_no) const
const Tensor< 1, spacedim > & boundary_form(const unsigned int i) const
static TableIndices< rank_ > unrolled_to_component_indices(const unsigned int i)
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 std::vector< Tensor< 4, spacedim > > & get_jacobian_pushed_forward_2nd_derivatives() const
::Tensor< 2, spacedim > hessian_type
const FEValuesBase< dim, spacedim > & fe_values
int single_nonzero_component
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.
::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
const FESubfaceValues< dim, spacedim > & get_present_fe_values() const
const FEValuesBase< dim, spacedim > & fe_values
bool is_nonzero_shape_function_component
const Tensor< 4, spacedim > & jacobian_pushed_forward_2nd_derivative(const unsigned int quadrature_point) const
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
const Tensor< 1, spacedim > & shape_grad(const unsigned int function_no, const unsigned int quadrature_point) const
std::vector< ShapeFunctionData > shape_function_data
const FEValuesBase< dim, spacedim > & fe_values
::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