15#ifndef dealii_fe_point_evaluation_h
16#define dealii_fe_point_evaluation_h
47 <<
"You are requesting information from an FEPointEvaluationBase "
48 <<
"object for which this kind of information has not been computed. "
49 <<
"What information these objects compute is determined by the update_* "
50 <<
"flags you pass to MappingInfo() in the Constructor. Here, "
51 <<
"the operation you are attempting requires the <" << arg1
52 <<
"> flag to be set, but it was apparently not specified "
53 <<
"upon initialization.");
63 typename Enable =
void>
69 typename ::internal::VectorizedArrayTrait<
78 n_components == spacedim,
90 const unsigned int component,
94 result[component] = vector_entry;
108 for (
unsigned int c = 0; c < n_components; ++c)
109 result_scalar[c] = result[c].sum();
111 return result_scalar;
119 return result[component].
sum();
124 const unsigned int vector_lane,
127 for (
unsigned int i = 0; i < n_components; ++i)
128 for (
unsigned int d = 0; d < dim; ++d)
131 value[d][i], vector_lane);
136 const unsigned int vector_lane,
139 for (
unsigned int i = 0; i < n_components; ++i)
140 for (
unsigned int d = 0; d < dim; ++d)
142 value[d][i], vector_lane) = result[i][d];
147 const unsigned int vector_lane,
150 for (
unsigned int i = 0; i < n_components; ++i)
151 for (
unsigned int d = 0; d < dim; ++d)
153 value[d][i], vector_lane) = result[i][d];
158 const unsigned int vector_lane)
160 for (
unsigned int i = 0; i < n_components; ++i)
161 for (
unsigned int d = 0; d < spacedim; ++d)
168 const unsigned int vector_lane,
171 for (
unsigned int i = 0; i < n_components; ++i)
172 result[i] =
value[i][vector_lane];
185 const unsigned int vector_lane,
188 for (
unsigned int i = 0; i < n_components; ++i)
189 value[i][vector_lane] = result[i];
203 for (
unsigned int i = 0; i < n_components; ++i)
210 const unsigned int vector_lane,
211 const unsigned int component,
215 vector_lane) += shape_value;
220 const unsigned int vector_lane,
221 const unsigned int component)
229 const unsigned int vector_lane,
230 const unsigned int component,
233 for (
unsigned int d = 0; d < spacedim; ++d)
241 const unsigned int vector_lane,
242 const unsigned int component)
245 for (
unsigned int d = 0; d < spacedim; ++d)
253 template <
int dim,
int spacedim,
typename Number>
259 typename ::internal::VectorizedArrayTrait<
276 result = vector_entry;
299 const unsigned int vector_lane,
302 for (
unsigned int d = 0; d < dim; ++d)
303 result[d] =
value[d][vector_lane];
316 const unsigned int vector_lane,
319 for (
unsigned int d = 0; d < dim; ++d)
320 value[d][vector_lane] = result[d];
333 const unsigned int vector_lane)
335 for (
unsigned int d = 0; d < spacedim; ++d)
342 const unsigned int vector_lane,
345 result =
value[vector_lane];
358 const unsigned int vector_lane,
361 value[vector_lane] = result;
380 const unsigned int vector_lane,
390 const unsigned int vector_lane,
398 const unsigned int vector_lane,
402 for (
unsigned int d = 0; d < spacedim; ++d)
409 const unsigned int vector_lane,
413 for (
unsigned int d = 0; d < spacedim; ++d)
420 template <
int dim,
typename Number>
425 std::enable_if_t<dim != 1>>
430 typename ::internal::VectorizedArrayTrait<
444 const unsigned int component,
448 result[component] = vector_entry;
462 for (
unsigned int c = 0; c < dim; ++c)
463 result_scalar[c] = result[c].sum();
465 return result_scalar;
473 return result[component].
sum();
478 const unsigned int vector_lane,
481 for (
unsigned int i = 0; i < dim; ++i)
482 for (
unsigned int d = 0; d < dim; ++d)
485 value[d][i], vector_lane);
490 const unsigned int vector_lane,
493 for (
unsigned int i = 0; i < dim; ++i)
494 for (
unsigned int d = 0; d < dim; ++d)
496 value[d][i], vector_lane) = result[i][d];
501 const unsigned int vector_lane)
503 for (
unsigned int i = 0; i < dim; ++i)
504 for (
unsigned int d = 0; d < dim; ++d)
511 const unsigned int vector_lane,
514 for (
unsigned int i = 0; i < dim; ++i)
515 result[i] =
value[i][vector_lane];
528 const unsigned int vector_lane,
531 for (
unsigned int i = 0; i < dim; ++i)
532 value[i][vector_lane] = result[i];
546 for (
unsigned int i = 0; i < dim; ++i)
553 const unsigned int vector_lane,
554 const unsigned int component,
558 vector_lane) += shape_value;
563 const unsigned int vector_lane,
564 const unsigned int component)
572 const unsigned int vector_lane,
573 const unsigned int component,
576 for (
unsigned int d = 0; d < dim; ++d)
584 const unsigned int vector_lane,
585 const unsigned int component)
588 for (
unsigned int d = 0; d < dim; ++d)
596 template <
int dim,
int spacedim>
599 const unsigned int base_element_number);
601 template <
int dim,
int spacedim>
605 template <
int dim,
int spacedim>
606 std::vector<Polynomials::Polynomial<double>>
619template <
int n_components_,
622 typename Number =
double>
635 using ETT =
typename internal::FEPointEvaluation::
636 EvaluatorTypeTraits<dim, spacedim, n_components, Number>;
642 typename ETT::interface_vectorized_unit_gradient_type;
666 const unsigned int first_selected_component = 0);
690 const unsigned int first_selected_component = 0,
782 Tensor<1, (dim == 2 ? 1 : dim), Number>
808 JxW(
const unsigned int point_index)
const;
874 setup(
const unsigned int first_selected_component);
881 template <
bool is_face,
bool is_linear>
914 std::vector<Polynomials::Polynomial<double>>
poly;
1050 std::unique_ptr<NonMatching::MappingInfo<dim, spacedim, Number>>
1130template <
int n_components_,
1133 typename Number =
double>
1147 using ETT =
typename internal::FEPointEvaluation::
1148 EvaluatorTypeTraits<dim, spacedim, n_components, Number>;
1155 typename ETT::interface_vectorized_unit_gradient_type;
1178 const unsigned int first_selected_component = 0);
1199 const unsigned int first_selected_component = 0);
1242 template <std::
size_t str
ide_view>
1285 template <std::
size_t str
ide_view>
1289 const bool sum_into_values =
false);
1316 const bool sum_into_values =
false);
1344 template <std::
size_t str
ide_view>
1349 const bool sum_into_values =
false);
1380 const bool sum_into_values =
false);
1418 template <
bool is_linear, std::
size_t str
ide_view>
1427 template <
bool is_linear, std::
size_t str
ide_view>
1432 const unsigned int n_shapes,
1433 const unsigned int qb,
1440 template <
bool is_linear, std::
size_t str
ide_view>
1449 template <std::
size_t str
ide_view>
1460 template <
bool is_linear>
1464 const unsigned int n_shapes,
1465 const unsigned int qb,
1475 template <
bool is_linear, std::
size_t str
ide_view>
1480 const bool sum_into_values);
1485 template <
bool do_JxW,
bool is_linear, std::
size_t str
ide_view>
1490 const bool sum_into_values);
1495 template <
bool do_JxW, std::
size_t str
ide_view>
1500 const bool sum_into_values);
1505 template <
bool do_JxW, std::
size_t str
ide_view>
1510 const bool sum_into_values);
1540template <
int n_components_,
1543 typename Number =
double>
1557 using ETT =
typename internal::FEPointEvaluation::
1558 EvaluatorTypeTraits<dim, spacedim, n_components, Number>;
1565 typename ETT::interface_vectorized_unit_gradient_type;
1574 const unsigned int first_selected_component = 0);
1588 reinit(
const unsigned int face_index);
1601 template <std::
size_t str
ide_view>
1644 template <std::
size_t str
ide_view>
1648 const bool sum_into_values =
false);
1675 const bool sum_into_values =
false);
1699 template <std::
size_t str
ide_view>
1704 const bool sum_into_values =
false);
1731 const bool sum_into_values =
false);
1739 template <
int str
ide_face_dof = VectorizedArrayType::size()>
1750 template <
int str
ide_face_dof = VectorizedArrayType::size()>
1754 const bool sum_into_values =
false);
1788 template <
bool is_linear, std::
size_t str
ide_view>
1794 template <
bool do_JxW,
bool is_linear, std::
size_t str
ide_view>
1799 const bool sum_into_values);
1805 template <
bool is_linear,
int str
ide_face_dof>
1814 template <
bool do_JxW,
bool is_linear,
int str
ide_face_dof>
1819 const bool sum_into_values);
1827template <
int n_components_,
int dim,
int spacedim,
typename Number>
1832 const unsigned int first_selected_component)
1833 : n_q_batches(
numbers::invalid_unsigned_int)
1834 , n_q_points(
numbers::invalid_unsigned_int)
1835 , n_q_points_scalar(
numbers::invalid_unsigned_int)
1839 , update_flags(update_flags)
1840 , mapping_info_on_the_fly(
1841 std::make_unique<
NonMatching::MappingInfo<dim, spacedim, Number>>(
1844 , mapping_info(mapping_info_on_the_fly.get())
1845 , current_cell_index(
numbers::invalid_unsigned_int)
1846 , current_face_number(
numbers::invalid_unsigned_int)
1847 , must_reinitialize_pointers(false)
1850 setup(first_selected_component);
1855template <
int n_components_,
int dim,
int spacedim,
typename Number>
1860 const unsigned int first_selected_component,
1861 const bool is_interior)
1862 : n_q_batches(
numbers::invalid_unsigned_int)
1863 , n_q_points(
numbers::invalid_unsigned_int)
1864 , n_q_points_scalar(
numbers::invalid_unsigned_int)
1865 , mapping(&mapping_info.get_mapping())
1868 , update_flags(mapping_info.get_update_flags())
1869 , mapping_info(&mapping_info)
1870 , current_cell_index(
numbers::invalid_unsigned_int)
1871 , current_face_number(
numbers::invalid_unsigned_int)
1872 , must_reinitialize_pointers(true)
1873 , is_interior(is_interior)
1875 setup(first_selected_component);
1880template <
int n_components_,
int dim,
int spacedim,
typename Number>
1884 : n_q_batches(other.n_q_batches)
1885 , n_q_points(other.n_q_points)
1886 , n_q_points_scalar(other.n_q_points_scalar)
1887 , mapping(other.mapping)
1890 , use_linear_path(other.use_linear_path)
1891 , renumber(other.renumber)
1892 , solution_renumbered(other.solution_renumbered)
1893 , solution_renumbered_vectorized(other.solution_renumbered_vectorized)
1894 , values(other.values)
1895 , gradients(other.gradients)
1896 , dofs_per_component(other.dofs_per_component)
1897 , dofs_per_component_face(other.dofs_per_component_face)
1898 , component_in_base_element(other.component_in_base_element)
1899 , nonzero_shape_function_component(other.nonzero_shape_function_component)
1900 , update_flags(other.update_flags)
1901 , fe_values(other.fe_values)
1902 , mapping_info_on_the_fly(
1903 other.mapping_info_on_the_fly ?
1908 , mapping_info(other.mapping_info)
1909 , current_cell_index(other.current_cell_index)
1910 , current_face_number(other.current_face_number)
1911 , fast_path(other.fast_path)
1912 , must_reinitialize_pointers(true)
1913 , is_interior(other.is_interior)
1918template <
int n_components_,
int dim,
int spacedim,
typename Number>
1924template <
int n_components_,
int dim,
int spacedim,
typename Number>
1927 const unsigned int first_selected_component)
1930 fe->n_components() + 1);
1932 shapes.reserve(100);
1934 bool same_base_element =
true;
1935 unsigned int base_element_number = 0;
1936 component_in_base_element = 0;
1937 unsigned int component = 0;
1938 for (; base_element_number < fe->n_base_elements(); ++base_element_number)
1939 if (component + fe->element_multiplicity(base_element_number) >
1940 first_selected_component)
1942 if (first_selected_component + n_components >
1943 component + fe->element_multiplicity(base_element_number))
1944 same_base_element =
false;
1945 component_in_base_element = first_selected_component - component;
1949 component += fe->element_multiplicity(base_element_number);
1953 *fe, base_element_number) &&
1956 shape_info.reinit(
QMidpoint<1>(), *fe, base_element_number);
1957 renumber = shape_info.lexicographic_numbering;
1958 dofs_per_component = shape_info.dofs_per_component_on_cell;
1959 dofs_per_component_face = shape_info.dofs_per_component_on_face;
1961 fe->base_element(base_element_number));
1963 bool is_lexicographic =
true;
1964 for (
unsigned int i = 0; i < renumber.size(); ++i)
1965 if (i != renumber[i])
1966 is_lexicographic =
false;
1968 if (is_lexicographic)
1971 use_linear_path = (poly.size() == 2 && poly[0].value(0.) == 1. &&
1972 poly[0].value(1.) == 0. && poly[1].value(0.) == 0. &&
1973 poly[1].value(1.) == 1.) &&
1974 (fe->n_components() == n_components);
1976 const unsigned int size_face = 3 * dofs_per_component_face * n_components;
1977 const unsigned int size_cell = dofs_per_component * n_components;
1978 scratch_data_scalar.resize(size_face + size_cell);
1980 solution_renumbered.resize(dofs_per_component);
1981 solution_renumbered_vectorized.resize(dofs_per_component);
1987 nonzero_shape_function_component.resize(fe->n_dofs_per_cell());
1988 for (
unsigned int d = 0; d < n_components; ++d)
1990 const unsigned int component = first_selected_component + d;
1991 for (
unsigned int i = 0; i < fe->n_dofs_per_cell(); ++i)
1993 const bool is_primitive =
1994 fe->is_primitive() || fe->is_primitive(i);
1996 nonzero_shape_function_component[i][d] =
1997 (component == fe->system_to_component_index(i).first);
1999 nonzero_shape_function_component[i][d] =
2000 (fe->get_nonzero_components(i)[component] ==
true);
2010template <
int n_components_,
int dim,
int spacedim,
typename Number>
2011template <
bool is_face,
bool is_linear>
2015 const unsigned int geometry_index =
2016 mapping_info->template compute_geometry_index_offset<is_face>(
2017 current_cell_index, current_face_number);
2019 cell_type = mapping_info->get_cell_type(geometry_index);
2021 const_cast<unsigned int &
>(n_q_points_scalar) =
2022 mapping_info->get_n_q_points_unvectorized(geometry_index);
2025 const_cast<unsigned int &
>(n_q_batches) =
2026 (n_q_points_scalar + n_lanes_internal - 1) / n_lanes_internal;
2028 const unsigned int n_q_points_before = n_q_points;
2030 const_cast<unsigned int &
>(n_q_points) =
2031 (stride == 1) ? n_q_batches : n_q_points_scalar;
2033 if (n_q_points != n_q_points_before)
2036 values.resize(n_q_points);
2038 gradients.resize(n_q_points);
2041 if (n_q_points == 0)
2045 const unsigned int unit_point_offset =
2046 mapping_info->compute_unit_point_index_offset(geometry_index);
2049 unit_point_faces_ptr =
2050 mapping_info->get_unit_point_faces(unit_point_offset);
2052 unit_point_ptr = mapping_info->get_unit_point(unit_point_offset);
2055 const unsigned int data_offset =
2056 mapping_info->compute_data_index_offset(geometry_index);
2057 const unsigned int compressed_data_offset =
2058 mapping_info->compute_compressed_data_index_offset(geometry_index);
2062 mapping_info->get_update_flags_mapping();
2064 real_point_ptr = mapping_info->get_real_point(data_offset);
2067 mapping_info->get_jacobian(compressed_data_offset, is_interior);
2069 inverse_jacobian_ptr =
2070 mapping_info->get_inverse_jacobian(compressed_data_offset,
2073 normal_ptr = mapping_info->get_normal_vector(data_offset);
2075 JxW_ptr = mapping_info->get_JxW(data_offset);
2079 real_point_ptr = mapping_info->get_real_point(data_offset);
2081 mapping_info->get_jacobian(compressed_data_offset, is_interior);
2082 inverse_jacobian_ptr =
2083 mapping_info->get_inverse_jacobian(compressed_data_offset, is_interior);
2084 normal_ptr = mapping_info->get_normal_vector(data_offset);
2085 JxW_ptr = mapping_info->get_JxW(data_offset);
2088 if (!is_linear && fast_path)
2090 const std::size_t n_shapes = poly.size();
2092 shapes_faces.resize_fast(n_q_batches * n_shapes);
2094 shapes.resize_fast(n_q_batches * n_shapes);
2096 for (
unsigned int qb = 0; qb < n_q_batches; ++qb)
2102 shapes_faces.data() + qb * n_shapes,
2104 unit_point_faces_ptr[qb],
2117 else if (qb + 1 < n_q_batches)
2122 shapes.data() + qb * n_shapes,
2125 unit_point_ptr[qb + 1]);
2141template <
int n_components_,
int dim,
int spacedim,
typename Number>
2145 Number>::value_type &
2147 const unsigned int point_index)
const
2150 return values[point_index];
2155template <
int n_components_,
int dim,
int spacedim,
typename Number>
2159 Number>::gradient_type &
2161 const unsigned int point_index)
const
2164 return gradients[point_index];
2169template <
int n_components_,
int dim,
int spacedim,
typename Number>
2172 const unsigned int point_index)
const
2174 static_assert(n_components == dim,
2175 "Only makes sense for a vector field with dim components");
2178 return trace(gradients[point_index]);
2183template <
int n_components_,
int dim,
int spacedim,
typename Number>
2187 const unsigned int point_index)
2190 values[point_index] = value;
2195template <
int n_components_,
int dim,
int spacedim,
typename Number>
2199 const unsigned int point_index)
2202 gradients[point_index] = gradient;
2207template <
int n_components_,
int dim,
int spacedim,
typename Number>
2210 const Number &value,
2211 const unsigned int point_index)
2213 static_assert(n_components == dim,
2214 "Only makes sense for a vector field with dim components");
2218 for (
unsigned int d = 0; d < dim; ++d)
2219 gradients[point_index][d][d] = value;
2224template <
int n_components_,
int dim,
int spacedim,
typename Number>
2225Tensor<1, (dim == 2 ? 1 : dim), Number>
2227 const unsigned int point_index)
const
2230 dim > 1 && n_components == dim,
2231 "Only makes sense for a vector field with dim components and dim > 1");
2234 Tensor<1, (dim == 2 ? 1 : dim), Number> curl;
2238 curl[0] = grad[1][0] - grad[0][1];
2241 curl[0] = grad[2][1] - grad[1][2];
2242 curl[1] = grad[0][2] - grad[2][0];
2243 curl[2] = grad[1][0] - grad[0][1];
2253template <
int n_components_,
int dim,
int spacedim,
typename Number>
2256 const unsigned int point_index)
const
2259 Assert(jacobian_ptr !=
nullptr,
2261 ExcFEPointEvaluationAccessToUninitializedMappingField(
2262 "update_jacobians"));
2263 return jacobian_ptr[cell_type <= ::internal::MatrixFreeFunctions::
2264 GeometryType::affine ?
2271template <
int n_components_,
int dim,
int spacedim,
typename Number>
2274 const unsigned int point_index)
const
2277 Assert(inverse_jacobian_ptr !=
nullptr,
2279 ExcFEPointEvaluationAccessToUninitializedMappingField(
2280 "update_inverse_jacobians"));
2281 return inverse_jacobian_ptr
2290template <
int n_components_,
int dim,
int spacedim,
typename Number>
2293 const unsigned int point_index)
const
2296 Assert(JxW_ptr !=
nullptr,
2298 ExcFEPointEvaluationAccessToUninitializedMappingField(
2299 "update_JxW_values"));
2300 return JxW_ptr[point_index];
2305template <
int n_components_,
int dim,
int spacedim,
typename Number>
2308 const unsigned int point_index)
const
2310 return quadrature_point(point_index);
2315template <
int n_components_,
int dim,
int spacedim,
typename Number>
2318 const unsigned int point_index)
const
2321 Assert(real_point_ptr !=
nullptr,
2323 ExcFEPointEvaluationAccessToUninitializedMappingField(
2324 "update_quadrature_points"));
2325 return real_point_ptr[point_index];
2330template <
int n_components_,
int dim,
int spacedim,
typename Number>
2333 const unsigned int point_index)
const
2336 Assert(unit_point_ptr !=
nullptr,
ExcMessage(
"unit_point_ptr is not set!"));
2338 for (
unsigned int d = 0; d < dim; ++d)
2340 unit_point_ptr[point_index / stride][d], point_index % stride);
2346template <
int n_components_,
int dim,
int spacedim,
typename Number>
2357template <
int n_components_,
int dim,
int spacedim,
typename Number>
2361 const unsigned int first_selected_component)
2365 first_selected_component)
2370template <
int n_components_,
int dim,
int spacedim,
typename Number>
2375 const unsigned int first_selected_component)
2380 first_selected_component)
2385template <
int n_components_,
int dim,
int spacedim,
typename Number>
2393 if (this->use_linear_path)
2394 this->
template do_reinit<false, true>();
2396 this->
template do_reinit<false, false>();
2401template <
int n_components_,
int dim,
int spacedim,
typename Number>
2405 internal_reinit_single_cell_state_mapping_info();
2406 this->must_reinitialize_pointers =
false;
2411template <
int n_components_,
int dim,
int spacedim,
typename Number>
2418 AssertThrow(this->mapping_info_on_the_fly.get() !=
nullptr,
2421 this->mapping_info_on_the_fly->reinit(cell, unit_points);
2422 this->must_reinitialize_pointers =
false;
2424 if (!this->fast_path)
2426 this->fe_values = std::make_shared<FEValues<dim, spacedim>>(
2430 std::vector<
Point<dim>>(unit_points.begin(), unit_points.end())),
2431 this->update_flags);
2432 this->fe_values->reinit(cell);
2435 if (this->use_linear_path)
2436 this->
template do_reinit<false, true>();
2438 this->
template do_reinit<false, false>();
2443template <
int n_components_,
int dim,
int spacedim,
typename Number>
2450 this->must_reinitialize_pointers =
false;
2452 if (this->use_linear_path)
2453 this->
template do_reinit<false, true>();
2455 this->
template do_reinit<false, false>();
2457 if (!this->fast_path)
2459 std::vector<Point<dim>> unit_points(this->n_q_points_scalar);
2461 for (
unsigned int v = 0; v < this->n_q_points_scalar; ++v)
2462 for (
unsigned int d = 0; d < dim; ++d)
2464 this->unit_point_ptr[v / n_lanes_internal][d][v % n_lanes_internal];
2466 this->fe_values = std::make_shared<FEValues<dim, spacedim>>(
2470 std::vector<
Point<dim>>(unit_points.begin(), unit_points.end())),
2471 this->update_flags);
2473 this->fe_values->reinit(
2474 this->mapping_info->get_cell_iterator(this->current_cell_index));
2480template <
int n_components_,
int dim,
int spacedim,
typename Number>
2481template <std::
size_t str
ide_view>
2493 if (this->must_reinitialize_pointers)
2494 internal_reinit_single_cell_state_mapping_info();
2496 if (this->n_q_points == 0)
2500 if (this->fast_path)
2502 if (this->use_linear_path)
2503 evaluate_fast<true>(solution_values, evaluation_flags);
2505 evaluate_fast<false>(solution_values, evaluation_flags);
2508 evaluate_slow(solution_values, evaluation_flags);
2513template <
int n_components_,
int dim,
int spacedim,
typename Number>
2520 solution_values.
size()),
2526template <
int n_components_,
int dim,
int spacedim,
typename Number>
2527template <std::
size_t str
ide_view>
2532 const bool sum_into_values)
2534 do_integrate<true>(solution_values, integration_flags, sum_into_values);
2539template <
int n_components_,
int dim,
int spacedim,
typename Number>
2544 const bool sum_into_values)
2547 solution_values.
size()),
2554template <
int n_components_,
int dim,
int spacedim,
typename Number>
2562 for (
const auto point_index : this->quadrature_point_indices())
2563 return_value += values[point_index] * this->JxW(point_index);
2565 return ETT::sum_value(return_value);
2570template <
int n_components_,
int dim,
int spacedim,
typename Number>
2577 "Calling this function only makes sense in fully vectorized mode."));
2578 if (q == n_q_batches - 1)
2580 const unsigned int n_filled_lanes =
2581 n_q_points_scalar & (n_lanes_user_interface - 1);
2582 if (n_filled_lanes == 0)
2583 return n_lanes_user_interface;
2585 return n_filled_lanes;
2588 return n_lanes_user_interface;
2593template <
int n_components_,
int dim,
int spacedim,
typename Number>
2594template <std::
size_t str
ide_view>
2599 const bool sum_into_values)
2601 do_integrate<false>(solution_values, integration_flags, sum_into_values);
2606template <
int n_components_,
int dim,
int spacedim,
typename Number>
2611 const bool sum_into_values)
2614 solution_values.
size()),
2621template <
int n_components_,
int dim,
int spacedim,
typename Number>
2622template <
bool is_linear, std::
size_t str
ide_view>
2627 const unsigned int dofs_per_comp =
2630 for (
unsigned int comp = 0; comp < n_components; ++comp)
2632 const std::size_t offset =
2633 (this->component_in_base_element + comp) * dofs_per_comp;
2635 if ((is_linear && n_components == 1) || this->renumber.empty())
2637 for (
unsigned int i = 0; i < dofs_per_comp; ++i)
2638 ETT::read_value(solution_values[i + offset],
2640 this->solution_renumbered[i]);
2644 const unsigned int *renumber_ptr = this->renumber.data() + offset;
2645 for (
unsigned int i = 0; i < dofs_per_comp; ++i)
2646 ETT::read_value(solution_values[renumber_ptr[i]],
2648 this->solution_renumbered[i]);
2655template <
int n_components_,
int dim,
int spacedim,
typename Number>
2656template <
bool is_linear, std::
size_t str
ide_view>
2661 const unsigned int n_shapes,
2662 const unsigned int qb,
2668 std::array<vectorized_value_type, dim + 1> result;
2669 if constexpr (is_linear)
2671 if constexpr (n_components == 1)
2678 stride_view>(solution_values.
data(), this->unit_point_ptr[qb]);
2682 this->solution_renumbered.data(), this->unit_point_ptr[qb]);
2690 false>(this->shapes.data() + qb * n_shapes,
2692 this->solution_renumbered.data());
2693 gradient[0] = result[0];
2695 gradient[1] = result[1];
2697 gradient[2] = result[2];
2698 value = result[dim];
2702 if constexpr (is_linear)
2704 if constexpr (n_components == 1)
2709 stride_view>(solution_values.
data(), this->unit_point_ptr[qb]);
2712 this->solution_renumbered.data(), this->unit_point_ptr[qb]);
2720 this->shapes.data() + qb * n_shapes,
2722 this->solution_renumbered.data());
2728template <
int n_components_,
int dim,
int spacedim,
typename Number>
2729template <
bool is_linear, std::
size_t str
ide_view>
2735 if (!(is_linear && n_components == 1))
2736 prepare_evaluate_fast<is_linear>(solution_values);
2739 const unsigned int n_shapes = is_linear ? 2 : this->poly.size();
2741 for (
unsigned int qb = 0; qb < this->n_q_batches; ++qb)
2746 compute_evaluate_fast<is_linear>(
2747 solution_values, evaluation_flags, n_shapes, qb, value, gradient);
2751 for (
unsigned int v = 0, offset = qb * stride;
2752 v < stride && (stride == 1 || offset < this->n_q_points_scalar);
2754 ETT::set_value(value, v, this->values[offset]);
2762 for (
unsigned int v = 0, offset = qb * stride;
2763 v < stride && (stride == 1 || offset < this->n_q_points_scalar);
2767 ETT::set_gradient(gradient, v, unit_gradient);
2768 this->gradients[offset] =
2772 this->inverse_jacobian_ptr[0].
transpose(), unit_gradient) :
2775 ->inverse_jacobian_ptr[this->cell_type <=
2777 GeometryType::affine ?
2789template <
int n_components_,
int dim,
int spacedim,
typename Number>
2790template <std::
size_t str
ide_view>
2797 Assert(this->fe_values.get() !=
nullptr,
2799 "Not initialized. Please call FEPointEvaluation::reinit()!"));
2801 const std::size_t n_points = this->fe_values->get_quadrature().size();
2805 this->values.resize(this->n_q_points);
2806 std::fill(this->values.begin(), this->values.end(),
value_type());
2807 for (
unsigned int i = 0; i < this->fe->n_dofs_per_cell(); ++i)
2810 for (
unsigned int d = 0; d < n_components; ++d)
2811 if (this->nonzero_shape_function_component[i][d] &&
2812 (this->fe->is_primitive(i) || this->fe->is_primitive()))
2813 for (
unsigned int qb = 0, q = 0; q < n_points;
2814 ++qb, q += n_lanes_user_interface)
2815 for (
unsigned int v = 0;
2816 v < n_lanes_user_interface && q + v < n_points;
2818 ETT::access(this->values[qb],
2821 this->fe_values->shape_value(i, q + v) * value);
2822 else if (this->nonzero_shape_function_component[i][d])
2823 for (
unsigned int qb = 0, q = 0; q < n_points;
2824 ++qb, q += n_lanes_user_interface)
2825 for (
unsigned int v = 0;
2826 v < n_lanes_user_interface && q + v < n_points;
2828 ETT::access(this->values[qb],
2831 this->fe_values->shape_value_component(i,
2840 this->gradients.resize(this->n_q_points);
2841 std::fill(this->gradients.begin(),
2842 this->gradients.end(),
2844 for (
unsigned int i = 0; i < this->fe->n_dofs_per_cell(); ++i)
2847 for (
unsigned int d = 0; d < n_components; ++d)
2848 if (this->nonzero_shape_function_component[i][d] &&
2849 (this->fe->is_primitive(i) || this->fe->is_primitive()))
2850 for (
unsigned int qb = 0, q = 0; q < n_points;
2851 ++qb, q += n_lanes_user_interface)
2852 for (
unsigned int v = 0;
2853 v < n_lanes_user_interface && q + v < n_points;
2855 ETT::access(this->gradients[qb],
2858 this->fe_values->shape_grad(i, q + v) * value);
2859 else if (this->nonzero_shape_function_component[i][d])
2860 for (
unsigned int qb = 0, q = 0; q < n_points;
2861 ++qb, q += n_lanes_user_interface)
2862 for (
unsigned int v = 0;
2863 v < n_lanes_user_interface && q + v < n_points;
2866 this->gradients[qb],
2869 this->fe_values->shape_grad_component(i, q + v, d) * value);
2876template <
int n_components_,
int dim,
int spacedim,
typename Number>
2877template <
bool is_linear>
2881 const unsigned int n_shapes,
2882 const unsigned int qb,
2897 solution_values_vectorized_linear :
2898 this->solution_renumbered_vectorized.data(),
2899 this->unit_point_ptr[qb],
2906 this->shapes.data() + qb * n_shapes,
2909 is_linear ? solution_values_vectorized_linear :
2910 this->solution_renumbered_vectorized.data(),
2911 this->unit_point_ptr[qb],
2917template <
int n_components_,
int dim,
int spacedim,
typename Number>
2918template <
bool is_linear, std::
size_t str
ide_view>
2923 const bool sum_into_values)
2925 if (!sum_into_values && this->fe->n_components() > n_components)
2926 for (
unsigned int i = 0; i < solution_values.
size(); ++i)
2927 solution_values[i] = 0;
2929 const unsigned int dofs_per_comp =
2932 for (
unsigned int comp = 0; comp < n_components; ++comp)
2934 const std::size_t offset =
2935 (this->component_in_base_element + comp) * dofs_per_comp;
2937 if (is_linear || this->renumber.empty())
2939 for (
unsigned int i = 0; i < dofs_per_comp; ++i)
2940 if (sum_into_values)
2941 solution_values[i + offset] +=
2942 ETT::sum_value(comp,
2944 *(solution_values_vectorized_linear + i) :
2945 this->solution_renumbered_vectorized[i]);
2947 solution_values[i + offset] =
2948 ETT::sum_value(comp,
2950 *(solution_values_vectorized_linear + i) :
2951 this->solution_renumbered_vectorized[i]);
2955 const unsigned int *renumber_ptr = this->renumber.
data() + offset;
2956 for (
unsigned int i = 0; i < dofs_per_comp; ++i)
2957 if (sum_into_values)
2958 solution_values[renumber_ptr[i]] +=
2959 ETT::sum_value(comp, this->solution_renumbered_vectorized[i]);
2961 solution_values[renumber_ptr[i]] =
2962 ETT::sum_value(comp, this->solution_renumbered_vectorized[i]);
2969template <
int n_components_,
int dim,
int spacedim,
typename Number>
2970template <
bool do_JxW,
bool is_linear, std::
size_t str
ide_view>
2975 const bool sum_into_values)
2978 if constexpr (stride == 1)
2979 if (
const unsigned int n_filled_lanes =
2980 this->n_q_points_scalar & (n_lanes_internal - 1);
2984 for (
unsigned int v = n_filled_lanes; v < n_lanes_internal; ++v)
2985 ETT::set_zero_value(this->values.back(), v);
2987 for (
unsigned int v = n_filled_lanes; v < n_lanes_internal; ++v)
2988 ETT::set_zero_gradient(this->gradients.back(), v);
2992 solution_values_vectorized_linear = {};
2995 const unsigned int n_shapes = is_linear ? 2 : this->poly.size();
2997 const bool cartesian_cell =
2999 const bool affine_cell =
3001 for (
unsigned int qb = 0; qb < this->n_q_batches; ++qb)
3007 for (
unsigned int v = 0, offset = qb * stride;
3008 v < stride && (stride == 1 || offset < this->n_q_points_scalar);
3010 ETT::get_value(value,
3012 do_JxW ? this->values[offset] * this->JxW_ptr[offset] :
3013 this->values[offset]);
3016 for (
unsigned int v = 0, offset = qb * stride;
3017 v < stride && (stride == 1 || offset < this->n_q_points_scalar);
3021 do_JxW ? this->gradients[offset] * this->JxW_ptr[offset] :
3022 this->gradients[offset];
3030 this->inverse_jacobian_ptr[affine_cell ? 0 : offset],
3034 compute_integrate_fast<is_linear>(
3040 solution_values_vectorized_linear.data());
3044 finish_integrate_fast<is_linear>(solution_values,
3045 solution_values_vectorized_linear.
data(),
3051template <
int n_components_,
int dim,
int spacedim,
typename Number>
3052template <
bool do_JxW, std::
size_t str
ide_view>
3057 const bool sum_into_values)
3060 Assert(this->fe_values.get() !=
nullptr,
3062 "Not initialized. Please call FEPointEvaluation::reinit()!"));
3063 if (!sum_into_values)
3064 for (
unsigned int i = 0; i < solution_values.
size(); ++i)
3065 solution_values[i] = 0;
3067 const std::size_t n_points = this->fe_values->get_quadrature().
size();
3072 for (
unsigned int i = 0; i < this->fe->n_dofs_per_cell(); ++i)
3074 for (
unsigned int d = 0; d < n_components; ++d)
3075 if (this->nonzero_shape_function_component[i][d] &&
3076 (this->fe->is_primitive(i) || this->fe->is_primitive()))
3077 for (
unsigned int qb = 0, q = 0; q < n_points;
3078 ++qb, q += n_lanes_user_interface)
3079 for (
unsigned int v = 0;
3080 v < n_lanes_user_interface && q + v < n_points;
3082 solution_values[i] +=
3083 this->fe_values->shape_value(i, q + v) *
3084 ETT::access(this->values[qb], v, d) *
3085 (do_JxW ? this->fe_values->JxW(q + v) : 1.);
3086 else if (this->nonzero_shape_function_component[i][d])
3087 for (
unsigned int qb = 0, q = 0; q < n_points;
3088 ++qb, q += n_lanes_user_interface)
3089 for (
unsigned int v = 0;
3090 v < n_lanes_user_interface && q + v < n_points;
3092 solution_values[i] +=
3093 this->fe_values->shape_value_component(i, q + v, d) *
3094 ETT::access(this->values[qb], v, d) *
3095 (do_JxW ? this->fe_values->JxW(q + v) : 1.);
3102 for (
unsigned int i = 0; i < this->fe->n_dofs_per_cell(); ++i)
3104 for (
unsigned int d = 0; d < n_components; ++d)
3105 if (this->nonzero_shape_function_component[i][d] &&
3106 (this->fe->is_primitive(i) || this->fe->is_primitive()))
3107 for (
unsigned int qb = 0, q = 0; q < n_points;
3108 ++qb, q += n_lanes_user_interface)
3109 for (
unsigned int v = 0;
3110 v < n_lanes_user_interface && q + v < n_points;
3112 solution_values[i] +=
3113 this->fe_values->shape_grad(i, q + v) *
3114 ETT::access(this->gradients[qb], v, d) *
3115 (do_JxW ? this->fe_values->JxW(q + v) : 1.);
3116 else if (this->nonzero_shape_function_component[i][d])
3117 for (
unsigned int qb = 0, q = 0; q < n_points;
3118 ++qb, q += n_lanes_user_interface)
3119 for (
unsigned int v = 0;
3120 v < n_lanes_user_interface && q + v < n_points;
3122 solution_values[i] +=
3123 this->fe_values->shape_grad_component(i, q + v, d) *
3124 ETT::access(this->gradients[qb], v, d) *
3125 (do_JxW ? this->fe_values->JxW(q + v) : 1.);
3132template <
int n_components_,
int dim,
int spacedim,
typename Number>
3133template <
bool do_JxW, std::
size_t str
ide_view>
3138 const bool sum_into_values)
3140 if (this->must_reinitialize_pointers)
3141 internal_reinit_single_cell_state_mapping_info();
3145 if (this->n_q_points == 0 ||
3147 (integration_flags &
3150 if (!sum_into_values)
3151 for (
unsigned int i = 0; i < solution_values.
size(); ++i)
3152 solution_values[i] = 0;
3157 !do_JxW || this->JxW_ptr !=
nullptr,
3159 "JxW pointer is not set! If you do not want to integrate() use test_and_sum()"));
3162 if (this->fast_path)
3164 if (this->use_linear_path)
3165 integrate_fast<do_JxW, true>(solution_values,
3169 integrate_fast<do_JxW, false>(solution_values,
3174 integrate_slow<do_JxW>(solution_values, integration_flags, sum_into_values);
3179template <
int n_components_,
int dim,
int spacedim,
typename Number>
3182 const unsigned int point_index)
const
3185 Assert(this->normal_ptr !=
nullptr,
3187 ExcFEPointEvaluationAccessToUninitializedMappingField(
3188 "update_normal_vectors"));
3189 return this->normal_ptr[point_index];
3194template <
int n_components_,
int dim,
int spacedim,
typename Number>
3200 const unsigned int point_index)
const
3205 if constexpr (n_components == 1)
3207 this->gradients[point_index] * normal_vector(point_index);
3209 for (
unsigned int comp = 0; comp < n_components; ++comp)
3210 normal_derivative[comp] =
3211 this->gradients[point_index][comp] * normal_vector(point_index);
3213 return normal_derivative;
3218template <
int n_components_,
int dim,
int spacedim,
typename Number>
3222 const unsigned int point_index)
3225 if constexpr (n_components == 1)
3226 this->gradients[point_index] = value * normal_vector(point_index);
3228 for (
unsigned int comp = 0; comp < n_components; ++comp)
3229 this->gradients[point_index][comp] =
3230 value[comp] * normal_vector(point_index);
3235template <
int n_components_,
int dim,
int spacedim,
typename Number>
3240 const bool is_interior,
3241 const unsigned int first_selected_component)
3245 first_selected_component,
3251template <
int n_components_,
int dim,
int spacedim,
typename Number>
3255 const unsigned int face_number)
3258 this->current_face_number = face_number;
3259 this->must_reinitialize_pointers =
false;
3261 if (this->use_linear_path)
3262 this->
template do_reinit<true, true>();
3264 this->
template do_reinit<true, false>();
3269template <
int n_components_,
int dim,
int spacedim,
typename Number>
3272 const unsigned int face_index)
3274 this->current_cell_index = face_index;
3275 this->current_face_number =
3276 this->mapping_info->get_face_number(face_index, this->is_interior);
3277 this->must_reinitialize_pointers =
false;
3279 if (this->use_linear_path)
3280 this->
template do_reinit<true, true>();
3282 this->
template do_reinit<true, false>();
3287template <
int n_components_,
int dim,
int spacedim,
typename Number>
3288template <std::
size_t str
ide_view>
3294 Assert(!this->must_reinitialize_pointers,
3295 ExcMessage(
"Object has not been reinitialized!"));
3297 if (this->n_q_points == 0)
3308 if (this->use_linear_path)
3309 do_evaluate<true>(solution_values, evaluation_flags);
3311 do_evaluate<false>(solution_values, evaluation_flags);
3316template <
int n_components_,
int dim,
int spacedim,
typename Number>
3323 solution_values.
size()),
3329template <
int n_components_,
int dim,
int spacedim,
typename Number>
3330template <
bool is_linear, std::
size_t str
ide_view>
3336 const unsigned int dofs_per_comp =
3340 if (stride_view == 1 && this->component_in_base_element == 0 &&
3341 (is_linear || this->renumber.empty()))
3342 input = solution_values.
data();
3345 for (
unsigned int comp = 0; comp < n_components; ++comp)
3347 const std::size_t offset =
3348 (this->component_in_base_element + comp) * dofs_per_comp;
3350 if (is_linear || this->renumber.empty())
3352 for (
unsigned int i = 0; i < dofs_per_comp; ++i)
3353 this->scratch_data_scalar[i + comp * dofs_per_comp] =
3354 solution_values[i + offset];
3358 const unsigned int *renumber_ptr = this->renumber.
data() + offset;
3359 for (
unsigned int i = 0; i < dofs_per_comp; ++i)
3360 this->scratch_data_scalar[i + comp * dofs_per_comp] =
3361 solution_values[renumber_ptr[i]];
3364 input = this->scratch_data_scalar.
data();
3368 this->scratch_data_scalar.begin() + dofs_per_comp * n_components;
3371 template interpolate<true, false>(n_components,
3376 this->current_face_number);
3378 do_evaluate_in_face<is_linear, 1>(output, evaluation_flags);
3383template <
int n_components_,
int dim,
int spacedim,
typename Number>
3384template <std::
size_t str
ide_view>
3389 const bool sum_into_values)
3391 Assert(!this->must_reinitialize_pointers,
3392 ExcMessage(
"Object has not been reinitialized!"));
3396 if (this->n_q_points == 0 ||
3398 (integration_flags &
3401 if (!sum_into_values)
3402 for (
unsigned int i = 0; i < solution_values.
size(); ++i)
3403 solution_values[i] = 0;
3409 if (this->use_linear_path)
3410 do_integrate<true, true>(solution_values,
3414 do_integrate<true, false>(solution_values,
3421template <
int n_components_,
int dim,
int spacedim,
typename Number>
3426 const bool sum_into_values)
3429 solution_values.
size()),
3436template <
int n_components_,
int dim,
int spacedim,
typename Number>
3437template <std::
size_t str
ide_view>
3442 const bool sum_into_values)
3444 Assert(!this->must_reinitialize_pointers,
3445 ExcMessage(
"Object has not been reinitialized!"));
3449 if (this->n_q_points == 0 ||
3451 (integration_flags &
3454 if (!sum_into_values)
3455 for (
unsigned int i = 0; i < solution_values.
size(); ++i)
3456 solution_values[i] = 0;
3462 if (this->use_linear_path)
3463 do_integrate<false, true>(solution_values,
3467 do_integrate<false, false>(solution_values,
3474template <
int n_components_,
int dim,
int spacedim,
typename Number>
3479 const bool sum_into_values)
3482 solution_values.
size()),
3489template <
int n_components_,
int dim,
int spacedim,
typename Number>
3490template <
bool do_JxW,
bool is_linear, std::
size_t str
ide_view>
3495 const bool sum_into_values)
3497 if (!sum_into_values && this->fe->n_components() > n_components)
3498 for (
unsigned int i = 0; i < solution_values.
size(); ++i)
3499 solution_values[i] = 0;
3501 do_integrate_in_face<do_JxW, is_linear, 1>(this->scratch_data_scalar.begin(),
3505 ScalarNumber *input = this->scratch_data_scalar.begin();
3507 if (stride_view == 1 && this->component_in_base_element == 0 &&
3508 (is_linear || this->renumber.empty()))
3510 if (sum_into_values)
3512 FEFaceNormalEvaluationImpl<dim, is_linear ? 1 : -1,
ScalarNumber>::
3513 template interpolate<false, true>(n_components,
3517 solution_values.
data(),
3518 this->current_face_number);
3521 FEFaceNormalEvaluationImpl<dim, is_linear ? 1 : -1,
ScalarNumber>::
3522 template interpolate<false, false>(n_components,
3526 solution_values.
data(),
3527 this->current_face_number);
3531 const unsigned int dofs_per_comp_face =
3532 is_linear ?
Utilities::pow(2, dim - 1) : this->dofs_per_component_face;
3534 const unsigned int size_input = 3 * dofs_per_comp_face * n_components;
3538 FEFaceNormalEvaluationImpl<dim, is_linear ? 1 : -1,
ScalarNumber>::
3539 template interpolate<false, false>(n_components,
3544 this->current_face_number);
3546 const unsigned int dofs_per_comp =
3549 for (
unsigned int comp = 0; comp < n_components; ++comp)
3551 const std::size_t offset =
3552 (this->component_in_base_element + comp) * dofs_per_comp;
3554 if (is_linear || this->renumber.empty())
3556 for (
unsigned int i = 0; i < dofs_per_comp; ++i)
3557 if (sum_into_values)
3558 solution_values[i + offset] +=
3559 output[i + comp * dofs_per_comp];
3561 solution_values[i + offset] =
3562 output[i + comp * dofs_per_comp];
3566 const unsigned int *renumber_ptr = this->renumber.
data() + offset;
3567 for (
unsigned int i = 0; i < dofs_per_comp; ++i)
3568 if (sum_into_values)
3569 solution_values[renumber_ptr[i]] +=
3570 output[i + comp * dofs_per_comp];
3572 solution_values[renumber_ptr[i]] =
3573 output[i + comp * dofs_per_comp];
3581template <
int n_components_,
int dim,
int spacedim,
typename Number>
3582template <
int str
ide_face_dof>
3588 if (this->use_linear_path)
3589 do_evaluate_in_face<true, stride_face_dof>(face_dof_values,
3592 do_evaluate_in_face<false, stride_face_dof>(face_dof_values,
3598template <
int n_components_,
int dim,
int spacedim,
typename Number>
3599template <
bool is_linear,
int str
ide_face_dof>
3606 if constexpr (n_components == 1)
3607 face_dof_values_ptr = face_dof_values;
3610 const unsigned int dofs_per_comp_face =
3611 is_linear ?
Utilities::pow(2, dim - 1) : this->dofs_per_component_face;
3612 for (
unsigned int comp = 0; comp < n_components; ++comp)
3613 for (
unsigned int i = 0; i < 2 * dofs_per_comp_face; ++i)
3614 ETT::read_value(face_dof_values[(i + comp * 3 * dofs_per_comp_face) *
3617 this->solution_renumbered[i]);
3619 face_dof_values_ptr = this->solution_renumbered.data();
3622 constexpr int stride_face_dof_actual =
3623 n_components == 1 ? stride_face_dof : 1;
3626 const unsigned int n_shapes = is_linear ? 2 : this->poly.size();
3628 for (
unsigned int qb = 0; qb < this->n_q_batches; ++qb)
3635 const std::array<vectorized_value_type, dim + 1> interpolated_value =
3642 stride_face_dof_actual>(face_dof_values_ptr,
3643 this->unit_point_faces_ptr[qb]) :
3650 stride_face_dof_actual>(this->shapes_faces.data() +
3653 face_dof_values_ptr);
3655 value = interpolated_value[dim - 1];
3658 if (this->current_face_number / 2 == 0)
3660 gradient[0] = interpolated_value[dim];
3662 gradient[1] = interpolated_value[0];
3664 gradient[2] = interpolated_value[1];
3666 else if (this->current_face_number / 2 == 1)
3669 gradient[1] = interpolated_value[dim];
3672 gradient[0] = interpolated_value[1];
3673 gradient[2] = interpolated_value[0];
3676 gradient[0] = interpolated_value[0];
3680 else if (this->current_face_number / 2 == 2)
3684 gradient[0] = interpolated_value[0];
3685 gradient[1] = interpolated_value[1];
3686 gradient[2] = interpolated_value[dim];
3701 stride_face_dof_actual>(face_dof_values_ptr,
3702 this->unit_point_faces_ptr[qb]) :
3708 stride_face_dof_actual>(this->shapes_faces.data() +
3711 face_dof_values_ptr);
3716 for (
unsigned int v = 0, offset = qb * stride;
3717 v < stride && (stride == 1 || offset < this->n_q_points_scalar);
3719 ETT::set_value(value, v, this->values[offset]);
3727 for (
unsigned int v = 0, offset = qb * stride;
3728 v < stride && (stride == 1 || offset < this->n_q_points_scalar);
3732 ETT::set_gradient(gradient, v, unit_gradient);
3733 this->gradients[offset] =
3737 this->inverse_jacobian_ptr[0].
transpose(), unit_gradient) :
3740 ->inverse_jacobian_ptr[this->cell_type <=
3742 GeometryType::affine ?
3754template <
int n_components_,
int dim,
int spacedim,
typename Number>
3755template <
int str
ide_face_dof>
3760 const bool sum_into_values)
3762 if (this->use_linear_path)
3763 do_integrate_in_face<true, true, stride_face_dof>(face_dof_values,
3767 do_integrate_in_face<true, false, stride_face_dof>(face_dof_values,
3774template <
int n_components_,
int dim,
int spacedim,
typename Number>
3775template <
bool do_JxW,
bool is_linear,
int str
ide_face_dof>
3781 const bool sum_into_values)
3784 if constexpr (stride == 1)
3785 if (
const unsigned int n_filled_lanes =
3786 this->n_q_points_scalar & (n_lanes_internal - 1);
3790 for (
unsigned int v = n_filled_lanes; v < n_lanes_internal; ++v)
3791 ETT::set_zero_value(this->values.back(), v);
3793 for (
unsigned int v = n_filled_lanes; v < n_lanes_internal; ++v)
3794 ETT::set_zero_gradient(this->gradients.back(), v);
3799 solution_values_vectorized_linear = {};
3802 const unsigned int n_shapes = is_linear ? 2 : this->poly.size();
3804 const bool cartesian_cell =
3806 const bool affine_cell =
3808 for (
unsigned int qb = 0; qb < this->n_q_batches; ++qb)
3814 for (
unsigned int v = 0, offset = qb * stride;
3815 v < stride && (stride == 1 || offset < this->n_q_points_scalar);
3817 ETT::get_value(value,
3819 do_JxW ? this->values[offset] * this->JxW_ptr[offset] :
3820 this->values[offset]);
3823 for (
unsigned int v = 0, offset = qb * stride;
3824 v < stride && (stride == 1 || offset < this->n_q_points_scalar);
3828 do_JxW ? this->gradients[offset] * this->JxW_ptr[offset] :
3829 this->gradients[offset];
3837 this->inverse_jacobian_ptr[affine_cell ? 0 : offset],
3843 std::array<vectorized_value_type, 2> value_face = {};
3846 value_face[0] = value;
3849 if (this->current_face_number / 2 == 0)
3851 value_face[1] = gradient[0];
3853 gradient_in_face[0] = gradient[1];
3855 gradient_in_face[1] = gradient[2];
3857 else if (this->current_face_number / 2 == 1)
3860 value_face[1] = gradient[1];
3863 gradient_in_face[0] = gradient[2];
3864 gradient_in_face[1] = gradient[0];
3867 gradient_in_face[0] = gradient[0];
3871 else if (this->current_face_number / 2 == 2)
3875 value_face[1] = gradient[2];
3876 gradient_in_face[0] = gradient[0];
3877 gradient_in_face[1] = gradient[1];
3890 2>(this->shapes_faces.data() + qb * n_shapes,
3894 is_linear ? solution_values_vectorized_linear.data() :
3895 this->solution_renumbered_vectorized.data(),
3896 this->unit_point_faces_ptr[qb],
3904 this->shapes_faces.data() + qb * n_shapes,
3907 is_linear ? solution_values_vectorized_linear.data() :
3908 this->solution_renumbered_vectorized.data(),
3909 this->unit_point_faces_ptr[qb],
3913 const unsigned int dofs_per_comp_face =
3914 is_linear ?
Utilities::pow(2, dim - 1) : this->dofs_per_component_face;
3916 for (
unsigned int comp = 0; comp < n_components; ++comp)
3917 for (
unsigned int i = 0; i < 2 * dofs_per_comp_face; ++i)
3918 if (sum_into_values)
3919 face_dof_values[(i + comp * 3 * dofs_per_comp_face) *
3921 ETT::sum_value(comp,
3923 *(solution_values_vectorized_linear.data() + i) :
3924 this->solution_renumbered_vectorized[i]);
3926 face_dof_values[(i + comp * 3 * dofs_per_comp_face) * stride_face_dof] =
3927 ETT::sum_value(comp,
3929 *(solution_values_vectorized_linear.data() + i) :
3930 this->solution_renumbered_vectorized[i]);
3935template <
int n_components_,
int dim,
int spacedim,
typename Number>
3938 const unsigned int point_index)
const
3941 Assert(this->normal_ptr !=
nullptr,
3943 ExcFEPointEvaluationAccessToUninitializedMappingField(
3944 "update_normal_vectors"));
3948 for (
unsigned int d = 0; d < dim; ++d)
3957 return this->normal_ptr[point_index];
3963template <
int n_components_,
int dim,
int spacedim,
typename Number>
3974 if constexpr (n_components == 1)
3976 this->gradients[point_index] * normal_vector(point_index);
3978 for (
unsigned int comp = 0; comp < n_components; ++comp)
3979 normal_derivative[comp] =
3980 this->gradients[point_index][comp] * normal_vector(point_index);
3982 return normal_derivative;
3987template <
int n_components_,
int dim,
int spacedim,
typename Number>
3991 const unsigned int point_index)
3994 if constexpr (n_components == 1)
3995 this->gradients[point_index] = value * normal_vector(point_index);
3997 for (
unsigned int comp = 0; comp < n_components; ++comp)
3998 this->gradients[point_index][comp] =
3999 value[comp] * normal_vector(point_index);
874 setup(
const unsigned int first_selected_component); {
…}
808 JxW(
const unsigned int point_index)
const; {
…}
value_type * data() const noexcept
typename ETT::vectorized_value_type vectorized_value_type
static constexpr std::size_t n_lanes_internal
void do_integrate_in_face(ScalarNumber *face_dof_values, const EvaluationFlags::EvaluationFlags &integration_flags, const bool sum_into_values)
void integrate_in_face(ScalarNumber *face_dof_values, const EvaluationFlags::EvaluationFlags &integration_flags, const bool sum_into_values=false)
typename ETT::scalar_value_type scalar_value_type
typename internal::VectorizedArrayTrait< Number >::value_type ScalarNumber
void evaluate(const StridedArrayView< const ScalarNumber, stride_view > &solution_values, const EvaluationFlags::EvaluationFlags &evaluation_flags)
void do_evaluate_in_face(const ScalarNumber *face_dof_values, const EvaluationFlags::EvaluationFlags &evaluation_flags)
void submit_normal_derivative(const value_type &, const unsigned int point_index)
typename internal::FEPointEvaluation::EvaluatorTypeTraits< dim, spacedim, n_components, Number > ETT
static constexpr std::size_t stride
typename ETT::interface_vectorized_unit_gradient_type interface_vectorized_unit_gradient_type
void do_evaluate(const StridedArrayView< const ScalarNumber, stride_view > &solution_values, const EvaluationFlags::EvaluationFlags &evaluation_flags)
static constexpr unsigned int n_components
static constexpr std::size_t n_lanes_user_interface
typename ETT::real_gradient_type gradient_type
FEFacePointEvaluation(const NonMatching::MappingInfo< dim, spacedim, Number > &mapping_info, const FiniteElement< dim, spacedim > &fe, const bool is_interior=true, const unsigned int first_selected_component=0)
void evaluate_in_face(const ScalarNumber *face_dof_values, const EvaluationFlags::EvaluationFlags &evaluation_flags)
void reinit(const unsigned int cell_index, const unsigned int face_number)
void test_and_sum(const StridedArrayView< ScalarNumber, stride_view > &solution_values, const EvaluationFlags::EvaluationFlags &integration_flags, const bool sum_into_values=false)
typename ETT::unit_gradient_type unit_gradient_type
const value_type get_normal_derivative(const unsigned int point_index) const
void do_integrate(const StridedArrayView< ScalarNumber, stride_view > &solution_values, const EvaluationFlags::EvaluationFlags &integration_flags, const bool sum_into_values)
void integrate(const StridedArrayView< ScalarNumber, stride_view > &solution_values, const EvaluationFlags::EvaluationFlags &integration_flags, const bool sum_into_values=false)
static constexpr unsigned int dimension
Tensor< 1, spacedim, Number > normal_vector(const unsigned int point_index) const
typename ETT::value_type value_type
typename ::internal::VectorizedArrayTrait< Number >::vectorized_value_type VectorizedArrayType
ObserverPointer< const Mapping< dim, spacedim > > mapping
const DerivativeForm< 1, dim, spacedim, Number > * jacobian_ptr
std::unique_ptr< NonMatching::MappingInfo< dim, spacedim, Number > > mapping_info_on_the_fly
std::vector< gradient_type > gradients
unsigned int dofs_per_component
typename ETT::interface_vectorized_unit_gradient_type interface_vectorized_unit_gradient_type
Number get_divergence(const unsigned int point_index) const
unsigned int dofs_per_component_face
Number JxW(const unsigned int point_index) const
const UpdateFlags update_flags
static constexpr std::size_t n_lanes_user_interface
static constexpr std::size_t stride
internal::MatrixFreeFunctions::GeometryType cell_type
std::vector< Polynomials::Polynomial< double > > poly
Point< spacedim, Number > real_point(const unsigned int point_index) const
Tensor< 1,(dim==2 ? 1 :dim), Number > get_curl(const unsigned int point_index) const
const unsigned int n_q_points_scalar
const Point< dim, VectorizedArrayType > * unit_point_ptr
FEPointEvaluationBase(FEPointEvaluationBase &&other) noexcept
AlignedVector< ScalarNumber > scratch_data_scalar
std_cxx20::ranges::iota_view< unsigned int, unsigned int > quadrature_point_indices() const
const value_type & get_value(const unsigned int point_index) const
const Point< dim - 1, VectorizedArrayType > * unit_point_faces_ptr
std::vector< scalar_value_type > solution_renumbered
const unsigned int n_q_batches
ObserverPointer< const FiniteElement< dim, spacedim > > fe
std::vector< std::array< bool, n_components > > nonzero_shape_function_component
Point< spacedim, Number > quadrature_point(const unsigned int point_index) const
unsigned int n_active_entries_per_quadrature_batch(unsigned int q)
const gradient_type & get_gradient(const unsigned int point_index) const
std::shared_ptr< FEValues< dim, spacedim > > fe_values
FEPointEvaluationBase(const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const UpdateFlags update_flags, const unsigned int first_selected_component=0)
std::vector< unsigned int > renumber
std::vector< value_type > values
AlignedVector<::ndarray< VectorizedArrayType, 2, dim - 1 > > shapes_faces
Point< dim, Number > unit_point(const unsigned int point_index) const
const unsigned int n_q_points
void submit_divergence(const Number &value, const unsigned int point_index)
bool must_reinitialize_pointers
const Tensor< 1, spacedim, Number > * normal_ptr
DerivativeForm< 1, spacedim, dim, Number > inverse_jacobian(const unsigned int point_index) const
AlignedVector< vectorized_value_type > solution_renumbered_vectorized
static constexpr std::size_t n_lanes_internal
ObserverPointer< const NonMatching::MappingInfo< dim, spacedim, Number > > mapping_info
const DerivativeForm< 1, spacedim, dim, Number > * inverse_jacobian_ptr
AlignedVector<::ndarray< VectorizedArrayType, 2, dim > > shapes
typename ETT::value_type value_type
unsigned int current_cell_index
scalar_value_type integrate_value() const
void submit_gradient(const gradient_type &, const unsigned int point_index)
void setup(const unsigned int first_selected_component)
typename ETT::scalar_value_type scalar_value_type
static constexpr unsigned int dimension
typename ::internal::VectorizedArrayTrait< Number >::vectorized_value_type VectorizedArrayType
typename ETT::vectorized_value_type vectorized_value_type
static constexpr unsigned int n_components
typename internal::FEPointEvaluation::EvaluatorTypeTraits< dim, spacedim, n_components, Number > ETT
unsigned int current_face_number
FEPointEvaluationBase(FEPointEvaluationBase &other) noexcept
DerivativeForm< 1, dim, spacedim, Number > jacobian(const unsigned int point_index) const
typename ETT::real_gradient_type gradient_type
const Point< spacedim, Number > * real_point_ptr
FEPointEvaluationBase(const NonMatching::MappingInfo< dim, spacedim, Number > &mapping_info, const FiniteElement< dim, spacedim > &fe, const unsigned int first_selected_component=0, const bool is_interior=true)
unsigned int component_in_base_element
internal::MatrixFreeFunctions::ShapeInfo< ScalarNumber > shape_info
void submit_value(const value_type &value, const unsigned int point_index)
typename internal::VectorizedArrayTrait< Number >::value_type ScalarNumber
static constexpr std::size_t stride
void integrate_slow(const StridedArrayView< ScalarNumber, stride_view > &solution_values, const EvaluationFlags::EvaluationFlags &integration_flags, const bool sum_into_values)
static constexpr unsigned int dimension
typename ETT::value_type value_type
void do_integrate(const StridedArrayView< ScalarNumber, stride_view > &solution_values, const EvaluationFlags::EvaluationFlags &integration_flags, const bool sum_into_values)
void evaluate_fast(const StridedArrayView< const ScalarNumber, stride_view > &solution_values, const EvaluationFlags::EvaluationFlags &evaluation_flags)
typename ETT::scalar_value_type scalar_value_type
typename ETT::unit_gradient_type unit_gradient_type
void evaluate(const StridedArrayView< const ScalarNumber, stride_view > &solution_values, const EvaluationFlags::EvaluationFlags &evaluation_flags)
const value_type get_normal_derivative(const unsigned int point_index) const
void prepare_evaluate_fast(const StridedArrayView< const ScalarNumber, stride_view > &solution_values)
void test_and_sum(const StridedArrayView< ScalarNumber, stride_view > &solution_values, const EvaluationFlags::EvaluationFlags &integration_flags, const bool sum_into_values=false)
void integrate(const StridedArrayView< ScalarNumber, stride_view > &solution_values, const EvaluationFlags::EvaluationFlags &integration_flags, const bool sum_into_values=false)
typename internal::VectorizedArrayTrait< Number >::value_type ScalarNumber
typename ETT::vectorized_value_type vectorized_value_type
void compute_integrate_fast(const EvaluationFlags::EvaluationFlags &integration_flags, const unsigned int n_shapes, const unsigned int qb, const vectorized_value_type value, const interface_vectorized_unit_gradient_type gradient, vectorized_value_type *solution_values_vectorized_linear)
typename ::internal::VectorizedArrayTrait< Number >::vectorized_value_type VectorizedArrayType
typename ETT::real_gradient_type gradient_type
FEPointEvaluation(const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const UpdateFlags update_flags, const unsigned int first_selected_component=0)
void submit_normal_derivative(const value_type &, const unsigned int point_index)
void evaluate_slow(const StridedArrayView< const ScalarNumber, stride_view > &solution_values, const EvaluationFlags::EvaluationFlags &evaluation_flags)
static constexpr std::size_t n_lanes_internal
static constexpr unsigned int n_components
void finish_integrate_fast(const StridedArrayView< ScalarNumber, stride_view > &solution_values, vectorized_value_type *solution_values_vectorized_linear, const bool sum_into_values)
typename internal::FEPointEvaluation::EvaluatorTypeTraits< dim, spacedim, n_components, Number > ETT
void internal_reinit_single_cell_state_mapping_info()
void compute_evaluate_fast(const StridedArrayView< const ScalarNumber, stride_view > &solution_values, const EvaluationFlags::EvaluationFlags &evaluation_flags, const unsigned int n_shapes, const unsigned int qb, vectorized_value_type &value, interface_vectorized_unit_gradient_type &gradient)
void integrate_fast(const StridedArrayView< ScalarNumber, stride_view > &solution_values, const EvaluationFlags::EvaluationFlags &integration_flags, const bool sum_into_values)
static constexpr std::size_t n_lanes_user_interface
typename ETT::interface_vectorized_unit_gradient_type interface_vectorized_unit_gradient_type
Tensor< 1, spacedim, Number > normal_vector(const unsigned int point_index) const
Abstract base class for mapping classes.
value_type * data() const noexcept
Tensor< rank, dim, Number > sum(const Tensor< rank, dim, Number > &local, const MPI_Comm mpi_communicator)
#define DEAL_II_DEPRECATED
#define DEAL_II_NAMESPACE_OPEN
constexpr bool running_in_debug_mode()
#define DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_ASSERT_UNREACHABLE()
#define DEAL_II_NOT_IMPLEMENTED()
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcFEPointEvaluationAccessToUninitializedMappingField(std::string arg1)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcNotInitialized()
#define DeclException1(Exception1, type1, outsequence)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
@ update_values
Shape function values.
@ update_normal_vectors
Normal vectors.
@ update_JxW_values
Transformed quadrature weights.
@ update_jacobians
Volume element.
@ update_inverse_jacobians
Volume element.
@ update_gradients
Shape function gradients.
@ update_quadrature_points
Transformed quadrature points.
EvaluationFlags
The EvaluationFlags enum.
constexpr T pow(const T base, const int iexp)
std::vector< Polynomials::Polynomial< double > > get_polynomial_space(const FiniteElement< dim, spacedim > &fe)
bool is_fast_path_supported(const FiniteElement< dim, spacedim > &fe, const unsigned int base_element_number)
std::array< typename ProductTypeNoPoint< Number, Number2 >::type, dim+n_values > evaluate_tensor_product_value_and_gradient_shapes(const ::ndarray< Number2, 2, dim > *shapes, const int n_shapes, const Number *values, const std::vector< unsigned int > &renumber={})
void integrate_tensor_product_value_and_gradient(const ::ndarray< Number, 2, dim > *shapes, const unsigned int n_shapes, const Number2 *value, const Tensor< 1, dim, Number2 > &gradient, Number2 *values, const Point< dim, Number > &p, const bool do_add)
void compute_values_of_array_in_pairs(::ndarray< Number, 2, dim > *shapes, const std::vector< Polynomials::Polynomial< double > > &poly, const Point< dim, Number > &p0, const Point< dim, Number > &p1)
ProductTypeNoPoint< Number, Number2 >::type evaluate_tensor_product_value_shapes(const ::ndarray< Number2, 2, dim > *shapes, const int n_shapes, const Number *values, const std::vector< unsigned int > &renumber={})
void compute_values_of_array(::ndarray< Number, 2, dim > *shapes, const std::vector< Polynomials::Polynomial< double > > &poly, const Point< dim, Number > &p, const unsigned int derivative=1)
ProductTypeNoPoint< Number, Number2 >::type evaluate_tensor_product_value_linear(const Number *values, const Point< dim, Number2 > &p)
std::array< typename ProductTypeNoPoint< Number, Number2 >::type, dim+n_values > evaluate_tensor_product_value_and_gradient_linear(const Number *values, const Point< dim, Number2 > &p)
void integrate_tensor_product_value(const ::ndarray< Number, 2, dim > *shapes, const unsigned int n_shapes, const Number2 &value, Number2 *values, const Point< dim, Number > &p, const bool do_add)
constexpr unsigned int invalid_unsigned_int
boost::integer_range< IncrementableType > iota_view
typename internal::ndarray::HelperArray< T, Ns... >::type ndarray
static void set_zero_gradient(unit_gradient_type &value, const unsigned int vector_lane)
static void access(value_type &value, const unsigned int vector_lane, const unsigned int component, const ScalarNumber &shape_value)
static void set_value(const vectorized_value_type &value, const unsigned int vector_lane, scalar_value_type &result)
static Tensor< 1, dim, ScalarNumber > access(const real_gradient_type &value, const unsigned int vector_lane, const unsigned int component)
static void get_value(vectorized_value_type &value, const unsigned int, const vectorized_value_type &result)
static void set_zero_value(value_type &value, const unsigned int vector_lane)
static void set_gradient(const interface_vectorized_unit_gradient_type &value, const unsigned int vector_lane, unit_gradient_type &result)
static void read_value(const ScalarNumber vector_entry, const unsigned int component, scalar_value_type &result)
static void access(real_gradient_type &value, const unsigned int vector_lane, const unsigned int component, const Tensor< 1, dim, ScalarNumber > &shape_gradient)
static scalar_value_type sum_value(const vectorized_value_type &result)
static scalar_value_type sum_value(const scalar_value_type &result)
static ScalarNumber sum_value(const unsigned int component, const vectorized_value_type &result)
static ScalarNumber access(const value_type &value, const unsigned int vector_lane, const unsigned int component)
static void get_value(vectorized_value_type &value, const unsigned int vector_lane, const scalar_value_type &result)
typename ::internal::VectorizedArrayTrait< Number >::vectorized_value_type VectorizedArrayType
static void get_gradient(interface_vectorized_unit_gradient_type &value, const unsigned int vector_lane, const unit_gradient_type &result)
static void set_value(const vectorized_value_type &value, const unsigned int, vectorized_value_type &result)
typename internal::VectorizedArrayTrait< Number >::value_type ScalarNumber
static ScalarNumber access(const value_type &value, const unsigned int vector_lane, const unsigned int)
static void get_gradient(vectorized_unit_gradient_type &value, const unsigned int, const vectorized_unit_gradient_type &result)
static scalar_value_type sum_value(const scalar_value_type &result)
static void get_value(vectorized_value_type &value, const unsigned int, const vectorized_value_type &result)
static void set_gradient(const vectorized_unit_gradient_type &value, const unsigned int vector_lane, scalar_unit_gradient_type &result)
static void set_zero_value(value_type &value, const unsigned int vector_lane)
static Tensor< 1, spacedim, ScalarNumber > access(const real_gradient_type &value, const unsigned int vector_lane, const unsigned int)
static void set_value(const vectorized_value_type &value, const unsigned int, vectorized_value_type &result)
VectorizedArrayType vectorized_value_type
typename internal::VectorizedArrayTrait< Number >::value_type ScalarNumber
static void get_value(vectorized_value_type &value, const unsigned int vector_lane, const scalar_value_type &result)
static scalar_value_type sum_value(const vectorized_value_type &result)
static void access(real_gradient_type &value, const unsigned int vector_lane, const unsigned int, const Tensor< 1, spacedim, ScalarNumber > &shape_gradient)
static void set_gradient(const vectorized_unit_gradient_type &value, const unsigned int, vectorized_unit_gradient_type &result)
ScalarNumber scalar_value_type
static void set_value(const vectorized_value_type &value, const unsigned int vector_lane, scalar_value_type &result)
static void set_zero_gradient(real_gradient_type &value, const unsigned int vector_lane)
static void access(value_type &value, const unsigned int vector_lane, const unsigned int, const ScalarNumber &shape_value)
static void get_gradient(vectorized_unit_gradient_type &value, const unsigned int vector_lane, const scalar_unit_gradient_type &result)
static ScalarNumber sum_value(const unsigned int, const vectorized_value_type &result)
static void read_value(const ScalarNumber vector_entry, const unsigned int, scalar_value_type &result)
typename ::internal::VectorizedArrayTrait< Number >::vectorized_value_type VectorizedArrayType
static Tensor< 1, spacedim, ScalarNumber > access(const real_gradient_type &value, const unsigned int vector_lane, const unsigned int component)
static void get_gradient(interface_vectorized_unit_gradient_type &value, const unsigned int vector_lane, const DerivativeForm< 1, dim, n_components, Number > &result)
static void get_value(vectorized_value_type &value, const unsigned int, const vectorized_value_type &result)
static void get_value(vectorized_value_type &value, const unsigned int vector_lane, const scalar_value_type &result)
static void set_value(const vectorized_value_type &value, const unsigned int, vectorized_value_type &result)
static scalar_value_type sum_value(const scalar_value_type &result)
typename ::internal::VectorizedArrayTrait< Number >::vectorized_value_type VectorizedArrayType
static scalar_value_type sum_value(const vectorized_value_type &result)
static void read_value(const ScalarNumber vector_entry, const unsigned int component, scalar_value_type &result)
Tensor< 1, n_components, ScalarNumber > scalar_value_type
static void set_value(const vectorized_value_type &value, const unsigned int vector_lane, scalar_value_type &result)
static ScalarNumber access(const value_type &value, const unsigned int vector_lane, const unsigned int component)
static ScalarNumber sum_value(const unsigned int component, const vectorized_value_type &result)
Tensor< 1, n_components, Tensor< 1, dim, VectorizedArrayType > > vectorized_unit_gradient_type
static void set_zero_value(value_type &value, const unsigned int vector_lane)
static void set_gradient(const interface_vectorized_unit_gradient_type &value, const unsigned int vector_lane, unit_gradient_type &result)
static void set_zero_gradient(real_gradient_type &value, const unsigned int vector_lane)
typename internal::VectorizedArrayTrait< Number >::value_type ScalarNumber
static void get_gradient(interface_vectorized_unit_gradient_type &value, const unsigned int vector_lane, const unit_gradient_type &result)
std::conditional_t< n_components==spacedim, Tensor< 2, spacedim, Number >, Tensor< 1, n_components, Tensor< 1, spacedim, Number > > > real_gradient_type
Tensor< 1, n_components, VectorizedArrayType > vectorized_value_type
static void access(real_gradient_type &value, const unsigned int vector_lane, const unsigned int component, const Tensor< 1, spacedim, ScalarNumber > &shape_gradient)
Tensor< 1, n_components, Tensor< 1, dim, Number > > unit_gradient_type
static void access(value_type &value, const unsigned int vector_lane, const unsigned int component, const ScalarNumber &shape_value)
static constexpr std::size_t width()
static constexpr std::size_t stride()
static value_type & get(value_type &value, unsigned int c)
DEAL_II_HOST constexpr Number trace(const SymmetricTensor< 2, dim2, Number > &)