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>
1921 : n_q_batches(other.n_q_batches)
1922 , n_q_points(other.n_q_points)
1923 , n_q_points_scalar(other.n_q_points_scalar)
1924 , mapping(other.mapping)
1927 , use_linear_path(other.use_linear_path)
1928 , renumber(other.renumber)
1929 , solution_renumbered(other.solution_renumbered)
1930 , solution_renumbered_vectorized(other.solution_renumbered_vectorized)
1931 , values(other.values)
1932 , gradients(other.gradients)
1933 , dofs_per_component(other.dofs_per_component)
1934 , dofs_per_component_face(other.dofs_per_component_face)
1935 , component_in_base_element(other.component_in_base_element)
1936 , nonzero_shape_function_component(other.nonzero_shape_function_component)
1937 , update_flags(other.update_flags)
1938 , fe_values(other.fe_values)
1939 , mapping_info_on_the_fly(std::move(other.mapping_info_on_the_fly))
1940 , mapping_info(other.mapping_info)
1941 , current_cell_index(other.current_cell_index)
1942 , current_face_number(other.current_face_number)
1943 , fast_path(other.fast_path)
1944 , must_reinitialize_pointers(other.must_reinitialize_pointers)
1945 , is_interior(other.is_interior)
1950template <
int n_components_,
int dim,
int spacedim,
typename Number>
1953 const unsigned int first_selected_component)
1956 fe->n_components() + 1);
1958 shapes.reserve(100);
1960 bool same_base_element =
true;
1961 unsigned int base_element_number = 0;
1962 component_in_base_element = 0;
1963 unsigned int component = 0;
1964 for (; base_element_number < fe->n_base_elements(); ++base_element_number)
1965 if (component + fe->element_multiplicity(base_element_number) >
1966 first_selected_component)
1968 if (first_selected_component + n_components >
1969 component + fe->element_multiplicity(base_element_number))
1970 same_base_element =
false;
1971 component_in_base_element = first_selected_component - component;
1975 component += fe->element_multiplicity(base_element_number);
1979 *fe, base_element_number) &&
1982 shape_info.reinit(
QMidpoint<1>(), *fe, base_element_number);
1983 renumber = shape_info.lexicographic_numbering;
1984 dofs_per_component = shape_info.dofs_per_component_on_cell;
1985 dofs_per_component_face = shape_info.dofs_per_component_on_face;
1987 fe->base_element(base_element_number));
1989 bool is_lexicographic =
true;
1990 for (
unsigned int i = 0; i < renumber.size(); ++i)
1991 if (i != renumber[i])
1992 is_lexicographic =
false;
1994 if (is_lexicographic)
1997 use_linear_path = (poly.size() == 2 && poly[0].value(0.) == 1. &&
1998 poly[0].value(1.) == 0. && poly[1].value(0.) == 0. &&
1999 poly[1].value(1.) == 1.) &&
2000 (fe->n_components() == n_components);
2002 const unsigned int size_face = 3 * dofs_per_component_face * n_components;
2003 const unsigned int size_cell = dofs_per_component * n_components;
2004 scratch_data_scalar.resize(size_face + size_cell);
2006 solution_renumbered.resize(dofs_per_component);
2007 solution_renumbered_vectorized.resize(dofs_per_component);
2013 nonzero_shape_function_component.resize(fe->n_dofs_per_cell());
2014 for (
unsigned int d = 0; d < n_components; ++d)
2016 const unsigned int component = first_selected_component + d;
2017 for (
unsigned int i = 0; i < fe->n_dofs_per_cell(); ++i)
2019 const bool is_primitive =
2020 fe->is_primitive() || fe->is_primitive(i);
2022 nonzero_shape_function_component[i][d] =
2023 (component == fe->system_to_component_index(i).first);
2025 nonzero_shape_function_component[i][d] =
2026 (fe->get_nonzero_components(i)[component] ==
true);
2036template <
int n_components_,
int dim,
int spacedim,
typename Number>
2037template <
bool is_face,
bool is_linear>
2041 const unsigned int geometry_index =
2042 mapping_info->template compute_geometry_index_offset<is_face>(
2043 current_cell_index, current_face_number);
2045 cell_type = mapping_info->get_cell_type(geometry_index);
2047 const_cast<unsigned int &
>(n_q_points_scalar) =
2048 mapping_info->get_n_q_points_unvectorized(geometry_index);
2051 const_cast<unsigned int &
>(n_q_batches) =
2052 (n_q_points_scalar + n_lanes_internal - 1) / n_lanes_internal;
2054 const unsigned int n_q_points_before = n_q_points;
2056 const_cast<unsigned int &
>(n_q_points) =
2057 (stride == 1) ? n_q_batches : n_q_points_scalar;
2059 if (n_q_points != n_q_points_before)
2062 values.resize(n_q_points);
2064 gradients.resize(n_q_points);
2067 if (n_q_points == 0)
2071 const unsigned int unit_point_offset =
2072 mapping_info->compute_unit_point_index_offset(geometry_index);
2075 unit_point_faces_ptr =
2076 mapping_info->get_unit_point_faces(unit_point_offset);
2078 unit_point_ptr = mapping_info->get_unit_point(unit_point_offset);
2081 const unsigned int data_offset =
2082 mapping_info->compute_data_index_offset(geometry_index);
2083 const unsigned int compressed_data_offset =
2084 mapping_info->compute_compressed_data_index_offset(geometry_index);
2087 mapping_info->get_update_flags_mapping();
2089 real_point_ptr = mapping_info->get_real_point(data_offset);
2092 mapping_info->get_jacobian(compressed_data_offset, is_interior);
2094 inverse_jacobian_ptr =
2095 mapping_info->get_inverse_jacobian(compressed_data_offset, is_interior);
2097 normal_ptr = mapping_info->get_normal_vector(data_offset);
2099 JxW_ptr = mapping_info->get_JxW(data_offset);
2101 real_point_ptr = mapping_info->get_real_point(data_offset);
2103 mapping_info->get_jacobian(compressed_data_offset, is_interior);
2104 inverse_jacobian_ptr =
2105 mapping_info->get_inverse_jacobian(compressed_data_offset, is_interior);
2106 normal_ptr = mapping_info->get_normal_vector(data_offset);
2107 JxW_ptr = mapping_info->get_JxW(data_offset);
2110 if (!is_linear && fast_path)
2112 const std::size_t n_shapes = poly.size();
2114 shapes_faces.resize_fast(n_q_batches * n_shapes);
2116 shapes.resize_fast(n_q_batches * n_shapes);
2118 for (
unsigned int qb = 0; qb < n_q_batches; ++qb)
2124 shapes_faces.data() + qb * n_shapes,
2126 unit_point_faces_ptr[qb],
2139 else if (qb + 1 < n_q_batches)
2144 shapes.data() + qb * n_shapes,
2147 unit_point_ptr[qb + 1]);
2163template <
int n_components_,
int dim,
int spacedim,
typename Number>
2167 Number>::value_type &
2169 const unsigned int point_index)
const
2172 return values[point_index];
2177template <
int n_components_,
int dim,
int spacedim,
typename Number>
2181 Number>::gradient_type &
2183 const unsigned int point_index)
const
2186 return gradients[point_index];
2191template <
int n_components_,
int dim,
int spacedim,
typename Number>
2194 const unsigned int point_index)
const
2196 static_assert(n_components == dim,
2197 "Only makes sense for a vector field with dim components");
2200 return trace(gradients[point_index]);
2205template <
int n_components_,
int dim,
int spacedim,
typename Number>
2209 const unsigned int point_index)
2212 values[point_index] = value;
2217template <
int n_components_,
int dim,
int spacedim,
typename Number>
2221 const unsigned int point_index)
2224 gradients[point_index] = gradient;
2229template <
int n_components_,
int dim,
int spacedim,
typename Number>
2232 const Number &value,
2233 const unsigned int point_index)
2235 static_assert(n_components == dim,
2236 "Only makes sense for a vector field with dim components");
2240 for (
unsigned int d = 0; d < dim; ++d)
2241 gradients[point_index][d][d] = value;
2246template <
int n_components_,
int dim,
int spacedim,
typename Number>
2247Tensor<1, (dim == 2 ? 1 : dim), Number>
2249 const unsigned int point_index)
const
2252 dim > 1 && n_components == dim,
2253 "Only makes sense for a vector field with dim components and dim > 1");
2256 Tensor<1, (dim == 2 ? 1 : dim), Number> curl;
2260 curl[0] = grad[1][0] - grad[0][1];
2263 curl[0] = grad[2][1] - grad[1][2];
2264 curl[1] = grad[0][2] - grad[2][0];
2265 curl[2] = grad[1][0] - grad[0][1];
2275template <
int n_components_,
int dim,
int spacedim,
typename Number>
2278 const unsigned int point_index)
const
2281 Assert(jacobian_ptr !=
nullptr,
2283 ExcFEPointEvaluationAccessToUninitializedMappingField(
2284 "update_jacobians"));
2285 return jacobian_ptr[cell_type <= ::internal::MatrixFreeFunctions::
2286 GeometryType::affine ?
2293template <
int n_components_,
int dim,
int spacedim,
typename Number>
2296 const unsigned int point_index)
const
2299 Assert(inverse_jacobian_ptr !=
nullptr,
2301 ExcFEPointEvaluationAccessToUninitializedMappingField(
2302 "update_inverse_jacobians"));
2303 return inverse_jacobian_ptr
2312template <
int n_components_,
int dim,
int spacedim,
typename Number>
2315 const unsigned int point_index)
const
2318 Assert(JxW_ptr !=
nullptr,
2320 ExcFEPointEvaluationAccessToUninitializedMappingField(
2321 "update_JxW_values"));
2322 return JxW_ptr[point_index];
2327template <
int n_components_,
int dim,
int spacedim,
typename Number>
2330 const unsigned int point_index)
const
2332 return quadrature_point(point_index);
2338template <
int n_components_,
int dim,
int spacedim,
typename Number>
2341 const unsigned int point_index)
const
2344 Assert(real_point_ptr !=
nullptr,
2346 ExcFEPointEvaluationAccessToUninitializedMappingField(
2347 "update_quadrature_points"));
2348 return real_point_ptr[point_index];
2353template <
int n_components_,
int dim,
int spacedim,
typename Number>
2356 const unsigned int point_index)
const
2359 Assert(unit_point_ptr !=
nullptr,
ExcMessage(
"unit_point_ptr is not set!"));
2361 for (
unsigned int d = 0; d < dim; ++d)
2363 unit_point_ptr[point_index / stride][d], point_index % stride);
2369template <
int n_components_,
int dim,
int spacedim,
typename Number>
2374 return {0U, n_q_points};
2379template <
int n_components_,
int dim,
int spacedim,
typename Number>
2383 const unsigned int first_selected_component)
2387 first_selected_component)
2392template <
int n_components_,
int dim,
int spacedim,
typename Number>
2397 const unsigned int first_selected_component)
2402 first_selected_component)
2407template <
int n_components_,
int dim,
int spacedim,
typename Number>
2415 if (this->use_linear_path)
2416 this->
template do_reinit<false, true>();
2418 this->
template do_reinit<false, false>();
2423template <
int n_components_,
int dim,
int spacedim,
typename Number>
2427 internal_reinit_single_cell_state_mapping_info();
2428 this->must_reinitialize_pointers =
false;
2433template <
int n_components_,
int dim,
int spacedim,
typename Number>
2440 AssertThrow(this->mapping_info_on_the_fly.get() !=
nullptr,
2443 this->mapping_info_on_the_fly->reinit(cell, unit_points);
2444 this->must_reinitialize_pointers =
false;
2446 if (!this->fast_path)
2448 this->fe_values = std::make_shared<FEValues<dim, spacedim>>(
2452 std::vector<
Point<dim>>(unit_points.begin(), unit_points.end())),
2453 this->update_flags);
2454 this->fe_values->reinit(cell);
2457 if (this->use_linear_path)
2458 this->
template do_reinit<false, true>();
2460 this->
template do_reinit<false, false>();
2465template <
int n_components_,
int dim,
int spacedim,
typename Number>
2472 this->must_reinitialize_pointers =
false;
2474 if (this->use_linear_path)
2475 this->
template do_reinit<false, true>();
2477 this->
template do_reinit<false, false>();
2479 if (!this->fast_path)
2481 std::vector<Point<dim>> unit_points(this->n_q_points_scalar);
2483 for (
unsigned int v = 0; v < this->n_q_points_scalar; ++v)
2484 for (
unsigned int d = 0; d < dim; ++d)
2486 this->unit_point_ptr[v / n_lanes_internal][d][v % n_lanes_internal];
2488 this->fe_values = std::make_shared<FEValues<dim, spacedim>>(
2492 std::vector<
Point<dim>>(unit_points.begin(), unit_points.end())),
2493 this->update_flags);
2495 this->fe_values->reinit(
2496 this->mapping_info->get_cell_iterator(this->current_cell_index));
2502template <
int n_components_,
int dim,
int spacedim,
typename Number>
2503template <std::
size_t str
ide_view>
2515 if (this->must_reinitialize_pointers)
2516 internal_reinit_single_cell_state_mapping_info();
2518 if (this->n_q_points == 0)
2522 if (this->fast_path)
2524 if (this->use_linear_path)
2525 evaluate_fast<true>(solution_values, evaluation_flags);
2527 evaluate_fast<false>(solution_values, evaluation_flags);
2530 evaluate_slow(solution_values, evaluation_flags);
2535template <
int n_components_,
int dim,
int spacedim,
typename Number>
2542 solution_values.
size()),
2548template <
int n_components_,
int dim,
int spacedim,
typename Number>
2549template <std::
size_t str
ide_view>
2554 const bool sum_into_values)
2556 do_integrate<true>(solution_values, integration_flags, sum_into_values);
2561template <
int n_components_,
int dim,
int spacedim,
typename Number>
2566 const bool sum_into_values)
2569 solution_values.
size()),
2576template <
int n_components_,
int dim,
int spacedim,
typename Number>
2584 for (
const auto point_index : this->quadrature_point_indices())
2585 return_value += values[point_index] * this->JxW(point_index);
2587 return ETT::sum_value(return_value);
2592template <
int n_components_,
int dim,
int spacedim,
typename Number>
2599 "Calling this function only makes sense in fully vectorized mode."));
2600 if (q == n_q_batches - 1)
2602 const unsigned int n_filled_lanes =
2603 n_q_points_scalar & (n_lanes_user_interface - 1);
2604 if (n_filled_lanes == 0)
2605 return n_lanes_user_interface;
2607 return n_filled_lanes;
2610 return n_lanes_user_interface;
2615template <
int n_components_,
int dim,
int spacedim,
typename Number>
2616template <std::
size_t str
ide_view>
2621 const bool sum_into_values)
2623 do_integrate<false>(solution_values, integration_flags, sum_into_values);
2628template <
int n_components_,
int dim,
int spacedim,
typename Number>
2633 const bool sum_into_values)
2636 solution_values.
size()),
2643template <
int n_components_,
int dim,
int spacedim,
typename Number>
2644template <
bool is_linear, std::
size_t str
ide_view>
2649 const unsigned int dofs_per_comp =
2652 for (
unsigned int comp = 0; comp < n_components; ++comp)
2654 const std::size_t offset =
2655 (this->component_in_base_element + comp) * dofs_per_comp;
2657 if ((is_linear && n_components == 1) || this->renumber.empty())
2659 for (
unsigned int i = 0; i < dofs_per_comp; ++i)
2660 ETT::read_value(solution_values[i + offset],
2662 this->solution_renumbered[i]);
2666 const unsigned int *renumber_ptr = this->renumber.data() + offset;
2667 for (
unsigned int i = 0; i < dofs_per_comp; ++i)
2668 ETT::read_value(solution_values[renumber_ptr[i]],
2670 this->solution_renumbered[i]);
2677template <
int n_components_,
int dim,
int spacedim,
typename Number>
2678template <
bool is_linear, std::
size_t str
ide_view>
2683 const unsigned int n_shapes,
2684 const unsigned int qb,
2690 std::array<vectorized_value_type, dim + 1> result;
2691 if constexpr (is_linear)
2693 if constexpr (n_components == 1)
2700 stride_view>(solution_values.
data(), this->unit_point_ptr[qb]);
2704 this->solution_renumbered.data(), this->unit_point_ptr[qb]);
2712 false>(this->shapes.data() + qb * n_shapes,
2714 this->solution_renumbered.data());
2715 gradient[0] = result[0];
2717 gradient[1] = result[1];
2719 gradient[2] = result[2];
2720 value = result[dim];
2724 if constexpr (is_linear)
2726 if constexpr (n_components == 1)
2731 stride_view>(solution_values.
data(), this->unit_point_ptr[qb]);
2734 this->solution_renumbered.data(), this->unit_point_ptr[qb]);
2742 this->shapes.data() + qb * n_shapes,
2744 this->solution_renumbered.data());
2750template <
int n_components_,
int dim,
int spacedim,
typename Number>
2751template <
bool is_linear, std::
size_t str
ide_view>
2757 if (!(is_linear && n_components == 1))
2758 prepare_evaluate_fast<is_linear>(solution_values);
2761 const unsigned int n_shapes = is_linear ? 2 : this->poly.size();
2763 for (
unsigned int qb = 0; qb < this->n_q_batches; ++qb)
2768 compute_evaluate_fast<is_linear>(
2769 solution_values, evaluation_flags, n_shapes, qb, value, gradient);
2773 for (
unsigned int v = 0, offset = qb * stride;
2774 v < stride && (stride == 1 || offset < this->n_q_points_scalar);
2776 ETT::set_value(value, v, this->values[offset]);
2784 for (
unsigned int v = 0, offset = qb * stride;
2785 v < stride && (stride == 1 || offset < this->n_q_points_scalar);
2789 ETT::set_gradient(gradient, v, unit_gradient);
2790 this->gradients[offset] =
2794 this->inverse_jacobian_ptr[0].
transpose(), unit_gradient) :
2797 ->inverse_jacobian_ptr[this->cell_type <=
2799 GeometryType::affine ?
2811template <
int n_components_,
int dim,
int spacedim,
typename Number>
2812template <std::
size_t str
ide_view>
2819 Assert(this->fe_values.get() !=
nullptr,
2821 "Not initialized. Please call FEPointEvaluation::reinit()!"));
2823 const std::size_t n_points = this->fe_values->get_quadrature().size();
2827 this->values.resize(this->n_q_points);
2828 std::fill(this->values.begin(), this->values.end(),
value_type());
2829 for (
unsigned int i = 0; i < this->fe->n_dofs_per_cell(); ++i)
2832 for (
unsigned int d = 0; d < n_components; ++d)
2833 if (this->nonzero_shape_function_component[i][d] &&
2834 (this->fe->is_primitive(i) || this->fe->is_primitive()))
2835 for (
unsigned int qb = 0, q = 0; q < n_points;
2836 ++qb, q += n_lanes_user_interface)
2837 for (
unsigned int v = 0;
2838 v < n_lanes_user_interface && q + v < n_points;
2840 ETT::access(this->values[qb],
2843 this->fe_values->shape_value(i, q + v) * value);
2844 else if (this->nonzero_shape_function_component[i][d])
2845 for (
unsigned int qb = 0, q = 0; q < n_points;
2846 ++qb, q += n_lanes_user_interface)
2847 for (
unsigned int v = 0;
2848 v < n_lanes_user_interface && q + v < n_points;
2850 ETT::access(this->values[qb],
2853 this->fe_values->shape_value_component(i,
2862 this->gradients.resize(this->n_q_points);
2863 std::fill(this->gradients.begin(),
2864 this->gradients.end(),
2866 for (
unsigned int i = 0; i < this->fe->n_dofs_per_cell(); ++i)
2869 for (
unsigned int d = 0; d < n_components; ++d)
2870 if (this->nonzero_shape_function_component[i][d] &&
2871 (this->fe->is_primitive(i) || this->fe->is_primitive()))
2872 for (
unsigned int qb = 0, q = 0; q < n_points;
2873 ++qb, q += n_lanes_user_interface)
2874 for (
unsigned int v = 0;
2875 v < n_lanes_user_interface && q + v < n_points;
2877 ETT::access(this->gradients[qb],
2880 this->fe_values->shape_grad(i, q + v) * value);
2881 else if (this->nonzero_shape_function_component[i][d])
2882 for (
unsigned int qb = 0, q = 0; q < n_points;
2883 ++qb, q += n_lanes_user_interface)
2884 for (
unsigned int v = 0;
2885 v < n_lanes_user_interface && q + v < n_points;
2888 this->gradients[qb],
2891 this->fe_values->shape_grad_component(i, q + v, d) * value);
2898template <
int n_components_,
int dim,
int spacedim,
typename Number>
2899template <
bool is_linear>
2903 const unsigned int n_shapes,
2904 const unsigned int qb,
2919 solution_values_vectorized_linear :
2920 this->solution_renumbered_vectorized.data(),
2921 this->unit_point_ptr[qb],
2928 this->shapes.data() + qb * n_shapes,
2931 is_linear ? solution_values_vectorized_linear :
2932 this->solution_renumbered_vectorized.data(),
2933 this->unit_point_ptr[qb],
2939template <
int n_components_,
int dim,
int spacedim,
typename Number>
2940template <
bool is_linear, std::
size_t str
ide_view>
2945 const bool sum_into_values)
2947 if (!sum_into_values && this->fe->n_components() > n_components)
2948 for (
unsigned int i = 0; i < solution_values.
size(); ++i)
2949 solution_values[i] = 0;
2951 const unsigned int dofs_per_comp =
2954 for (
unsigned int comp = 0; comp < n_components; ++comp)
2956 const std::size_t offset =
2957 (this->component_in_base_element + comp) * dofs_per_comp;
2959 if (is_linear || this->renumber.empty())
2961 for (
unsigned int i = 0; i < dofs_per_comp; ++i)
2962 if (sum_into_values)
2963 solution_values[i + offset] +=
2964 ETT::sum_value(comp,
2966 *(solution_values_vectorized_linear + i) :
2967 this->solution_renumbered_vectorized[i]);
2969 solution_values[i + offset] =
2970 ETT::sum_value(comp,
2972 *(solution_values_vectorized_linear + i) :
2973 this->solution_renumbered_vectorized[i]);
2977 const unsigned int *renumber_ptr = this->renumber.
data() + offset;
2978 for (
unsigned int i = 0; i < dofs_per_comp; ++i)
2979 if (sum_into_values)
2980 solution_values[renumber_ptr[i]] +=
2981 ETT::sum_value(comp, this->solution_renumbered_vectorized[i]);
2983 solution_values[renumber_ptr[i]] =
2984 ETT::sum_value(comp, this->solution_renumbered_vectorized[i]);
2991template <
int n_components_,
int dim,
int spacedim,
typename Number>
2992template <
bool do_JxW,
bool is_linear, std::
size_t str
ide_view>
2997 const bool sum_into_values)
3000 if constexpr (stride == 1)
3001 if (
const unsigned int n_filled_lanes =
3002 this->n_q_points_scalar & (n_lanes_internal - 1);
3006 for (
unsigned int v = n_filled_lanes; v < n_lanes_internal; ++v)
3007 ETT::set_zero_value(this->values.back(), v);
3009 for (
unsigned int v = n_filled_lanes; v < n_lanes_internal; ++v)
3010 ETT::set_zero_gradient(this->gradients.back(), v);
3014 solution_values_vectorized_linear = {};
3017 const unsigned int n_shapes = is_linear ? 2 : this->poly.size();
3019 const bool cartesian_cell =
3021 const bool affine_cell =
3023 for (
unsigned int qb = 0; qb < this->n_q_batches; ++qb)
3029 for (
unsigned int v = 0, offset = qb * stride;
3030 v < stride && (stride == 1 || offset < this->n_q_points_scalar);
3032 ETT::get_value(value,
3034 do_JxW ? this->values[offset] * this->JxW_ptr[offset] :
3035 this->values[offset]);
3038 for (
unsigned int v = 0, offset = qb * stride;
3039 v < stride && (stride == 1 || offset < this->n_q_points_scalar);
3043 do_JxW ? this->gradients[offset] * this->JxW_ptr[offset] :
3044 this->gradients[offset];
3052 this->inverse_jacobian_ptr[affine_cell ? 0 : offset],
3056 compute_integrate_fast<is_linear>(
3062 solution_values_vectorized_linear.data());
3066 finish_integrate_fast<is_linear>(solution_values,
3067 solution_values_vectorized_linear.
data(),
3073template <
int n_components_,
int dim,
int spacedim,
typename Number>
3074template <
bool do_JxW, std::
size_t str
ide_view>
3079 const bool sum_into_values)
3082 Assert(this->fe_values.get() !=
nullptr,
3084 "Not initialized. Please call FEPointEvaluation::reinit()!"));
3085 if (!sum_into_values)
3086 for (
unsigned int i = 0; i < solution_values.
size(); ++i)
3087 solution_values[i] = 0;
3089 const std::size_t n_points = this->fe_values->get_quadrature().
size();
3094 for (
unsigned int i = 0; i < this->fe->n_dofs_per_cell(); ++i)
3096 for (
unsigned int d = 0; d < n_components; ++d)
3097 if (this->nonzero_shape_function_component[i][d] &&
3098 (this->fe->is_primitive(i) || this->fe->is_primitive()))
3099 for (
unsigned int qb = 0, q = 0; q < n_points;
3100 ++qb, q += n_lanes_user_interface)
3101 for (
unsigned int v = 0;
3102 v < n_lanes_user_interface && q + v < n_points;
3104 solution_values[i] +=
3105 this->fe_values->shape_value(i, q + v) *
3106 ETT::access(this->values[qb], v, d) *
3107 (do_JxW ? this->fe_values->JxW(q + v) : 1.);
3108 else if (this->nonzero_shape_function_component[i][d])
3109 for (
unsigned int qb = 0, q = 0; q < n_points;
3110 ++qb, q += n_lanes_user_interface)
3111 for (
unsigned int v = 0;
3112 v < n_lanes_user_interface && q + v < n_points;
3114 solution_values[i] +=
3115 this->fe_values->shape_value_component(i, q + v, d) *
3116 ETT::access(this->values[qb], v, d) *
3117 (do_JxW ? this->fe_values->JxW(q + v) : 1.);
3124 for (
unsigned int i = 0; i < this->fe->n_dofs_per_cell(); ++i)
3126 for (
unsigned int d = 0; d < n_components; ++d)
3127 if (this->nonzero_shape_function_component[i][d] &&
3128 (this->fe->is_primitive(i) || this->fe->is_primitive()))
3129 for (
unsigned int qb = 0, q = 0; q < n_points;
3130 ++qb, q += n_lanes_user_interface)
3131 for (
unsigned int v = 0;
3132 v < n_lanes_user_interface && q + v < n_points;
3134 solution_values[i] +=
3135 this->fe_values->shape_grad(i, q + v) *
3136 ETT::access(this->gradients[qb], v, d) *
3137 (do_JxW ? this->fe_values->JxW(q + v) : 1.);
3138 else if (this->nonzero_shape_function_component[i][d])
3139 for (
unsigned int qb = 0, q = 0; q < n_points;
3140 ++qb, q += n_lanes_user_interface)
3141 for (
unsigned int v = 0;
3142 v < n_lanes_user_interface && q + v < n_points;
3144 solution_values[i] +=
3145 this->fe_values->shape_grad_component(i, q + v, d) *
3146 ETT::access(this->gradients[qb], v, d) *
3147 (do_JxW ? this->fe_values->JxW(q + v) : 1.);
3154template <
int n_components_,
int dim,
int spacedim,
typename Number>
3155template <
bool do_JxW, std::
size_t str
ide_view>
3160 const bool sum_into_values)
3162 if (this->must_reinitialize_pointers)
3163 internal_reinit_single_cell_state_mapping_info();
3167 if (this->n_q_points == 0 ||
3169 (integration_flags &
3172 if (!sum_into_values)
3173 for (
unsigned int i = 0; i < solution_values.
size(); ++i)
3174 solution_values[i] = 0;
3179 !do_JxW || this->JxW_ptr !=
nullptr,
3181 "JxW pointer is not set! If you do not want to integrate() use test_and_sum()"));
3184 if (this->fast_path)
3186 if (this->use_linear_path)
3187 integrate_fast<do_JxW, true>(solution_values,
3191 integrate_fast<do_JxW, false>(solution_values,
3196 integrate_slow<do_JxW>(solution_values, integration_flags, sum_into_values);
3201template <
int n_components_,
int dim,
int spacedim,
typename Number>
3204 const unsigned int point_index)
const
3207 Assert(this->normal_ptr !=
nullptr,
3209 ExcFEPointEvaluationAccessToUninitializedMappingField(
3210 "update_normal_vectors"));
3211 return this->normal_ptr[point_index];
3216template <
int n_components_,
int dim,
int spacedim,
typename Number>
3222 const unsigned int point_index)
const
3225 return this->gradients[point_index] * normal_vector(point_index);
3230template <
int n_components_,
int dim,
int spacedim,
typename Number>
3234 const unsigned int point_index)
3237 if constexpr (n_components == 1)
3238 this->gradients[point_index] = value * normal_vector(point_index);
3240 this->gradients[point_index] =
3246template <
int n_components_,
int dim,
int spacedim,
typename Number>
3251 const bool is_interior,
3252 const unsigned int first_selected_component)
3256 first_selected_component,
3262template <
int n_components_,
int dim,
int spacedim,
typename Number>
3266 const unsigned int face_number)
3269 this->current_face_number = face_number;
3270 this->must_reinitialize_pointers =
false;
3272 if (this->use_linear_path)
3273 this->
template do_reinit<true, true>();
3275 this->
template do_reinit<true, false>();
3280template <
int n_components_,
int dim,
int spacedim,
typename Number>
3283 const unsigned int face_index)
3285 this->current_cell_index = face_index;
3286 this->current_face_number =
3287 this->mapping_info->get_face_number(face_index, this->is_interior);
3288 this->must_reinitialize_pointers =
false;
3290 if (this->use_linear_path)
3291 this->
template do_reinit<true, true>();
3293 this->
template do_reinit<true, false>();
3298template <
int n_components_,
int dim,
int spacedim,
typename Number>
3299template <std::
size_t str
ide_view>
3305 Assert(!this->must_reinitialize_pointers,
3306 ExcMessage(
"Object has not been reinitialized!"));
3308 if (this->n_q_points == 0)
3319 if (this->use_linear_path)
3320 do_evaluate<true>(solution_values, evaluation_flags);
3322 do_evaluate<false>(solution_values, evaluation_flags);
3327template <
int n_components_,
int dim,
int spacedim,
typename Number>
3334 solution_values.
size()),
3340template <
int n_components_,
int dim,
int spacedim,
typename Number>
3341template <
bool is_linear, std::
size_t str
ide_view>
3347 const unsigned int dofs_per_comp =
3351 if (stride_view == 1 && this->component_in_base_element == 0 &&
3352 (is_linear || this->renumber.empty()))
3353 input = solution_values.
data();
3356 for (
unsigned int comp = 0; comp < n_components; ++comp)
3358 const std::size_t offset =
3359 (this->component_in_base_element + comp) * dofs_per_comp;
3361 if (is_linear || this->renumber.empty())
3363 for (
unsigned int i = 0; i < dofs_per_comp; ++i)
3364 this->scratch_data_scalar[i + comp * dofs_per_comp] =
3365 solution_values[i + offset];
3369 const unsigned int *renumber_ptr = this->renumber.
data() + offset;
3370 for (
unsigned int i = 0; i < dofs_per_comp; ++i)
3371 this->scratch_data_scalar[i + comp * dofs_per_comp] =
3372 solution_values[renumber_ptr[i]];
3375 input = this->scratch_data_scalar.
data();
3379 this->scratch_data_scalar.begin() + dofs_per_comp * n_components;
3382 template interpolate<true, false>(n_components,
3387 this->current_face_number);
3389 do_evaluate_in_face<is_linear, 1>(output, evaluation_flags);
3394template <
int n_components_,
int dim,
int spacedim,
typename Number>
3395template <std::
size_t str
ide_view>
3400 const bool sum_into_values)
3402 Assert(!this->must_reinitialize_pointers,
3403 ExcMessage(
"Object has not been reinitialized!"));
3407 if (this->n_q_points == 0 ||
3409 (integration_flags &
3412 if (!sum_into_values)
3413 for (
unsigned int i = 0; i < solution_values.
size(); ++i)
3414 solution_values[i] = 0;
3420 if (this->use_linear_path)
3421 do_integrate<true, true>(solution_values,
3425 do_integrate<true, false>(solution_values,
3432template <
int n_components_,
int dim,
int spacedim,
typename Number>
3437 const bool sum_into_values)
3440 solution_values.
size()),
3447template <
int n_components_,
int dim,
int spacedim,
typename Number>
3448template <std::
size_t str
ide_view>
3453 const bool sum_into_values)
3455 Assert(!this->must_reinitialize_pointers,
3456 ExcMessage(
"Object has not been reinitialized!"));
3460 if (this->n_q_points == 0 ||
3462 (integration_flags &
3465 if (!sum_into_values)
3466 for (
unsigned int i = 0; i < solution_values.
size(); ++i)
3467 solution_values[i] = 0;
3473 if (this->use_linear_path)
3474 do_integrate<false, true>(solution_values,
3478 do_integrate<false, false>(solution_values,
3485template <
int n_components_,
int dim,
int spacedim,
typename Number>
3490 const bool sum_into_values)
3493 solution_values.
size()),
3500template <
int n_components_,
int dim,
int spacedim,
typename Number>
3501template <
bool do_JxW,
bool is_linear, std::
size_t str
ide_view>
3506 const bool sum_into_values)
3508 if (!sum_into_values && this->fe->n_components() > n_components)
3509 for (
unsigned int i = 0; i < solution_values.
size(); ++i)
3510 solution_values[i] = 0;
3512 do_integrate_in_face<do_JxW, is_linear, 1>(this->scratch_data_scalar.begin(),
3516 ScalarNumber *input = this->scratch_data_scalar.begin();
3518 if (stride_view == 1 && this->component_in_base_element == 0 &&
3519 (is_linear || this->renumber.empty()))
3521 if (sum_into_values)
3523 FEFaceNormalEvaluationImpl<dim, is_linear ? 1 : -1,
ScalarNumber>::
3524 template interpolate<false, true>(n_components,
3528 solution_values.
data(),
3529 this->current_face_number);
3532 FEFaceNormalEvaluationImpl<dim, is_linear ? 1 : -1,
ScalarNumber>::
3533 template interpolate<false, false>(n_components,
3537 solution_values.
data(),
3538 this->current_face_number);
3542 const unsigned int dofs_per_comp_face =
3543 is_linear ?
Utilities::pow(2, dim - 1) : this->dofs_per_component_face;
3545 const unsigned int size_input = 3 * dofs_per_comp_face * n_components;
3549 FEFaceNormalEvaluationImpl<dim, is_linear ? 1 : -1,
ScalarNumber>::
3550 template interpolate<false, false>(n_components,
3555 this->current_face_number);
3557 const unsigned int dofs_per_comp =
3560 for (
unsigned int comp = 0; comp < n_components; ++comp)
3562 const std::size_t offset =
3563 (this->component_in_base_element + comp) * dofs_per_comp;
3565 if (is_linear || this->renumber.empty())
3567 for (
unsigned int i = 0; i < dofs_per_comp; ++i)
3568 if (sum_into_values)
3569 solution_values[i + offset] +=
3570 output[i + comp * dofs_per_comp];
3572 solution_values[i + offset] =
3573 output[i + comp * dofs_per_comp];
3577 const unsigned int *renumber_ptr = this->renumber.
data() + offset;
3578 for (
unsigned int i = 0; i < dofs_per_comp; ++i)
3579 if (sum_into_values)
3580 solution_values[renumber_ptr[i]] +=
3581 output[i + comp * dofs_per_comp];
3583 solution_values[renumber_ptr[i]] =
3584 output[i + comp * dofs_per_comp];
3592template <
int n_components_,
int dim,
int spacedim,
typename Number>
3593template <
int str
ide_face_dof>
3599 if (this->use_linear_path)
3600 do_evaluate_in_face<true, stride_face_dof>(face_dof_values,
3603 do_evaluate_in_face<false, stride_face_dof>(face_dof_values,
3609template <
int n_components_,
int dim,
int spacedim,
typename Number>
3610template <
bool is_linear,
int str
ide_face_dof>
3617 if constexpr (n_components == 1)
3618 face_dof_values_ptr = face_dof_values;
3621 const unsigned int dofs_per_comp_face =
3622 is_linear ?
Utilities::pow(2, dim - 1) : this->dofs_per_component_face;
3623 for (
unsigned int comp = 0; comp < n_components; ++comp)
3624 for (
unsigned int i = 0; i < 2 * dofs_per_comp_face; ++i)
3625 ETT::read_value(face_dof_values[(i + comp * 3 * dofs_per_comp_face) *
3628 this->solution_renumbered[i]);
3630 face_dof_values_ptr = this->solution_renumbered.data();
3633 constexpr int stride_face_dof_actual =
3634 n_components == 1 ? stride_face_dof : 1;
3637 const unsigned int n_shapes = is_linear ? 2 : this->poly.size();
3639 for (
unsigned int qb = 0; qb < this->n_q_batches; ++qb)
3646 const std::array<vectorized_value_type, dim + 1> interpolated_value =
3653 stride_face_dof_actual>(face_dof_values_ptr,
3654 this->unit_point_faces_ptr[qb]) :
3661 stride_face_dof_actual>(this->shapes_faces.data() +
3664 face_dof_values_ptr);
3666 value = interpolated_value[dim - 1];
3669 if (this->current_face_number / 2 == 0)
3671 gradient[0] = interpolated_value[dim];
3673 gradient[1] = interpolated_value[0];
3675 gradient[2] = interpolated_value[1];
3677 else if (this->current_face_number / 2 == 1)
3680 gradient[1] = interpolated_value[dim];
3683 gradient[0] = interpolated_value[1];
3684 gradient[2] = interpolated_value[0];
3687 gradient[0] = interpolated_value[0];
3691 else if (this->current_face_number / 2 == 2)
3695 gradient[0] = interpolated_value[0];
3696 gradient[1] = interpolated_value[1];
3697 gradient[2] = interpolated_value[dim];
3712 stride_face_dof_actual>(face_dof_values_ptr,
3713 this->unit_point_faces_ptr[qb]) :
3719 stride_face_dof_actual>(this->shapes_faces.data() +
3722 face_dof_values_ptr);
3727 for (
unsigned int v = 0, offset = qb * stride;
3728 v < stride && (stride == 1 || offset < this->n_q_points_scalar);
3730 ETT::set_value(value, v, this->values[offset]);
3738 for (
unsigned int v = 0, offset = qb * stride;
3739 v < stride && (stride == 1 || offset < this->n_q_points_scalar);
3743 ETT::set_gradient(gradient, v, unit_gradient);
3744 this->gradients[offset] =
3748 this->inverse_jacobian_ptr[0].
transpose(), unit_gradient) :
3751 ->inverse_jacobian_ptr[this->cell_type <=
3753 GeometryType::affine ?
3765template <
int n_components_,
int dim,
int spacedim,
typename Number>
3766template <
int str
ide_face_dof>
3771 const bool sum_into_values)
3773 if (this->use_linear_path)
3774 do_integrate_in_face<true, true, stride_face_dof>(face_dof_values,
3778 do_integrate_in_face<true, false, stride_face_dof>(face_dof_values,
3785template <
int n_components_,
int dim,
int spacedim,
typename Number>
3786template <
bool do_JxW,
bool is_linear,
int str
ide_face_dof>
3792 const bool sum_into_values)
3795 if constexpr (stride == 1)
3796 if (
const unsigned int n_filled_lanes =
3797 this->n_q_points_scalar & (n_lanes_internal - 1);
3801 for (
unsigned int v = n_filled_lanes; v < n_lanes_internal; ++v)
3802 ETT::set_zero_value(this->values.back(), v);
3804 for (
unsigned int v = n_filled_lanes; v < n_lanes_internal; ++v)
3805 ETT::set_zero_gradient(this->gradients.back(), v);
3810 solution_values_vectorized_linear = {};
3813 const unsigned int n_shapes = is_linear ? 2 : this->poly.size();
3815 const bool cartesian_cell =
3817 const bool affine_cell =
3819 for (
unsigned int qb = 0; qb < this->n_q_batches; ++qb)
3825 for (
unsigned int v = 0, offset = qb * stride;
3826 v < stride && (stride == 1 || offset < this->n_q_points_scalar);
3828 ETT::get_value(value,
3830 do_JxW ? this->values[offset] * this->JxW_ptr[offset] :
3831 this->values[offset]);
3834 for (
unsigned int v = 0, offset = qb * stride;
3835 v < stride && (stride == 1 || offset < this->n_q_points_scalar);
3839 do_JxW ? this->gradients[offset] * this->JxW_ptr[offset] :
3840 this->gradients[offset];
3848 this->inverse_jacobian_ptr[affine_cell ? 0 : offset],
3854 std::array<vectorized_value_type, 2> value_face = {};
3857 value_face[0] = value;
3860 if (this->current_face_number / 2 == 0)
3862 value_face[1] = gradient[0];
3864 gradient_in_face[0] = gradient[1];
3866 gradient_in_face[1] = gradient[2];
3868 else if (this->current_face_number / 2 == 1)
3871 value_face[1] = gradient[1];
3874 gradient_in_face[0] = gradient[2];
3875 gradient_in_face[1] = gradient[0];
3878 gradient_in_face[0] = gradient[0];
3882 else if (this->current_face_number / 2 == 2)
3886 value_face[1] = gradient[2];
3887 gradient_in_face[0] = gradient[0];
3888 gradient_in_face[1] = gradient[1];
3901 2>(this->shapes_faces.data() + qb * n_shapes,
3905 is_linear ? solution_values_vectorized_linear.data() :
3906 this->solution_renumbered_vectorized.data(),
3907 this->unit_point_faces_ptr[qb],
3915 this->shapes_faces.data() + qb * n_shapes,
3918 is_linear ? solution_values_vectorized_linear.data() :
3919 this->solution_renumbered_vectorized.data(),
3920 this->unit_point_faces_ptr[qb],
3924 const unsigned int dofs_per_comp_face =
3925 is_linear ?
Utilities::pow(2, dim - 1) : this->dofs_per_component_face;
3927 for (
unsigned int comp = 0; comp < n_components; ++comp)
3928 for (
unsigned int i = 0; i < 2 * dofs_per_comp_face; ++i)
3929 if (sum_into_values)
3930 face_dof_values[(i + comp * 3 * dofs_per_comp_face) *
3932 ETT::sum_value(comp,
3934 *(solution_values_vectorized_linear.data() + i) :
3935 this->solution_renumbered_vectorized[i]);
3937 face_dof_values[(i + comp * 3 * dofs_per_comp_face) * stride_face_dof] =
3938 ETT::sum_value(comp,
3940 *(solution_values_vectorized_linear.data() + i) :
3941 this->solution_renumbered_vectorized[i]);
3946template <
int n_components_,
int dim,
int spacedim,
typename Number>
3949 const unsigned int point_index)
const
3952 Assert(this->normal_ptr !=
nullptr,
3954 ExcFEPointEvaluationAccessToUninitializedMappingField(
3955 "update_normal_vectors"));
3959 for (
unsigned int d = 0; d < dim; ++d)
3968 return this->normal_ptr[point_index];
3974template <
int n_components_,
int dim,
int spacedim,
typename Number>
3983 return this->gradients[point_index] * normal_vector(point_index);
3988template <
int n_components_,
int dim,
int spacedim,
typename Number>
3992 const unsigned int point_index)
3995 if constexpr (n_components == 1)
3996 this->gradients[point_index] = value * normal_vector(point_index);
3998 this->gradients[point_index] =
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
#define DEAL_II_NAMESPACE_CLOSE
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.
#define DEAL_II_ASSERT_UNREACHABLE()
#define DEAL_II_NOT_IMPLEMENTED()
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)
static const 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 SymmetricTensor< 4, dim, Number > outer_product(const SymmetricTensor< 2, dim, Number > &t1, const SymmetricTensor< 2, dim, Number > &t2)
DEAL_II_HOST constexpr Number trace(const SymmetricTensor< 2, dim2, Number > &)