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,
776 JxW(
const unsigned int point_index)
const;
842 setup(
const unsigned int first_selected_component);
849 template <
bool is_face,
bool is_linear>
882 std::vector<Polynomials::Polynomial<double>>
poly;
1018 std::unique_ptr<NonMatching::MappingInfo<dim, spacedim, Number>>
1100template <
int n_components_,
1103 typename Number =
double>
1117 using ETT =
typename internal::FEPointEvaluation::
1118 EvaluatorTypeTraits<dim, spacedim, n_components, Number>;
1125 typename ETT::interface_vectorized_unit_gradient_type;
1148 const unsigned int first_selected_component = 0);
1169 const unsigned int first_selected_component = 0);
1212 template <std::
size_t str
ide_view>
1255 template <std::
size_t str
ide_view>
1259 const bool sum_into_values =
false);
1286 const bool sum_into_values =
false);
1314 template <std::
size_t str
ide_view>
1319 const bool sum_into_values =
false);
1350 const bool sum_into_values =
false);
1372 template <
bool is_linear, std::
size_t str
ide_view>
1381 template <
bool is_linear, std::
size_t str
ide_view>
1386 const unsigned int n_shapes,
1387 const unsigned int qb,
1394 template <
bool is_linear, std::
size_t str
ide_view>
1403 template <std::
size_t str
ide_view>
1414 template <
bool is_linear>
1418 const unsigned int n_shapes,
1419 const unsigned int qb,
1429 template <
bool is_linear, std::
size_t str
ide_view>
1434 const bool sum_into_values);
1439 template <
bool do_JxW,
bool is_linear, std::
size_t str
ide_view>
1444 const bool sum_into_values);
1449 template <
bool do_JxW, std::
size_t str
ide_view>
1454 const bool sum_into_values);
1459 template <
bool do_JxW, std::
size_t str
ide_view>
1464 const bool sum_into_values);
1486template <
int n_components_,
1489 typename Number =
double>
1503 using ETT =
typename internal::FEPointEvaluation::
1504 EvaluatorTypeTraits<dim, spacedim, n_components, Number>;
1511 typename ETT::interface_vectorized_unit_gradient_type;
1520 const unsigned int first_selected_component = 0);
1534 reinit(
const unsigned int face_index);
1547 template <std::
size_t str
ide_view>
1590 template <std::
size_t str
ide_view>
1594 const bool sum_into_values =
false);
1621 const bool sum_into_values =
false);
1645 template <std::
size_t str
ide_view>
1650 const bool sum_into_values =
false);
1677 const bool sum_into_values =
false);
1685 template <
int str
ide_face_dof = VectorizedArrayType::size()>
1696 template <
int str
ide_face_dof = VectorizedArrayType::size()>
1700 const bool sum_into_values =
false);
1718 template <
bool is_linear, std::
size_t str
ide_view>
1724 template <
bool do_JxW,
bool is_linear, std::
size_t str
ide_view>
1729 const bool sum_into_values);
1735 template <
bool is_linear,
int str
ide_face_dof>
1744 template <
bool do_JxW,
bool is_linear,
int str
ide_face_dof>
1749 const bool sum_into_values);
1757template <
int n_components_,
int dim,
int spacedim,
typename Number>
1762 const unsigned int first_selected_component)
1763 : n_q_batches(
numbers::invalid_unsigned_int)
1764 , n_q_points(
numbers::invalid_unsigned_int)
1765 , n_q_points_scalar(
numbers::invalid_unsigned_int)
1769 , update_flags(update_flags)
1770 , mapping_info_on_the_fly(
1771 std::make_unique<
NonMatching::MappingInfo<dim, spacedim, Number>>(
1774 , mapping_info(mapping_info_on_the_fly.get())
1775 , current_cell_index(
numbers::invalid_unsigned_int)
1776 , current_face_number(
numbers::invalid_unsigned_int)
1777 , is_reinitialized(false)
1780 setup(first_selected_component);
1785template <
int n_components_,
int dim,
int spacedim,
typename Number>
1790 const unsigned int first_selected_component,
1791 const bool is_interior)
1792 : n_q_batches(
numbers::invalid_unsigned_int)
1793 , n_q_points(
numbers::invalid_unsigned_int)
1794 , n_q_points_scalar(
numbers::invalid_unsigned_int)
1795 , mapping(&mapping_info.get_mapping())
1798 , update_flags(mapping_info.get_update_flags())
1799 , mapping_info(&mapping_info)
1800 , current_cell_index(
numbers::invalid_unsigned_int)
1801 , current_face_number(
numbers::invalid_unsigned_int)
1802 , is_reinitialized(false)
1803 , is_interior(is_interior)
1805 setup(first_selected_component);
1812template <
int n_components_,
int dim,
int spacedim,
typename Number>
1816 : n_q_batches(other.n_q_batches)
1817 , n_q_points(other.n_q_points)
1818 , n_q_points_scalar(other.n_q_points_scalar)
1819 , mapping(other.mapping)
1822 , use_linear_path(other.use_linear_path)
1823 , renumber(other.renumber)
1824 , solution_renumbered(other.solution_renumbered)
1825 , solution_renumbered_vectorized(other.solution_renumbered_vectorized)
1826 , values(other.values)
1827 , gradients(other.gradients)
1828 , dofs_per_component(other.dofs_per_component)
1829 , dofs_per_component_face(other.dofs_per_component_face)
1830 , component_in_base_element(other.component_in_base_element)
1831 , nonzero_shape_function_component(other.nonzero_shape_function_component)
1832 , update_flags(other.update_flags)
1833 , fe_values(other.fe_values)
1834 , mapping_info_on_the_fly(
1835 other.mapping_info_on_the_fly ?
1840 , mapping_info(other.mapping_info)
1841 , current_cell_index(other.current_cell_index)
1842 , current_face_number(other.current_face_number)
1843 , fast_path(other.fast_path)
1844 , is_reinitialized(false)
1845 , is_interior(other.is_interior)
1847 connection_is_reinitialized = mapping_info->connect_is_reinitialized(
1848 [
this]() { this->is_reinitialized =
false; });
1853template <
int n_components_,
int dim,
int spacedim,
typename Number>
1856 : n_q_batches(other.n_q_batches)
1857 , n_q_points(other.n_q_points)
1858 , n_q_points_scalar(other.n_q_points_scalar)
1859 , mapping(other.mapping)
1862 , use_linear_path(other.use_linear_path)
1863 , renumber(other.renumber)
1864 , solution_renumbered(other.solution_renumbered)
1865 , solution_renumbered_vectorized(other.solution_renumbered_vectorized)
1866 , values(other.values)
1867 , gradients(other.gradients)
1868 , dofs_per_component(other.dofs_per_component)
1869 , dofs_per_component_face(other.dofs_per_component_face)
1870 , component_in_base_element(other.component_in_base_element)
1871 , nonzero_shape_function_component(other.nonzero_shape_function_component)
1872 , update_flags(other.update_flags)
1873 , fe_values(other.fe_values)
1874 , mapping_info_on_the_fly(std::move(other.mapping_info_on_the_fly))
1875 , mapping_info(other.mapping_info)
1876 , current_cell_index(other.current_cell_index)
1877 , current_face_number(other.current_face_number)
1878 , fast_path(other.fast_path)
1879 , is_reinitialized(
false)
1880 , is_interior(other.is_interior)
1882 connection_is_reinitialized = mapping_info->connect_is_reinitialized(
1883 [
this]() { this->is_reinitialized =
false; });
1888template <
int n_components_,
int dim,
int spacedim,
typename Number>
1892 connection_is_reinitialized.disconnect();
1897template <
int n_components_,
int dim,
int spacedim,
typename Number>
1900 const unsigned int first_selected_component)
1903 fe->n_components() + 1);
1905 shapes.reserve(100);
1907 bool same_base_element =
true;
1908 unsigned int base_element_number = 0;
1909 component_in_base_element = 0;
1910 unsigned int component = 0;
1911 for (; base_element_number < fe->n_base_elements(); ++base_element_number)
1912 if (component + fe->element_multiplicity(base_element_number) >
1913 first_selected_component)
1915 if (first_selected_component + n_components >
1916 component + fe->element_multiplicity(base_element_number))
1917 same_base_element =
false;
1918 component_in_base_element = first_selected_component - component;
1922 component += fe->element_multiplicity(base_element_number);
1926 *fe, base_element_number) &&
1929 shape_info.reinit(
QMidpoint<1>(), *fe, base_element_number);
1930 renumber = shape_info.lexicographic_numbering;
1931 dofs_per_component = shape_info.dofs_per_component_on_cell;
1932 dofs_per_component_face = shape_info.dofs_per_component_on_face;
1934 fe->base_element(base_element_number));
1936 bool is_lexicographic =
true;
1937 for (
unsigned int i = 0; i < renumber.size(); ++i)
1938 if (i != renumber[i])
1939 is_lexicographic =
false;
1941 if (is_lexicographic)
1944 use_linear_path = (poly.size() == 2 && poly[0].value(0.) == 1. &&
1945 poly[0].value(1.) == 0. && poly[1].value(0.) == 0. &&
1946 poly[1].value(1.) == 1.) &&
1947 (fe->n_components() == n_components);
1949 const unsigned int size_face = 3 * dofs_per_component_face * n_components;
1950 const unsigned int size_cell = dofs_per_component * n_components;
1951 scratch_data_scalar.resize(size_face + size_cell);
1953 solution_renumbered.resize(dofs_per_component);
1954 solution_renumbered_vectorized.resize(dofs_per_component);
1960 nonzero_shape_function_component.resize(fe->n_dofs_per_cell());
1961 for (
unsigned int d = 0; d < n_components; ++d)
1963 const unsigned int component = first_selected_component + d;
1964 for (
unsigned int i = 0; i < fe->n_dofs_per_cell(); ++i)
1966 const bool is_primitive =
1967 fe->is_primitive() || fe->is_primitive(i);
1969 nonzero_shape_function_component[i][d] =
1970 (component == fe->system_to_component_index(i).first);
1972 nonzero_shape_function_component[i][d] =
1973 (fe->get_nonzero_components(i)[component] ==
true);
1983template <
int n_components_,
int dim,
int spacedim,
typename Number>
1984template <
bool is_face,
bool is_linear>
1988 const unsigned int geometry_index =
1989 mapping_info->template compute_geometry_index_offset<is_face>(
1990 current_cell_index, current_face_number);
1992 cell_type = mapping_info->get_cell_type(geometry_index);
1994 const_cast<unsigned int &
>(n_q_points_scalar) =
1995 mapping_info->get_n_q_points_unvectorized(geometry_index);
1998 const_cast<unsigned int &
>(n_q_batches) =
1999 (n_q_points_scalar + n_lanes_internal - 1) / n_lanes_internal;
2001 const unsigned int n_q_points_before = n_q_points;
2003 const_cast<unsigned int &
>(n_q_points) =
2004 (stride == 1) ? n_q_batches : n_q_points_scalar;
2006 if (n_q_points != n_q_points_before)
2009 values.resize(n_q_points);
2011 gradients.resize(n_q_points);
2014 if (n_q_points == 0)
2016 is_reinitialized =
true;
2021 const unsigned int unit_point_offset =
2022 mapping_info->compute_unit_point_index_offset(geometry_index);
2025 unit_point_faces_ptr =
2026 mapping_info->get_unit_point_faces(unit_point_offset);
2028 unit_point_ptr = mapping_info->get_unit_point(unit_point_offset);
2031 const unsigned int data_offset =
2032 mapping_info->compute_data_index_offset(geometry_index);
2033 const unsigned int compressed_data_offset =
2034 mapping_info->compute_compressed_data_index_offset(geometry_index);
2037 mapping_info->get_update_flags_mapping();
2039 real_point_ptr = mapping_info->get_real_point(data_offset);
2042 mapping_info->get_jacobian(compressed_data_offset, is_interior);
2044 inverse_jacobian_ptr =
2045 mapping_info->get_inverse_jacobian(compressed_data_offset, is_interior);
2047 normal_ptr = mapping_info->get_normal_vector(data_offset);
2049 JxW_ptr = mapping_info->get_JxW(data_offset);
2051 real_point_ptr = mapping_info->get_real_point(data_offset);
2053 mapping_info->get_jacobian(compressed_data_offset, is_interior);
2054 inverse_jacobian_ptr =
2055 mapping_info->get_inverse_jacobian(compressed_data_offset, is_interior);
2056 normal_ptr = mapping_info->get_normal_vector(data_offset);
2057 JxW_ptr = mapping_info->get_JxW(data_offset);
2060 if (!is_linear && fast_path)
2062 const std::size_t n_shapes = poly.size();
2064 shapes_faces.resize_fast(n_q_batches * n_shapes);
2066 shapes.resize_fast(n_q_batches * n_shapes);
2068 for (
unsigned int qb = 0; qb < n_q_batches; ++qb)
2074 shapes_faces.data() + qb * n_shapes,
2076 unit_point_faces_ptr[qb],
2089 else if (qb + 1 < n_q_batches)
2094 shapes.data() + qb * n_shapes,
2097 unit_point_ptr[qb + 1]);
2110 is_reinitialized =
true;
2115template <
int n_components_,
int dim,
int spacedim,
typename Number>
2119 Number>::value_type &
2121 const unsigned int point_index)
const
2124 return values[point_index];
2129template <
int n_components_,
int dim,
int spacedim,
typename Number>
2133 Number>::gradient_type &
2135 const unsigned int point_index)
const
2138 return gradients[point_index];
2143template <
int n_components_,
int dim,
int spacedim,
typename Number>
2147 const unsigned int point_index)
2150 values[point_index] = value;
2155template <
int n_components_,
int dim,
int spacedim,
typename Number>
2159 const unsigned int point_index)
2162 gradients[point_index] = gradient;
2167template <
int n_components_,
int dim,
int spacedim,
typename Number>
2170 const unsigned int point_index)
const
2173 Assert(jacobian_ptr !=
nullptr,
2175 ExcFEPointEvaluationAccessToUninitializedMappingField(
2176 "update_jacobians"));
2177 return jacobian_ptr[cell_type <= ::internal::MatrixFreeFunctions::
2178 GeometryType::affine ?
2185template <
int n_components_,
int dim,
int spacedim,
typename Number>
2188 const unsigned int point_index)
const
2191 Assert(inverse_jacobian_ptr !=
nullptr,
2193 ExcFEPointEvaluationAccessToUninitializedMappingField(
2194 "update_inverse_jacobians"));
2195 return inverse_jacobian_ptr
2204template <
int n_components_,
int dim,
int spacedim,
typename Number>
2207 const unsigned int point_index)
const
2210 Assert(JxW_ptr !=
nullptr,
2212 ExcFEPointEvaluationAccessToUninitializedMappingField(
2213 "update_JxW_values"));
2214 return JxW_ptr[point_index];
2219template <
int n_components_,
int dim,
int spacedim,
typename Number>
2222 const unsigned int point_index)
const
2224 return quadrature_point(point_index);
2230template <
int n_components_,
int dim,
int spacedim,
typename Number>
2233 const unsigned int point_index)
const
2236 Assert(real_point_ptr !=
nullptr,
2238 ExcFEPointEvaluationAccessToUninitializedMappingField(
2239 "update_quadrature_points"));
2240 return real_point_ptr[point_index];
2245template <
int n_components_,
int dim,
int spacedim,
typename Number>
2248 const unsigned int point_index)
const
2251 Assert(unit_point_ptr !=
nullptr,
ExcMessage(
"unit_point_ptr is not set!"));
2253 for (
unsigned int d = 0; d < dim; ++d)
2255 unit_point_ptr[point_index / stride][d], point_index % stride);
2261template <
int n_components_,
int dim,
int spacedim,
typename Number>
2266 return {0U, n_q_points};
2271template <
int n_components_,
int dim,
int spacedim,
typename Number>
2275 const unsigned int first_selected_component)
2279 first_selected_component)
2284template <
int n_components_,
int dim,
int spacedim,
typename Number>
2289 const unsigned int first_selected_component)
2294 first_selected_component)
2299template <
int n_components_,
int dim,
int spacedim,
typename Number>
2306 if (this->use_linear_path)
2307 this->
template do_reinit<false, true>();
2309 this->
template do_reinit<false, false>();
2314template <
int n_components_,
int dim,
int spacedim,
typename Number>
2321 AssertThrow(this->mapping_info_on_the_fly.get() !=
nullptr,
2324 this->mapping_info->reinit(cell, unit_points);
2326 if (!this->fast_path)
2328 this->fe_values = std::make_shared<FEValues<dim, spacedim>>(
2332 std::vector<
Point<dim>>(unit_points.begin(), unit_points.end())),
2333 this->update_flags);
2334 this->fe_values->reinit(cell);
2337 if (this->use_linear_path)
2338 this->
template do_reinit<false, true>();
2340 this->
template do_reinit<false, false>();
2345template <
int n_components_,
int dim,
int spacedim,
typename Number>
2353 if (this->use_linear_path)
2354 this->
template do_reinit<false, true>();
2356 this->
template do_reinit<false, false>();
2358 if (!this->fast_path)
2360 std::vector<Point<dim>> unit_points(this->n_q_points_scalar);
2362 for (
unsigned int v = 0; v < this->n_q_points_scalar; ++v)
2363 for (
unsigned int d = 0; d < dim; ++d)
2365 this->unit_point_ptr[v / n_lanes_internal][d][v % n_lanes_internal];
2367 this->fe_values = std::make_shared<FEValues<dim, spacedim>>(
2371 std::vector<
Point<dim>>(unit_points.begin(), unit_points.end())),
2372 this->update_flags);
2374 this->fe_values->reinit(
2375 this->mapping_info->get_cell_iterator(this->current_cell_index));
2381template <
int n_components_,
int dim,
int spacedim,
typename Number>
2382template <std::
size_t str
ide_view>
2388 if (!this->is_reinitialized)
2391 if (this->n_q_points == 0)
2401 if (this->fast_path)
2403 if (this->use_linear_path)
2404 evaluate_fast<true>(solution_values, evaluation_flags);
2406 evaluate_fast<false>(solution_values, evaluation_flags);
2409 evaluate_slow(solution_values, evaluation_flags);
2414template <
int n_components_,
int dim,
int spacedim,
typename Number>
2421 solution_values.
size()),
2427template <
int n_components_,
int dim,
int spacedim,
typename Number>
2428template <std::
size_t str
ide_view>
2433 const bool sum_into_values)
2435 do_integrate<true>(solution_values, integration_flags, sum_into_values);
2440template <
int n_components_,
int dim,
int spacedim,
typename Number>
2445 const bool sum_into_values)
2448 solution_values.
size()),
2455template <
int n_components_,
int dim,
int spacedim,
typename Number>
2463 for (
const auto point_index : this->quadrature_point_indices())
2464 return_value += values[point_index] * this->JxW(point_index);
2466 return ETT::sum_value(return_value);
2471template <
int n_components_,
int dim,
int spacedim,
typename Number>
2478 "Calling this function only makes sense in fully vectorized mode."));
2479 if (q == n_q_batches - 1)
2481 const unsigned int n_filled_lanes =
2482 n_q_points_scalar & (n_lanes_user_interface - 1);
2483 if (n_filled_lanes == 0)
2484 return n_lanes_user_interface;
2486 return n_filled_lanes;
2489 return n_lanes_user_interface;
2494template <
int n_components_,
int dim,
int spacedim,
typename Number>
2495template <std::
size_t str
ide_view>
2500 const bool sum_into_values)
2502 do_integrate<false>(solution_values, integration_flags, sum_into_values);
2507template <
int n_components_,
int dim,
int spacedim,
typename Number>
2512 const bool sum_into_values)
2515 solution_values.
size()),
2522template <
int n_components_,
int dim,
int spacedim,
typename Number>
2523template <
bool is_linear, std::
size_t str
ide_view>
2528 const unsigned int dofs_per_comp =
2531 for (
unsigned int comp = 0; comp < n_components; ++comp)
2533 const std::size_t offset =
2534 (this->component_in_base_element + comp) * dofs_per_comp;
2536 if ((is_linear && n_components == 1) || this->renumber.empty())
2538 for (
unsigned int i = 0; i < dofs_per_comp; ++i)
2539 ETT::read_value(solution_values[i + offset],
2541 this->solution_renumbered[i]);
2545 const unsigned int *renumber_ptr = this->renumber.data() + offset;
2546 for (
unsigned int i = 0; i < dofs_per_comp; ++i)
2547 ETT::read_value(solution_values[renumber_ptr[i]],
2549 this->solution_renumbered[i]);
2556template <
int n_components_,
int dim,
int spacedim,
typename Number>
2557template <
bool is_linear, std::
size_t str
ide_view>
2562 const unsigned int n_shapes,
2563 const unsigned int qb,
2569 std::array<vectorized_value_type, dim + 1> result;
2570 if constexpr (is_linear)
2572 if constexpr (n_components == 1)
2579 stride_view>(solution_values.
data(), this->unit_point_ptr[qb]);
2583 this->solution_renumbered.data(), this->unit_point_ptr[qb]);
2591 false>(this->shapes.data() + qb * n_shapes,
2593 this->solution_renumbered.data());
2594 gradient[0] = result[0];
2596 gradient[1] = result[1];
2598 gradient[2] = result[2];
2599 value = result[dim];
2603 if constexpr (is_linear)
2605 if constexpr (n_components == 1)
2610 stride_view>(solution_values.
data(), this->unit_point_ptr[qb]);
2613 this->solution_renumbered.data(), this->unit_point_ptr[qb]);
2621 this->shapes.data() + qb * n_shapes,
2623 this->solution_renumbered.data());
2629template <
int n_components_,
int dim,
int spacedim,
typename Number>
2630template <
bool is_linear, std::
size_t str
ide_view>
2636 if (!(is_linear && n_components == 1))
2637 prepare_evaluate_fast<is_linear>(solution_values);
2640 const unsigned int n_shapes = is_linear ? 2 : this->poly.size();
2642 for (
unsigned int qb = 0; qb < this->n_q_batches; ++qb)
2647 compute_evaluate_fast<is_linear>(
2648 solution_values, evaluation_flags, n_shapes, qb, value, gradient);
2652 for (
unsigned int v = 0, offset = qb * stride;
2653 v < stride && (stride == 1 || offset < this->n_q_points_scalar);
2655 ETT::set_value(value, v, this->values[offset]);
2663 for (
unsigned int v = 0, offset = qb * stride;
2664 v < stride && (stride == 1 || offset < this->n_q_points_scalar);
2668 ETT::set_gradient(gradient, v, unit_gradient);
2669 this->gradients[offset] =
2673 this->inverse_jacobian_ptr[0].
transpose(), unit_gradient) :
2676 ->inverse_jacobian_ptr[this->cell_type <=
2678 GeometryType::affine ?
2690template <
int n_components_,
int dim,
int spacedim,
typename Number>
2691template <std::
size_t str
ide_view>
2698 Assert(this->fe_values.get() !=
nullptr,
2700 "Not initialized. Please call FEPointEvaluation::reinit()!"));
2702 const std::size_t n_points = this->fe_values->get_quadrature().size();
2706 this->values.resize(this->n_q_points);
2707 std::fill(this->values.begin(), this->values.end(),
value_type());
2708 for (
unsigned int i = 0; i < this->fe->n_dofs_per_cell(); ++i)
2711 for (
unsigned int d = 0; d < n_components; ++d)
2712 if (this->nonzero_shape_function_component[i][d] &&
2713 (this->fe->is_primitive(i) || this->fe->is_primitive()))
2714 for (
unsigned int qb = 0, q = 0; q < n_points;
2715 ++qb, q += n_lanes_user_interface)
2716 for (
unsigned int v = 0;
2717 v < n_lanes_user_interface && q + v < n_points;
2719 ETT::access(this->values[qb],
2722 this->fe_values->shape_value(i, q + v) * value);
2723 else if (this->nonzero_shape_function_component[i][d])
2724 for (
unsigned int qb = 0, q = 0; q < n_points;
2725 ++qb, q += n_lanes_user_interface)
2726 for (
unsigned int v = 0;
2727 v < n_lanes_user_interface && q + v < n_points;
2729 ETT::access(this->values[qb],
2732 this->fe_values->shape_value_component(i,
2741 this->gradients.resize(this->n_q_points);
2742 std::fill(this->gradients.begin(),
2743 this->gradients.end(),
2745 for (
unsigned int i = 0; i < this->fe->n_dofs_per_cell(); ++i)
2748 for (
unsigned int d = 0; d < n_components; ++d)
2749 if (this->nonzero_shape_function_component[i][d] &&
2750 (this->fe->is_primitive(i) || this->fe->is_primitive()))
2751 for (
unsigned int qb = 0, q = 0; q < n_points;
2752 ++qb, q += n_lanes_user_interface)
2753 for (
unsigned int v = 0;
2754 v < n_lanes_user_interface && q + v < n_points;
2756 ETT::access(this->gradients[qb],
2759 this->fe_values->shape_grad(i, q + v) * value);
2760 else if (this->nonzero_shape_function_component[i][d])
2761 for (
unsigned int qb = 0, q = 0; q < n_points;
2762 ++qb, q += n_lanes_user_interface)
2763 for (
unsigned int v = 0;
2764 v < n_lanes_user_interface && q + v < n_points;
2767 this->gradients[qb],
2770 this->fe_values->shape_grad_component(i, q + v, d) * value);
2777template <
int n_components_,
int dim,
int spacedim,
typename Number>
2778template <
bool is_linear>
2782 const unsigned int n_shapes,
2783 const unsigned int qb,
2798 solution_values_vectorized_linear :
2799 this->solution_renumbered_vectorized.data(),
2800 this->unit_point_ptr[qb],
2807 this->shapes.data() + qb * n_shapes,
2810 is_linear ? solution_values_vectorized_linear :
2811 this->solution_renumbered_vectorized.data(),
2812 this->unit_point_ptr[qb],
2818template <
int n_components_,
int dim,
int spacedim,
typename Number>
2819template <
bool is_linear, std::
size_t str
ide_view>
2824 const bool sum_into_values)
2826 if (!sum_into_values && this->fe->n_components() > n_components)
2827 for (
unsigned int i = 0; i < solution_values.
size(); ++i)
2828 solution_values[i] = 0;
2830 const unsigned int dofs_per_comp =
2833 for (
unsigned int comp = 0; comp < n_components; ++comp)
2835 const std::size_t offset =
2836 (this->component_in_base_element + comp) * dofs_per_comp;
2838 if (is_linear || this->renumber.empty())
2840 for (
unsigned int i = 0; i < dofs_per_comp; ++i)
2841 if (sum_into_values)
2842 solution_values[i + offset] +=
2843 ETT::sum_value(comp,
2845 *(solution_values_vectorized_linear + i) :
2846 this->solution_renumbered_vectorized[i]);
2848 solution_values[i + offset] =
2849 ETT::sum_value(comp,
2851 *(solution_values_vectorized_linear + i) :
2852 this->solution_renumbered_vectorized[i]);
2856 const unsigned int *renumber_ptr = this->renumber.
data() + offset;
2857 for (
unsigned int i = 0; i < dofs_per_comp; ++i)
2858 if (sum_into_values)
2859 solution_values[renumber_ptr[i]] +=
2860 ETT::sum_value(comp, this->solution_renumbered_vectorized[i]);
2862 solution_values[renumber_ptr[i]] =
2863 ETT::sum_value(comp, this->solution_renumbered_vectorized[i]);
2870template <
int n_components_,
int dim,
int spacedim,
typename Number>
2871template <
bool do_JxW,
bool is_linear, std::
size_t str
ide_view>
2876 const bool sum_into_values)
2879 if constexpr (stride == 1)
2880 if (
const unsigned int n_filled_lanes =
2881 this->n_q_points_scalar & (n_lanes_internal - 1);
2885 for (
unsigned int v = n_filled_lanes; v < n_lanes_internal; ++v)
2886 ETT::set_zero_value(this->values.back(), v);
2888 for (
unsigned int v = n_filled_lanes; v < n_lanes_internal; ++v)
2889 ETT::set_zero_gradient(this->gradients.back(), v);
2893 solution_values_vectorized_linear = {};
2896 const unsigned int n_shapes = is_linear ? 2 : this->poly.size();
2898 const bool cartesian_cell =
2900 const bool affine_cell =
2902 for (
unsigned int qb = 0; qb < this->n_q_batches; ++qb)
2908 for (
unsigned int v = 0, offset = qb * stride;
2909 v < stride && (stride == 1 || offset < this->n_q_points_scalar);
2911 ETT::get_value(value,
2913 do_JxW ? this->values[offset] * this->JxW_ptr[offset] :
2914 this->values[offset]);
2917 for (
unsigned int v = 0, offset = qb * stride;
2918 v < stride && (stride == 1 || offset < this->n_q_points_scalar);
2922 do_JxW ? this->gradients[offset] * this->JxW_ptr[offset] :
2923 this->gradients[offset];
2931 this->inverse_jacobian_ptr[affine_cell ? 0 : offset],
2935 compute_integrate_fast<is_linear>(
2941 solution_values_vectorized_linear.data());
2945 finish_integrate_fast<is_linear>(solution_values,
2946 solution_values_vectorized_linear.
data(),
2952template <
int n_components_,
int dim,
int spacedim,
typename Number>
2953template <
bool do_JxW, std::
size_t str
ide_view>
2958 const bool sum_into_values)
2961 Assert(this->fe_values.get() !=
nullptr,
2963 "Not initialized. Please call FEPointEvaluation::reinit()!"));
2964 if (!sum_into_values)
2965 for (
unsigned int i = 0; i < solution_values.
size(); ++i)
2966 solution_values[i] = 0;
2968 const std::size_t n_points = this->fe_values->get_quadrature().
size();
2973 for (
unsigned int i = 0; i < this->fe->n_dofs_per_cell(); ++i)
2975 for (
unsigned int d = 0; d < n_components; ++d)
2976 if (this->nonzero_shape_function_component[i][d] &&
2977 (this->fe->is_primitive(i) || this->fe->is_primitive()))
2978 for (
unsigned int qb = 0, q = 0; q < n_points;
2979 ++qb, q += n_lanes_user_interface)
2980 for (
unsigned int v = 0;
2981 v < n_lanes_user_interface && q + v < n_points;
2983 solution_values[i] +=
2984 this->fe_values->shape_value(i, q + v) *
2985 ETT::access(this->values[qb], v, d) *
2986 (do_JxW ? this->fe_values->JxW(q + v) : 1.);
2987 else if (this->nonzero_shape_function_component[i][d])
2988 for (
unsigned int qb = 0, q = 0; q < n_points;
2989 ++qb, q += n_lanes_user_interface)
2990 for (
unsigned int v = 0;
2991 v < n_lanes_user_interface && q + v < n_points;
2993 solution_values[i] +=
2994 this->fe_values->shape_value_component(i, q + v, d) *
2995 ETT::access(this->values[qb], v, d) *
2996 (do_JxW ? this->fe_values->JxW(q + v) : 1.);
3003 for (
unsigned int i = 0; i < this->fe->n_dofs_per_cell(); ++i)
3005 for (
unsigned int d = 0; d < n_components; ++d)
3006 if (this->nonzero_shape_function_component[i][d] &&
3007 (this->fe->is_primitive(i) || this->fe->is_primitive()))
3008 for (
unsigned int qb = 0, q = 0; q < n_points;
3009 ++qb, q += n_lanes_user_interface)
3010 for (
unsigned int v = 0;
3011 v < n_lanes_user_interface && q + v < n_points;
3013 solution_values[i] +=
3014 this->fe_values->shape_grad(i, q + v) *
3015 ETT::access(this->gradients[qb], v, d) *
3016 (do_JxW ? this->fe_values->JxW(q + v) : 1.);
3017 else if (this->nonzero_shape_function_component[i][d])
3018 for (
unsigned int qb = 0, q = 0; q < n_points;
3019 ++qb, q += n_lanes_user_interface)
3020 for (
unsigned int v = 0;
3021 v < n_lanes_user_interface && q + v < n_points;
3023 solution_values[i] +=
3024 this->fe_values->shape_grad_component(i, q + v, d) *
3025 ETT::access(this->gradients[qb], v, d) *
3026 (do_JxW ? this->fe_values->JxW(q + v) : 1.);
3033template <
int n_components_,
int dim,
int spacedim,
typename Number>
3034template <
bool do_JxW, std::
size_t str
ide_view>
3039 const bool sum_into_values)
3041 if (!this->is_reinitialized)
3046 if (this->n_q_points == 0 ||
3048 (integration_flags &
3051 if (!sum_into_values)
3052 for (
unsigned int i = 0; i < solution_values.
size(); ++i)
3053 solution_values[i] = 0;
3058 !do_JxW || this->JxW_ptr !=
nullptr,
3060 "JxW pointer is not set! If you do not want to integrate() use test_and_sum()"));
3063 if (this->fast_path)
3065 if (this->use_linear_path)
3066 integrate_fast<do_JxW, true>(solution_values,
3070 integrate_fast<do_JxW, false>(solution_values,
3075 integrate_slow<do_JxW>(solution_values, integration_flags, sum_into_values);
3080template <
int n_components_,
int dim,
int spacedim,
typename Number>
3083 const unsigned int point_index)
const
3086 Assert(this->normal_ptr !=
nullptr,
3088 ExcFEPointEvaluationAccessToUninitializedMappingField(
3089 "update_normal_vectors"));
3090 if (this->is_interior)
3091 return this->normal_ptr[point_index];
3093 return -this->normal_ptr[point_index];
3098template <
int n_components_,
int dim,
int spacedim,
typename Number>
3103 const bool is_interior,
3104 const unsigned int first_selected_component)
3108 first_selected_component,
3114template <
int n_components_,
int dim,
int spacedim,
typename Number>
3118 const unsigned int face_number)
3121 this->current_face_number = face_number;
3123 if (this->use_linear_path)
3124 this->
template do_reinit<true, true>();
3126 this->
template do_reinit<true, false>();
3131template <
int n_components_,
int dim,
int spacedim,
typename Number>
3134 const unsigned int face_index)
3136 this->current_cell_index = face_index;
3137 this->current_face_number =
3138 this->mapping_info->get_face_number(face_index, this->is_interior);
3140 if (this->use_linear_path)
3141 this->
template do_reinit<true, true>();
3143 this->
template do_reinit<true, false>();
3148template <
int n_components_,
int dim,
int spacedim,
typename Number>
3149template <std::
size_t str
ide_view>
3157 if (this->n_q_points == 0)
3168 if (this->use_linear_path)
3169 do_evaluate<true>(solution_values, evaluation_flags);
3171 do_evaluate<false>(solution_values, evaluation_flags);
3176template <
int n_components_,
int dim,
int spacedim,
typename Number>
3183 solution_values.
size()),
3189template <
int n_components_,
int dim,
int spacedim,
typename Number>
3190template <
bool is_linear, std::
size_t str
ide_view>
3196 const unsigned int dofs_per_comp =
3200 if (stride_view == 1 && this->component_in_base_element == 0 &&
3201 (is_linear || this->renumber.empty()))
3202 input = solution_values.
data();
3205 for (
unsigned int comp = 0; comp < n_components; ++comp)
3207 const std::size_t offset =
3208 (this->component_in_base_element + comp) * dofs_per_comp;
3210 if (is_linear || this->renumber.empty())
3212 for (
unsigned int i = 0; i < dofs_per_comp; ++i)
3213 this->scratch_data_scalar[i + comp * dofs_per_comp] =
3214 solution_values[i + offset];
3218 const unsigned int *renumber_ptr = this->renumber.
data() + offset;
3219 for (
unsigned int i = 0; i < dofs_per_comp; ++i)
3220 this->scratch_data_scalar[i + comp * dofs_per_comp] =
3221 solution_values[renumber_ptr[i]];
3224 input = this->scratch_data_scalar.
data();
3228 this->scratch_data_scalar.begin() + dofs_per_comp * n_components;
3231 template interpolate<true, false>(n_components,
3236 this->current_face_number);
3238 do_evaluate_in_face<is_linear, 1>(output, evaluation_flags);
3243template <
int n_components_,
int dim,
int spacedim,
typename Number>
3244template <std::
size_t str
ide_view>
3249 const bool sum_into_values)
3255 if (this->n_q_points == 0 ||
3257 (integration_flags &
3260 if (!sum_into_values)
3261 for (
unsigned int i = 0; i < solution_values.
size(); ++i)
3262 solution_values[i] = 0;
3268 if (this->use_linear_path)
3269 do_integrate<true, true>(solution_values,
3273 do_integrate<true, false>(solution_values,
3280template <
int n_components_,
int dim,
int spacedim,
typename Number>
3285 const bool sum_into_values)
3288 solution_values.
size()),
3295template <
int n_components_,
int dim,
int spacedim,
typename Number>
3296template <std::
size_t str
ide_view>
3301 const bool sum_into_values)
3307 if (this->n_q_points == 0 ||
3309 (integration_flags &
3312 if (!sum_into_values)
3313 for (
unsigned int i = 0; i < solution_values.
size(); ++i)
3314 solution_values[i] = 0;
3320 if (this->use_linear_path)
3321 do_integrate<false, true>(solution_values,
3325 do_integrate<false, false>(solution_values,
3332template <
int n_components_,
int dim,
int spacedim,
typename Number>
3337 const bool sum_into_values)
3340 solution_values.
size()),
3347template <
int n_components_,
int dim,
int spacedim,
typename Number>
3348template <
bool do_JxW,
bool is_linear, std::
size_t str
ide_view>
3353 const bool sum_into_values)
3355 if (!sum_into_values && this->fe->n_components() > n_components)
3356 for (
unsigned int i = 0; i < solution_values.
size(); ++i)
3357 solution_values[i] = 0;
3359 do_integrate_in_face<do_JxW, is_linear, 1>(this->scratch_data_scalar.begin(),
3363 ScalarNumber *input = this->scratch_data_scalar.begin();
3365 if (stride_view == 1 && this->component_in_base_element == 0 &&
3366 (is_linear || this->renumber.empty()))
3368 if (sum_into_values)
3370 FEFaceNormalEvaluationImpl<dim, is_linear ? 1 : -1,
ScalarNumber>::
3371 template interpolate<false, true>(n_components,
3375 solution_values.
data(),
3376 this->current_face_number);
3379 FEFaceNormalEvaluationImpl<dim, is_linear ? 1 : -1,
ScalarNumber>::
3380 template interpolate<false, false>(n_components,
3384 solution_values.
data(),
3385 this->current_face_number);
3389 const unsigned int dofs_per_comp_face =
3390 is_linear ?
Utilities::pow(2, dim - 1) : this->dofs_per_component_face;
3392 const unsigned int size_input = 3 * dofs_per_comp_face * n_components;
3396 FEFaceNormalEvaluationImpl<dim, is_linear ? 1 : -1,
ScalarNumber>::
3397 template interpolate<false, false>(n_components,
3402 this->current_face_number);
3404 const unsigned int dofs_per_comp =
3407 for (
unsigned int comp = 0; comp < n_components; ++comp)
3409 const std::size_t offset =
3410 (this->component_in_base_element + comp) * dofs_per_comp;
3412 if (is_linear || this->renumber.empty())
3414 for (
unsigned int i = 0; i < dofs_per_comp; ++i)
3415 if (sum_into_values)
3416 solution_values[i + offset] +=
3417 output[i + comp * dofs_per_comp];
3419 solution_values[i + offset] =
3420 output[i + comp * dofs_per_comp];
3424 const unsigned int *renumber_ptr = this->renumber.
data() + offset;
3425 for (
unsigned int i = 0; i < dofs_per_comp; ++i)
3426 if (sum_into_values)
3427 solution_values[renumber_ptr[i]] +=
3428 output[i + comp * dofs_per_comp];
3430 solution_values[renumber_ptr[i]] =
3431 output[i + comp * dofs_per_comp];
3439template <
int n_components_,
int dim,
int spacedim,
typename Number>
3440template <
int str
ide_face_dof>
3446 if (this->use_linear_path)
3447 do_evaluate_in_face<true, stride_face_dof>(face_dof_values,
3450 do_evaluate_in_face<false, stride_face_dof>(face_dof_values,
3456template <
int n_components_,
int dim,
int spacedim,
typename Number>
3457template <
bool is_linear,
int str
ide_face_dof>
3464 if constexpr (n_components == 1)
3465 face_dof_values_ptr = face_dof_values;
3468 const unsigned int dofs_per_comp_face =
3469 is_linear ?
Utilities::pow(2, dim - 1) : this->dofs_per_component_face;
3470 for (
unsigned int comp = 0; comp < n_components; ++comp)
3471 for (
unsigned int i = 0; i < 2 * dofs_per_comp_face; ++i)
3472 ETT::read_value(face_dof_values[(i + comp * 3 * dofs_per_comp_face) *
3475 this->solution_renumbered[i]);
3477 face_dof_values_ptr = this->solution_renumbered.data();
3480 constexpr int stride_face_dof_actual =
3481 n_components == 1 ? stride_face_dof : 1;
3484 const unsigned int n_shapes = is_linear ? 2 : this->poly.size();
3486 for (
unsigned int qb = 0; qb < this->n_q_batches; ++qb)
3493 const std::array<vectorized_value_type, dim + 1> interpolated_value =
3500 stride_face_dof_actual>(face_dof_values_ptr,
3501 this->unit_point_faces_ptr[qb]) :
3508 stride_face_dof_actual>(this->shapes_faces.data() +
3511 face_dof_values_ptr);
3513 value = interpolated_value[dim - 1];
3516 if (this->current_face_number / 2 == 0)
3518 gradient[0] = interpolated_value[dim];
3520 gradient[1] = interpolated_value[0];
3522 gradient[2] = interpolated_value[1];
3524 else if (this->current_face_number / 2 == 1)
3527 gradient[1] = interpolated_value[dim];
3530 gradient[0] = interpolated_value[1];
3531 gradient[2] = interpolated_value[0];
3534 gradient[0] = interpolated_value[0];
3538 else if (this->current_face_number / 2 == 2)
3542 gradient[0] = interpolated_value[0];
3543 gradient[1] = interpolated_value[1];
3544 gradient[2] = interpolated_value[dim];
3559 stride_face_dof_actual>(face_dof_values_ptr,
3560 this->unit_point_faces_ptr[qb]) :
3566 stride_face_dof_actual>(this->shapes_faces.data() +
3569 face_dof_values_ptr);
3574 for (
unsigned int v = 0, offset = qb * stride;
3575 v < stride && (stride == 1 || offset < this->n_q_points_scalar);
3577 ETT::set_value(value, v, this->values[offset]);
3585 for (
unsigned int v = 0, offset = qb * stride;
3586 v < stride && (stride == 1 || offset < this->n_q_points_scalar);
3590 ETT::set_gradient(gradient, v, unit_gradient);
3591 this->gradients[offset] =
3595 this->inverse_jacobian_ptr[0].
transpose(), unit_gradient) :
3598 ->inverse_jacobian_ptr[this->cell_type <=
3600 GeometryType::affine ?
3612template <
int n_components_,
int dim,
int spacedim,
typename Number>
3613template <
int str
ide_face_dof>
3618 const bool sum_into_values)
3620 if (this->use_linear_path)
3621 do_integrate_in_face<true, true, stride_face_dof>(face_dof_values,
3625 do_integrate_in_face<true, false, stride_face_dof>(face_dof_values,
3632template <
int n_components_,
int dim,
int spacedim,
typename Number>
3633template <
bool do_JxW,
bool is_linear,
int str
ide_face_dof>
3639 const bool sum_into_values)
3642 if constexpr (stride == 1)
3643 if (
const unsigned int n_filled_lanes =
3644 this->n_q_points_scalar & (n_lanes_internal - 1);
3648 for (
unsigned int v = n_filled_lanes; v < n_lanes_internal; ++v)
3649 ETT::set_zero_value(this->values.back(), v);
3651 for (
unsigned int v = n_filled_lanes; v < n_lanes_internal; ++v)
3652 ETT::set_zero_gradient(this->gradients.back(), v);
3657 solution_values_vectorized_linear = {};
3660 const unsigned int n_shapes = is_linear ? 2 : this->poly.size();
3662 const bool cartesian_cell =
3664 const bool affine_cell =
3666 for (
unsigned int qb = 0; qb < this->n_q_batches; ++qb)
3672 for (
unsigned int v = 0, offset = qb * stride;
3673 v < stride && (stride == 1 || offset < this->n_q_points_scalar);
3675 ETT::get_value(value,
3677 do_JxW ? this->values[offset] * this->JxW_ptr[offset] :
3678 this->values[offset]);
3681 for (
unsigned int v = 0, offset = qb * stride;
3682 v < stride && (stride == 1 || offset < this->n_q_points_scalar);
3686 do_JxW ? this->gradients[offset] * this->JxW_ptr[offset] :
3687 this->gradients[offset];
3695 this->inverse_jacobian_ptr[affine_cell ? 0 : offset],
3701 std::array<vectorized_value_type, 2> value_face = {};
3704 value_face[0] = value;
3707 if (this->current_face_number / 2 == 0)
3709 value_face[1] = gradient[0];
3711 gradient_in_face[0] = gradient[1];
3713 gradient_in_face[1] = gradient[2];
3715 else if (this->current_face_number / 2 == 1)
3718 value_face[1] = gradient[1];
3721 gradient_in_face[0] = gradient[2];
3722 gradient_in_face[1] = gradient[0];
3725 gradient_in_face[0] = gradient[0];
3729 else if (this->current_face_number / 2 == 2)
3733 value_face[1] = gradient[2];
3734 gradient_in_face[0] = gradient[0];
3735 gradient_in_face[1] = gradient[1];
3748 2>(this->shapes_faces.data() + qb * n_shapes,
3752 is_linear ? solution_values_vectorized_linear.data() :
3753 this->solution_renumbered_vectorized.data(),
3754 this->unit_point_faces_ptr[qb],
3762 this->shapes_faces.data() + qb * n_shapes,
3765 is_linear ? solution_values_vectorized_linear.data() :
3766 this->solution_renumbered_vectorized.data(),
3767 this->unit_point_faces_ptr[qb],
3771 const unsigned int dofs_per_comp_face =
3772 is_linear ?
Utilities::pow(2, dim - 1) : this->dofs_per_component_face;
3774 for (
unsigned int comp = 0; comp < n_components; ++comp)
3775 for (
unsigned int i = 0; i < 2 * dofs_per_comp_face; ++i)
3776 if (sum_into_values)
3777 face_dof_values[(i + comp * 3 * dofs_per_comp_face) *
3779 ETT::sum_value(comp,
3781 *(solution_values_vectorized_linear.data() + i) :
3782 this->solution_renumbered_vectorized[i]);
3784 face_dof_values[(i + comp * 3 * dofs_per_comp_face) * stride_face_dof] =
3785 ETT::sum_value(comp,
3787 *(solution_values_vectorized_linear.data() + i) :
3788 this->solution_renumbered_vectorized[i]);
3793template <
int n_components_,
int dim,
int spacedim,
typename Number>
3796 const unsigned int point_index)
const
3799 Assert(this->normal_ptr !=
nullptr,
3801 ExcFEPointEvaluationAccessToUninitializedMappingField(
3802 "update_normal_vectors"));
3806 for (
unsigned int d = 0; d < dim; ++d)
3810 if (this->is_interior)
3817 if (this->is_interior)
3818 return this->normal_ptr[point_index];
3820 return -(this->normal_ptr[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)
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
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)
FEFacePointEvaluation(NonMatching::MappingInfo< dim, spacedim, Number > &mapping_info, const FiniteElement< dim, spacedim > &fe, const bool is_interior=true, const unsigned int first_selected_component=0)
typename ETT::unit_gradient_type unit_gradient_type
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
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
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
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
SmartPointer< NonMatching::MappingInfo< dim, spacedim, Number > > mapping_info
const unsigned int n_q_batches
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)
SmartPointer< const FiniteElement< dim, spacedim > > fe
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
boost::signals2::connection connection_is_reinitialized
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
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)
FEPointEvaluationBase(NonMatching::MappingInfo< dim, spacedim, Number > &mapping_info, const FiniteElement< dim, spacedim > &fe, const unsigned int first_selected_component=0, const bool is_interior=true)
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
unsigned int component_in_base_element
SmartPointer< const Mapping< dim, spacedim > > mapping
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)
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 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 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_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()
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)