15#ifndef dealii_fe_values_views_h
16#define dealii_fe_values_views_h
38template <
int dim,
int spacedim = dim>
48 template <
int dim,
typename NumberType =
double>
57 template <
typename NumberType>
69 template <
typename NumberType>
81 template <
typename NumberType>
124 template <
int dim,
int spacedim = dim>
162 template <
typename Number>
171 template <
typename Number>
181 template <
typename Number>
191 template <
typename Number>
201 template <
typename Number>
289 value(const
unsigned int shape_function, const
unsigned int q_point) const;
303 const
unsigned int q_point) const;
317 const
unsigned int q_point) const;
331 const
unsigned int q_point) const;
350 template <typename Number>
389 template <class InputVector>
392 const InputVector &dof_values,
413 template <typename Number>
425 template <class InputVector>
428 const InputVector &dof_values,
449 template <typename Number>
461 template <class InputVector>
464 const InputVector &dof_values,
487 template <typename Number>
499 template <class InputVector>
502 const InputVector &dof_values,
525 template <typename Number>
538 template <class InputVector>
541 const InputVector &dof_values,
544 &third_derivatives) const;
596 template <
int dim,
int spacedim = dim>
642 using curl_type = typename ::internal::CurlType<spacedim>::type;
664 template <
typename Number>
673 template <
typename Number>
683 template <
typename Number>
693 template <
typename Number>
703 template <
typename Number>
713 template <
typename Number>
722 template <
typename Number>
732 template <
typename Number>
789 const unsigned int first_vector_component);
839 value(
const unsigned int shape_function,
const unsigned int q_point)
const;
856 const unsigned int q_point)
const;
875 const unsigned int q_point)
const;
889 const unsigned int q_point)
const;
912 curl(
const unsigned int shape_function,
const unsigned int q_point)
const;
926 const unsigned int q_point)
const;
940 const unsigned int q_point)
const;
959 template <
typename Number>
998 template <
class InputVector>
1001 const InputVector &dof_values,
1022 template <
typename Number>
1034 template <
class InputVector>
1037 const InputVector &dof_values,
1064 template <
typename Number>
1066 get_function_symmetric_gradients(
1069 &symmetric_gradients)
const;
1077 template <
class InputVector>
1079 get_function_symmetric_gradients_from_local_dof_values(
1080 const InputVector &dof_values,
1083 &symmetric_gradients)
const;
1103 template <
typename Number>
1105 get_function_divergences(
1115 template <
class InputVector>
1117 get_function_divergences_from_local_dof_values(
1118 const InputVector &dof_values,
1120 &divergences)
const;
1140 template <
typename Number>
1151 template <
class InputVector>
1153 get_function_curls_from_local_dof_values(
1154 const InputVector &dof_values,
1175 template <
typename Number>
1187 template <
class InputVector>
1190 const InputVector &dof_values,
1212 template <
typename Number>
1224 template <
class InputVector>
1227 const InputVector &dof_values,
1249 template <
typename Number>
1262 template <
class InputVector>
1265 const InputVector &dof_values,
1268 &third_derivatives)
const;
1289 template <
int rank,
int dim,
int spacedim = dim>
1314 template <
int dim,
int spacedim>
1343 template <
typename Number>
1352 template <
typename Number>
1371 bool is_nonzero_shape_function_component
1372 [value_type::n_independent_components];
1383 unsigned int row_index[value_type::n_independent_components];
1416 const unsigned int first_tensor_component);
1461 value(const
unsigned int shape_function, const
unsigned int q_point) const;
1477 divergence(const
unsigned int shape_function,
1478 const
unsigned int q_point) const;
1497 template <typename Number>
1499 get_function_values(const
ReadVector<Number> &fe_function,
1536 template <class InputVector>
1538 get_function_values_from_local_dof_values(
1539 const InputVector &dof_values,
1564 template <typename Number>
1566 get_function_divergences(
1576 template <class InputVector>
1578 get_function_divergences_from_local_dof_values(
1579 const InputVector &dof_values,
1581 &divergences) const;
1593 unsigned int first_tensor_component;
1602 template <
int rank,
int dim,
int spacedim = dim>
1623 template <
int dim,
int spacedim>
1650 template <
typename Number>
1659 template <
typename Number>
1669 template <
typename Number>
1688 bool is_nonzero_shape_function_component
1689 [value_type::n_independent_components];
1700 unsigned int row_index[value_type::n_independent_components];
1750 const unsigned int first_tensor_component);
1783 value(
const unsigned int shape_function,
const unsigned int q_point)
const;
1800 const unsigned int q_point)
const;
1817 const unsigned int q_point)
const;
1836 template <
typename Number>
1875 template <
class InputVector>
1877 get_function_values_from_local_dof_values(
1878 const InputVector &dof_values,
1903 template <
typename Number>
1905 get_function_divergences(
1915 template <
class InputVector>
1917 get_function_divergences_from_local_dof_values(
1918 const InputVector &dof_values,
1920 &divergences)
const;
1938 template <
typename Number>
1940 get_function_gradients(
1950 template <
class InputVector>
1952 get_function_gradients_from_local_dof_values(
1953 const InputVector &dof_values,
1986 template <
int dim,
int spacedim,
typename Extractor>
1997 template <
int dim,
int spacedim>
2000 using type = typename ::FEValuesViews::Scalar<dim, spacedim>;
2010 template <
int dim,
int spacedim>
2013 using type = typename ::FEValuesViews::Vector<dim, spacedim>;
2023 template <
int dim,
int spacedim,
int rank>
2026 using type = typename ::FEValuesViews::Tensor<rank, dim, spacedim>;
2036 template <
int dim,
int spacedim,
int rank>
2040 typename ::FEValuesViews::SymmetricTensor<rank, dim, spacedim>;
2050 template <
int dim,
int spacedim>
2057 std::vector<Lazy<::FEValuesViews::Scalar<dim, spacedim>>>
scalars;
2058 std::vector<Lazy<::FEValuesViews::Vector<dim, spacedim>>>
vectors;
2062 std::vector<Lazy<::FEValuesViews::Tensor<2, dim, spacedim>>>
2079 template <
int dim,
int spacedim,
typename Extractor>
2080 using View = typename ::internal::FEValuesViews::
2081 ViewType<dim, spacedim, Extractor>::type;
2090 template <
int dim,
int spacedim>
2091 inline typename Scalar<dim, spacedim>::value_type
2092 Scalar<dim, spacedim>::value(
const unsigned int shape_function,
2093 const unsigned int q_point)
const
2099 "update_values"))));
2104 if (shape_function_data[shape_function].is_nonzero_shape_function_component)
2105 return fe_values->finite_element_output.shape_values(
2106 shape_function_data[shape_function].row_index, q_point);
2113 template <
int dim,
int spacedim>
2114 inline typename Scalar<dim, spacedim>::gradient_type
2115 Scalar<dim, spacedim>::gradient(
const unsigned int shape_function,
2116 const unsigned int q_point)
const
2121 "update_gradients")));
2126 if (shape_function_data[shape_function].is_nonzero_shape_function_component)
2127 return fe_values->finite_element_output
2128 .shape_gradients[shape_function_data[shape_function].row_index]
2131 return gradient_type();
2136 template <
int dim,
int spacedim>
2137 inline typename Scalar<dim, spacedim>::hessian_type
2138 Scalar<dim, spacedim>::hessian(
const unsigned int shape_function,
2139 const unsigned int q_point)
const
2144 "update_hessians")));
2149 if (shape_function_data[shape_function].is_nonzero_shape_function_component)
2150 return fe_values->finite_element_output
2151 .shape_hessians[shape_function_data[shape_function].row_index][q_point];
2153 return hessian_type();
2158 template <
int dim,
int spacedim>
2159 inline typename Scalar<dim, spacedim>::third_derivative_type
2160 Scalar<dim, spacedim>::third_derivative(
const unsigned int shape_function,
2161 const unsigned int q_point)
const
2166 "update_3rd_derivatives")));
2171 if (shape_function_data[shape_function].is_nonzero_shape_function_component)
2172 return fe_values->finite_element_output
2173 .shape_3rd_derivatives[shape_function_data[shape_function].row_index]
2176 return third_derivative_type();
2181 template <
int dim,
int spacedim>
2184 const unsigned int q_point)
const
2193 shape_function_data[shape_function].single_nonzero_component;
2195 return value_type();
2198 value_type return_value;
2199 return_value[shape_function_data[shape_function]
2200 .single_nonzero_component_index] =
2201 fe_values->finite_element_output.shape_values(snc, q_point);
2202 return return_value;
2206 value_type return_value;
2207 for (
unsigned int d = 0;
d < dim; ++
d)
2208 if (shape_function_data[shape_function]
2209 .is_nonzero_shape_function_component[d])
2210 return_value[
d] = fe_values->finite_element_output.shape_values(
2211 shape_function_data[shape_function].row_index[d], q_point);
2213 return return_value;
2219 template <
int dim,
int spacedim>
2222 const unsigned int q_point)
const
2227 "update_gradients")));
2231 shape_function_data[shape_function].single_nonzero_component;
2233 return gradient_type();
2236 gradient_type return_value;
2237 return_value[shape_function_data[shape_function]
2238 .single_nonzero_component_index] =
2239 fe_values->finite_element_output.shape_gradients[snc][q_point];
2240 return return_value;
2244 gradient_type return_value;
2245 for (
unsigned int d = 0;
d < dim; ++
d)
2246 if (shape_function_data[shape_function]
2247 .is_nonzero_shape_function_component[d])
2249 fe_values->finite_element_output.shape_gradients
2250 [shape_function_data[shape_function].row_index[
d]][q_point];
2252 return return_value;
2258 template <
int dim,
int spacedim>
2261 const unsigned int q_point)
const
2267 "update_gradients")));
2271 shape_function_data[shape_function].single_nonzero_component;
2273 return divergence_type();
2275 return fe_values->finite_element_output
2276 .shape_gradients[snc][q_point][shape_function_data[shape_function]
2277 .single_nonzero_component_index];
2280 divergence_type return_value = 0;
2281 for (
unsigned int d = 0;
d < dim; ++
d)
2282 if (shape_function_data[shape_function]
2283 .is_nonzero_shape_function_component[d])
2285 fe_values->finite_element_output.shape_gradients
2286 [shape_function_data[shape_function].row_index[
d]][q_point][
d];
2288 return return_value;
2294 template <
int dim,
int spacedim>
2297 const unsigned int q_point)
const
2304 "update_gradients")));
2307 shape_function_data[shape_function].single_nonzero_component;
2319 "Computing the curl in 1d is not a useful operation"));
2327 curl_type return_value;
2330 if (shape_function_data[shape_function]
2331 .single_nonzero_component_index == 0)
2333 -1.0 * fe_values->finite_element_output
2334 .shape_gradients[snc][q_point][1];
2336 return_value[0] = fe_values->finite_element_output
2337 .shape_gradients[snc][q_point][0];
2339 return return_value;
2344 curl_type return_value;
2346 return_value[0] = 0.0;
2348 if (shape_function_data[shape_function]
2349 .is_nonzero_shape_function_component[0])
2351 fe_values->finite_element_output
2352 .shape_gradients[shape_function_data[shape_function]
2353 .row_index[0]][q_point][1];
2355 if (shape_function_data[shape_function]
2356 .is_nonzero_shape_function_component[1])
2358 fe_values->finite_element_output
2359 .shape_gradients[shape_function_data[shape_function]
2360 .row_index[1]][q_point][0];
2362 return return_value;
2370 curl_type return_value;
2372 switch (shape_function_data[shape_function]
2373 .single_nonzero_component_index)
2377 return_value[0] = 0;
2378 return_value[1] = fe_values->finite_element_output
2379 .shape_gradients[snc][q_point][2];
2381 -1.0 * fe_values->finite_element_output
2382 .shape_gradients[snc][q_point][1];
2383 return return_value;
2389 -1.0 * fe_values->finite_element_output
2390 .shape_gradients[snc][q_point][2];
2391 return_value[1] = 0;
2392 return_value[2] = fe_values->finite_element_output
2393 .shape_gradients[snc][q_point][0];
2394 return return_value;
2399 return_value[0] = fe_values->finite_element_output
2400 .shape_gradients[snc][q_point][1];
2402 -1.0 * fe_values->finite_element_output
2403 .shape_gradients[snc][q_point][0];
2404 return_value[2] = 0;
2405 return return_value;
2412 curl_type return_value;
2414 for (
unsigned int i = 0; i < dim; ++i)
2415 return_value[i] = 0.0;
2417 if (shape_function_data[shape_function]
2418 .is_nonzero_shape_function_component[0])
2421 fe_values->finite_element_output
2422 .shape_gradients[shape_function_data[shape_function]
2423 .row_index[0]][q_point][2];
2425 fe_values->finite_element_output
2426 .shape_gradients[shape_function_data[shape_function]
2427 .row_index[0]][q_point][1];
2430 if (shape_function_data[shape_function]
2431 .is_nonzero_shape_function_component[1])
2434 fe_values->finite_element_output
2435 .shape_gradients[shape_function_data[shape_function]
2436 .row_index[1]][q_point][2];
2438 fe_values->finite_element_output
2439 .shape_gradients[shape_function_data[shape_function]
2440 .row_index[1]][q_point][0];
2443 if (shape_function_data[shape_function]
2444 .is_nonzero_shape_function_component[2])
2447 fe_values->finite_element_output
2448 .shape_gradients[shape_function_data[shape_function]
2449 .row_index[2]][q_point][1];
2451 fe_values->finite_element_output
2452 .shape_gradients[shape_function_data[shape_function]
2453 .row_index[2]][q_point][0];
2456 return return_value;
2467 template <
int dim,
int spacedim>
2470 const unsigned int q_point)
const
2476 "update_hessians")));
2480 shape_function_data[shape_function].single_nonzero_component;
2482 return hessian_type();
2485 hessian_type return_value;
2486 return_value[shape_function_data[shape_function]
2487 .single_nonzero_component_index] =
2488 fe_values->finite_element_output.shape_hessians[snc][q_point];
2489 return return_value;
2493 hessian_type return_value;
2494 for (
unsigned int d = 0;
d < dim; ++
d)
2495 if (shape_function_data[shape_function]
2496 .is_nonzero_shape_function_component[d])
2498 fe_values->finite_element_output.shape_hessians
2499 [shape_function_data[shape_function].row_index[
d]][q_point];
2501 return return_value;
2507 template <
int dim,
int spacedim>
2510 const unsigned int q_point)
const
2516 "update_3rd_derivatives")));
2520 shape_function_data[shape_function].single_nonzero_component;
2522 return third_derivative_type();
2525 third_derivative_type return_value;
2526 return_value[shape_function_data[shape_function]
2527 .single_nonzero_component_index] =
2528 fe_values->finite_element_output.shape_3rd_derivatives[snc][q_point];
2529 return return_value;
2533 third_derivative_type return_value;
2534 for (
unsigned int d = 0;
d < dim; ++
d)
2535 if (shape_function_data[shape_function]
2536 .is_nonzero_shape_function_component[d])
2538 fe_values->finite_element_output.shape_3rd_derivatives
2539 [shape_function_data[shape_function].row_index[
d]][q_point];
2541 return return_value;
2553 inline ::SymmetricTensor<2, 1>
2554 symmetrize_single_row(
const unsigned int n, const ::Tensor<1, 1> &t)
2564 inline ::SymmetricTensor<2, 2>
2565 symmetrize_single_row(
const unsigned int n, const ::Tensor<1, 2> &t)
2571 return {{t[0], 0, t[1] / 2}};
2575 return {{0, t[1], t[0] / 2}};
2587 inline ::SymmetricTensor<2, 3>
2588 symmetrize_single_row(
const unsigned int n, const ::Tensor<1, 3> &t)
2594 return {{t[0], 0, 0, t[1] / 2, t[2] / 2, 0}};
2598 return {{0, t[1], 0, t[0] / 2, 0, t[2] / 2}};
2602 return {{0, 0, t[2], 0, t[0] / 2, t[1] / 2}};
2615 template <
int dim,
int spacedim>
2618 const unsigned int q_point)
const
2623 "update_gradients")));
2627 shape_function_data[shape_function].single_nonzero_component;
2629 return symmetric_gradient_type();
2631 return internal::symmetrize_single_row(
2632 shape_function_data[shape_function].single_nonzero_component_index,
2633 fe_values->finite_element_output.shape_gradients[snc][q_point]);
2636 gradient_type return_value;
2637 for (
unsigned int d = 0;
d < dim; ++
d)
2638 if (shape_function_data[shape_function]
2639 .is_nonzero_shape_function_component[d])
2641 fe_values->finite_element_output.shape_gradients
2642 [shape_function_data[shape_function].row_index[
d]][q_point];
2650 template <
int dim,
int spacedim>
2653 const unsigned int q_point)
const
2663 shape_function_data[shape_function].single_nonzero_component;
2668 return value_type();
2672 value_type return_value;
2673 const unsigned int comp =
2674 shape_function_data[shape_function].single_nonzero_component_index;
2675 return_value[value_type::unrolled_to_component_indices(comp)] =
2676 fe_values->finite_element_output.shape_values(snc, q_point);
2677 return return_value;
2681 value_type return_value;
2682 for (
unsigned int d = 0;
d < value_type::n_independent_components; ++
d)
2683 if (shape_function_data[shape_function]
2684 .is_nonzero_shape_function_component[d])
2685 return_value[value_type::unrolled_to_component_indices(d)] =
2686 fe_values->finite_element_output.shape_values(
2687 shape_function_data[shape_function].row_index[d], q_point);
2688 return return_value;
2694 template <
int dim,
int spacedim>
2697 const unsigned int shape_function,
2698 const unsigned int q_point)
const
2703 "update_gradients")));
2706 shape_function_data[shape_function].single_nonzero_component;
2711 return divergence_type();
2734 const unsigned int comp =
2735 shape_function_data[shape_function].single_nonzero_component_index;
2736 const unsigned int ii =
2737 value_type::unrolled_to_component_indices(comp)[0];
2738 const unsigned int jj =
2739 value_type::unrolled_to_component_indices(comp)[1];
2752 const ::Tensor<1, spacedim> &phi_grad =
2753 fe_values->finite_element_output.shape_gradients[snc][q_point];
2755 divergence_type return_value;
2756 return_value[ii] = phi_grad[jj];
2759 return_value[jj] = phi_grad[ii];
2761 return return_value;
2766 divergence_type return_value;
2767 return return_value;
2773 template <
int dim,
int spacedim>
2776 const unsigned int q_point)
const
2786 shape_function_data[shape_function].single_nonzero_component;
2791 return value_type();
2795 value_type return_value;
2796 const unsigned int comp =
2797 shape_function_data[shape_function].single_nonzero_component_index;
2800 return_value[indices] =
2801 fe_values->finite_element_output.shape_values(snc, q_point);
2802 return return_value;
2806 value_type return_value;
2807 for (
unsigned int d = 0;
d < dim * dim; ++
d)
2808 if (shape_function_data[shape_function]
2809 .is_nonzero_shape_function_component[d])
2813 return_value[indices] =
2814 fe_values->finite_element_output.shape_values(
2815 shape_function_data[shape_function].row_index[d], q_point);
2817 return return_value;
2823 template <
int dim,
int spacedim>
2826 const unsigned int q_point)
const
2831 "update_gradients")));
2834 shape_function_data[shape_function].single_nonzero_component;
2839 return divergence_type();
2853 const unsigned int comp =
2854 shape_function_data[shape_function].single_nonzero_component_index;
2857 const unsigned int ii = indices[0];
2858 const unsigned int jj = indices[1];
2860 const ::Tensor<1, spacedim> &phi_grad =
2861 fe_values->finite_element_output.shape_gradients[snc][q_point];
2863 divergence_type return_value;
2865 return_value[ii] = phi_grad[jj];
2867 return return_value;
2872 divergence_type return_value;
2873 return return_value;
2879 template <
int dim,
int spacedim>
2882 const unsigned int q_point)
const
2887 "update_gradients")));
2890 shape_function_data[shape_function].single_nonzero_component;
2895 return gradient_type();
2909 const unsigned int comp =
2910 shape_function_data[shape_function].single_nonzero_component_index;
2913 const unsigned int ii = indices[0];
2914 const unsigned int jj = indices[1];
2916 const ::Tensor<1, spacedim> &phi_grad =
2917 fe_values->finite_element_output.shape_gradients[snc][q_point];
2919 gradient_type return_value;
2920 return_value[ii][jj] = phi_grad;
2922 return return_value;
2927 gradient_type return_value;
2928 return return_value;
void get_function_third_derivatives_from_local_dof_values(const InputVector &dof_values, std::vector< solution_third_derivative_type< typename InputVector::value_type > > &third_derivatives) const
Scalar & operator=(const Scalar< dim, spacedim > &)=delete
void get_function_laplacians(const ReadVector< Number > &fe_function, std::vector< solution_laplacian_type< Number > > &laplacians) const
value_type value(const unsigned int shape_function, const unsigned int q_point) const
typename ProductType< Number, value_type >::type solution_laplacian_type
void get_function_hessians_from_local_dof_values(const InputVector &dof_values, std::vector< solution_hessian_type< typename InputVector::value_type > > &hessians) const
void get_function_values_from_local_dof_values(const InputVector &dof_values, std::vector< solution_value_type< typename InputVector::value_type > > &values) const
gradient_type gradient(const unsigned int shape_function, const unsigned int q_point) const
std::vector< ShapeFunctionData > shape_function_data
void get_function_gradients(const ReadVector< Number > &fe_function, std::vector< solution_gradient_type< Number > > &gradients) const
Scalar & operator=(Scalar< dim, spacedim > &&) noexcept=default
third_derivative_type third_derivative(const unsigned int shape_function, const unsigned int q_point) const
typename ProductType< Number, gradient_type >::type solution_gradient_type
typename ProductType< Number, value_type >::type solution_value_type
void get_function_values(const ReadVector< Number > &fe_function, std::vector< solution_value_type< Number > > &values) const
void get_function_laplacians_from_local_dof_values(const InputVector &dof_values, std::vector< solution_laplacian_type< typename InputVector::value_type > > &laplacians) const
SmartPointer< const FEValuesBase< dim, spacedim > > fe_values
void get_function_gradients_from_local_dof_values(const InputVector &dof_values, std::vector< solution_gradient_type< typename InputVector::value_type > > &gradients) const
Scalar(Scalar< dim, spacedim > &&)=default
typename ProductType< Number, hessian_type >::type solution_hessian_type
Scalar(const Scalar< dim, spacedim > &)=delete
hessian_type hessian(const unsigned int shape_function, const unsigned int q_point) const
void get_function_third_derivatives(const ReadVector< Number > &fe_function, std::vector< solution_third_derivative_type< Number > > &third_derivatives) const
typename ProductType< Number, third_derivative_type >::type solution_third_derivative_type
void get_function_hessians(const ReadVector< Number > &fe_function, std::vector< solution_hessian_type< Number > > &hessians) const
SymmetricTensor(const SymmetricTensor< 2, dim, spacedim > &)=delete
SymmetricTensor(SymmetricTensor< 2, dim, spacedim > &&)=default
SymmetricTensor & operator=(SymmetricTensor< 2, dim, spacedim > &&) noexcept=default
typename ProductType< Number, value_type >::type solution_value_type
SymmetricTensor & operator=(const SymmetricTensor< 2, dim, spacedim > &)=delete
typename ProductType< Number, divergence_type >::type solution_divergence_type
typename ProductType< Number, value_type >::type solution_value_type
divergence_type divergence(const unsigned int shape_function, const unsigned int q_point) const
unsigned int first_tensor_component
typename ProductType< Number, divergence_type >::type solution_divergence_type
typename ProductType< Number, gradient_type >::type solution_gradient_type
Tensor(const Tensor< 2, dim, spacedim > &)=delete
SmartPointer< const FEValuesBase< dim, spacedim > > fe_values
Tensor & operator=(const Tensor< 2, dim, spacedim > &)=delete
Tensor & operator=(Tensor< 2, dim, spacedim > &&)=default
value_type value(const unsigned int shape_function, const unsigned int q_point) const
gradient_type gradient(const unsigned int shape_function, const unsigned int q_point) const
Tensor(Tensor< 2, dim, spacedim > &&)=default
std::vector< ShapeFunctionData > shape_function_data
Vector(const Vector< dim, spacedim > &)=delete
unsigned int first_vector_component
hessian_type hessian(const unsigned int shape_function, const unsigned int q_point) const
typename ProductType< Number, value_type >::type solution_laplacian_type
typename ProductType< Number, third_derivative_type >::type solution_third_derivative_type
Vector & operator=(Vector< dim, spacedim > &&)=default
typename ProductType< Number, hessian_type >::type solution_hessian_type
typename ProductType< Number, divergence_type >::type solution_divergence_type
gradient_type gradient(const unsigned int shape_function, const unsigned int q_point) const
divergence_type divergence(const unsigned int shape_function, const unsigned int q_point) const
Vector & operator=(const Vector< dim, spacedim > &)=delete
third_derivative_type third_derivative(const unsigned int shape_function, const unsigned int q_point) const
typename ProductType< Number, value_type >::type solution_value_type
symmetric_gradient_type symmetric_gradient(const unsigned int shape_function, const unsigned int q_point) const
Vector(Vector< dim, spacedim > &&)=default
typename ProductType< Number, gradient_type >::type solution_gradient_type
typename ::internal::CurlType< spacedim >::type curl_type
typename ProductType< Number, curl_type >::type solution_curl_type
typename ProductType< Number, symmetric_gradient_type >::type solution_symmetric_gradient_type
std::vector< ShapeFunctionData > shape_function_data
value_type value(const unsigned int shape_function, const unsigned int q_point) const
curl_type curl(const unsigned int shape_function, const unsigned int q_point) const
SmartPointer< const FEValuesBase< dim, spacedim > > fe_values
std::conditional_t< rank_==1, Number, Tensor< rank_ - 1, dim, Number > > value_type
static constexpr TableIndices< rank_ > unrolled_to_component_indices(const unsigned int i)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcMessage(std::string arg1)
@ update_hessians
Second derivatives of shape functions.
@ update_values
Shape function values.
@ update_3rd_derivatives
Third derivatives of shape functions.
@ update_gradients
Shape function gradients.
#define DEAL_II_ASSERT_UNREACHABLE()
#define DEAL_II_NOT_IMPLEMENTED()
typename ::internal::FEValuesViews:: ViewType< dim, spacedim, Extractor >::type View
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
bool is_nonzero_shape_function_component
int single_nonzero_component
unsigned int single_nonzero_component_index
unsigned int single_nonzero_component_index
int single_nonzero_component
int single_nonzero_component
unsigned int single_nonzero_component_index
typename internal::ProductTypeImpl< std::decay_t< T >, std::decay_t< U > >::type type
std::vector< Lazy<::FEValuesViews::Vector< dim, spacedim > > > vectors
std::vector< Lazy<::FEValuesViews::SymmetricTensor< 2, dim, spacedim > > > symmetric_second_order_tensors
std::vector< Lazy<::FEValuesViews::Tensor< 2, dim, spacedim > > > second_order_tensors
std::vector< Lazy<::FEValuesViews::Scalar< dim, spacedim > > > scalars
DEAL_II_HOST constexpr SymmetricTensor< 2, dim, Number > symmetrize(const Tensor< 2, dim, Number > &t)