36#include <boost/container/small_vector.hpp>
47 template <
int dim,
int spacedim>
48 inline std::vector<unsigned int>
51 std::vector<unsigned int> shape_function_to_row_table(
60 unsigned int nth_nonzero_component = 0;
65 row + nth_nonzero_component;
66 ++nth_nonzero_component;
71 return shape_function_to_row_table;
78 template <
typename Number,
typename T =
void>
93 template <
typename Number>
96 std::enable_if_t<Differentiation::AD::is_ad_number<Number>::value>>
99 value(
const Number & )
110template <
int dim,
int spacedim>
118template <
int dim,
int spacedim>
126template <
int dim,
int spacedim>
134template <
int dim,
int spacedim>
138 return cell.has_value();
143template <
int dim,
int spacedim>
152 [](
auto &cell_iterator) ->
154 return cell_iterator;
161template <
int dim,
int spacedim>
167 switch (cell.value().index())
170 return std::get<1>(cell.value())->get_dof_handler().n_dofs();
172 return std::get<2>(cell.value())->get_dof_handler().n_dofs();
174 Assert(
false, ExcNeedsDoFHandler());
181template <
int dim,
int spacedim>
182template <
typename Number>
190 switch (cell.value().index())
193 std::get<1>(cell.value())->get_interpolated_dof_values(in, out);
197 std::get<2>(cell.value())->get_interpolated_dof_values(in, out);
201 Assert(
false, ExcNeedsDoFHandler());
211template <
int dim,
int spacedim>
213 const unsigned int n_q_points,
222 ,
fe(&
fe, typeid(*this).name())
228 ExcMessage(
"There is nothing useful you can do with an FEValues "
229 "object when using a quadrature formula with zero "
230 "quadrature points!"));
236template <
int dim,
int spacedim>
245template <
int dim,
int spacedim>
264 template <
typename Number,
typename Number2>
267 const ::Table<2, double> &shape_values,
268 std::vector<Number> &values)
271 const unsigned int dofs_per_cell = shape_values.n_rows();
272 const unsigned int n_quadrature_points = values.size();
275 std::fill_n(values.begin(),
286 for (
unsigned int shape_func = 0; shape_func < dofs_per_cell; ++shape_func)
288 const Number2
value = dof_values[shape_func];
296 const double *shape_value_ptr = &shape_values(shape_func, 0);
297 for (
unsigned int point = 0; point < n_quadrature_points; ++point)
298 values[point] +=
value * (*shape_value_ptr++);
304 template <
int dim,
int spacedim,
typename VectorType>
308 const ::Table<2, double> &shape_values,
310 const std::vector<unsigned int> &shape_function_to_row_table,
312 const bool quadrature_points_fastest =
false,
313 const unsigned int component_multiple = 1)
315 using Number =
typename VectorType::value_type;
317 for (
unsigned int i = 0; i < values.size(); ++i)
318 std::fill_n(values[i].begin(),
320 typename VectorType::value_type());
325 if (dofs_per_cell == 0)
328 const unsigned int n_quadrature_points =
329 quadrature_points_fastest ? values[0].size() : values.size();
333 const unsigned result_components = n_components * component_multiple;
334 (void)result_components;
335 if (quadrature_points_fastest)
338 for (
unsigned int i = 0; i < values.size(); ++i)
344 for (
unsigned int i = 0; i < values.size(); ++i)
351 for (
unsigned int mc = 0; mc < component_multiple; ++mc)
352 for (
unsigned int shape_func = 0; shape_func < dofs_per_cell;
355 const Number &
value = dof_values[shape_func + mc * dofs_per_cell];
359 if (::internal::CheckForZero<Number>::value(
value) ==
true)
364 const unsigned int comp =
367 const unsigned int row =
368 shape_function_to_row_table[shape_func * n_components + comp];
370 const double *shape_value_ptr = &shape_values(row, 0);
372 if (quadrature_points_fastest)
374 VectorType &values_comp = values[comp];
375 for (
unsigned int point = 0; point < n_quadrature_points;
377 values_comp[point] +=
value * (*shape_value_ptr++);
380 for (
unsigned int point = 0; point < n_quadrature_points;
382 values[point][comp] +=
value * (*shape_value_ptr++);
385 for (
unsigned int c = 0; c < n_components; ++c)
390 const unsigned int row =
391 shape_function_to_row_table[shape_func * n_components + c];
393 const double *shape_value_ptr = &shape_values(row, 0);
394 const unsigned int comp = c + mc * n_components;
396 if (quadrature_points_fastest)
398 VectorType &values_comp = values[comp];
399 for (
unsigned int point = 0; point < n_quadrature_points;
401 values_comp[point] +=
value * (*shape_value_ptr++);
404 for (
unsigned int point = 0; point < n_quadrature_points;
406 values[point][comp] +=
value * (*shape_value_ptr++);
415 template <
int order,
int spacedim,
typename Number>
417 do_function_derivatives(
422 const unsigned int dofs_per_cell = shape_derivatives.size()[0];
423 const unsigned int n_quadrature_points = derivatives.size();
426 std::fill_n(derivatives.begin(),
437 for (
unsigned int shape_func = 0; shape_func < dofs_per_cell; ++shape_func)
439 const Number &
value = dof_values[shape_func];
443 if (::internal::CheckForZero<Number>::value(
value) ==
true)
447 &shape_derivatives[shape_func][0];
448 for (
unsigned int point = 0; point < n_quadrature_points; ++point)
449 derivatives[point] +=
value * (*shape_derivative_ptr++);
455 template <
int order,
int dim,
int spacedim,
typename Number>
457 do_function_derivatives(
461 const std::vector<unsigned int> &shape_function_to_row_table,
463 const bool quadrature_points_fastest =
false,
464 const unsigned int component_multiple = 1)
467 for (
unsigned int i = 0; i < derivatives.size(); ++i)
468 std::fill_n(derivatives[i].begin(),
469 derivatives[i].size(),
475 if (dofs_per_cell == 0)
479 const unsigned int n_quadrature_points =
480 quadrature_points_fastest ? derivatives[0].size() : derivatives.size();
484 const unsigned result_components = n_components * component_multiple;
485 (void)result_components;
486 if (quadrature_points_fastest)
489 for (
unsigned int i = 0; i < derivatives.size(); ++i)
495 for (
unsigned int i = 0; i < derivatives.size(); ++i)
502 for (
unsigned int mc = 0; mc < component_multiple; ++mc)
503 for (
unsigned int shape_func = 0; shape_func < dofs_per_cell;
506 const Number &value = dof_values[shape_func + mc * dofs_per_cell];
510 if (::internal::CheckForZero<Number>::value(value) ==
true)
515 const unsigned int comp =
518 const unsigned int row =
519 shape_function_to_row_table[shape_func * n_components + comp];
522 &shape_derivatives[row][0];
524 if (quadrature_points_fastest)
525 for (
unsigned int point = 0; point < n_quadrature_points;
527 derivatives[comp][point] += value * (*shape_derivative_ptr++);
529 for (
unsigned int point = 0; point < n_quadrature_points;
531 derivatives[point][comp] += value * (*shape_derivative_ptr++);
534 for (
unsigned int c = 0; c < n_components; ++c)
539 const unsigned int row =
540 shape_function_to_row_table[shape_func * n_components + c];
543 &shape_derivatives[row][0];
544 const unsigned int comp = c + mc * n_components;
546 if (quadrature_points_fastest)
547 for (
unsigned int point = 0; point < n_quadrature_points;
549 derivatives[comp][point] +=
550 value * (*shape_derivative_ptr++);
552 for (
unsigned int point = 0; point < n_quadrature_points;
554 derivatives[point][comp] +=
555 value * (*shape_derivative_ptr++);
562 template <
int spacedim,
typename Number,
typename Number2>
564 do_function_laplacians(
567 std::vector<Number> &laplacians)
569 const unsigned int dofs_per_cell = shape_hessians.size()[0];
570 const unsigned int n_quadrature_points = laplacians.size();
573 std::fill_n(laplacians.begin(),
580 for (
unsigned int shape_func = 0; shape_func < dofs_per_cell; ++shape_func)
582 const Number2 value = dof_values[shape_func];
591 &shape_hessians[shape_func][0];
592 for (
unsigned int point = 0; point < n_quadrature_points; ++point)
593 laplacians[point] += value *
trace(*shape_hessian_ptr++);
599 template <
int dim,
int spacedim,
typename VectorType,
typename Number>
601 do_function_laplacians(
605 const std::vector<unsigned int> &shape_function_to_row_table,
606 std::vector<VectorType> &laplacians,
607 const bool quadrature_points_fastest =
false,
608 const unsigned int component_multiple = 1)
611 for (
unsigned int i = 0; i < laplacians.size(); ++i)
612 std::fill_n(laplacians[i].begin(),
613 laplacians[i].size(),
614 typename VectorType::value_type());
619 if (dofs_per_cell == 0)
623 const unsigned int n_quadrature_points = laplacians.size();
627 const unsigned result_components = n_components * component_multiple;
628 (void)result_components;
629 if (quadrature_points_fastest)
632 for (
unsigned int i = 0; i < laplacians.size(); ++i)
638 for (
unsigned int i = 0; i < laplacians.size(); ++i)
645 for (
unsigned int mc = 0; mc < component_multiple; ++mc)
646 for (
unsigned int shape_func = 0; shape_func < dofs_per_cell;
649 const Number &value = dof_values[shape_func + mc * dofs_per_cell];
653 if (::internal::CheckForZero<Number>::value(value) ==
true)
658 const unsigned int comp =
661 const unsigned int row =
662 shape_function_to_row_table[shape_func * n_components + comp];
665 &shape_hessians[row][0];
666 if (quadrature_points_fastest)
668 VectorType &laplacians_comp = laplacians[comp];
669 for (
unsigned int point = 0; point < n_quadrature_points;
671 laplacians_comp[point] +=
672 value *
trace(*shape_hessian_ptr++);
675 for (
unsigned int point = 0; point < n_quadrature_points;
677 laplacians[point][comp] +=
678 value *
trace(*shape_hessian_ptr++);
681 for (
unsigned int c = 0; c < n_components; ++c)
686 const unsigned int row =
687 shape_function_to_row_table[shape_func * n_components + c];
690 &shape_hessians[row][0];
691 const unsigned int comp = c + mc * n_components;
693 if (quadrature_points_fastest)
695 VectorType &laplacians_comp = laplacians[comp];
696 for (
unsigned int point = 0; point < n_quadrature_points;
698 laplacians_comp[point] +=
699 value *
trace(*shape_hessian_ptr++);
702 for (
unsigned int point = 0; point < n_quadrature_points;
704 laplacians[point][comp] +=
705 value *
trace(*shape_hessian_ptr++);
713template <
int dim,
int spacedim>
714template <
typename Number>
718 std::vector<Number> &values)
const
731 this->finite_element_output.shape_values,
737template <
int dim,
int spacedim>
738template <
typename Number>
743 std::vector<Number> &values)
const
750 boost::container::small_vector<Number, 200> dof_values(
dofs_per_cell);
760template <
int dim,
int spacedim>
761template <
typename Number>
778 this->finite_element_output.shape_values,
780 this->finite_element_output.shape_function_to_row_table,
786template <
int dim,
int spacedim>
787template <
typename Number>
808 this->finite_element_output.shape_function_to_row_table,
816template <
int dim,
int spacedim>
817template <
typename Number>
823 const bool quadrature_points_fastest)
const
833 boost::container::small_vector<Number, 200> dof_values(
dofs_per_cell);
840 this->finite_element_output.shape_function_to_row_table,
842 quadrature_points_fastest,
848template <
int dim,
int spacedim>
849template <
typename Number>
866 this->finite_element_output.shape_gradients,
872template <
int dim,
int spacedim>
873template <
typename Number>
885 boost::container::small_vector<Number, 200> dof_values(
dofs_per_cell);
895template <
int dim,
int spacedim>
896template <
typename Number>
912 this->finite_element_output.shape_gradients,
914 this->finite_element_output.shape_function_to_row_table,
920template <
int dim,
int spacedim>
921template <
typename Number>
927 const bool quadrature_points_fastest)
const
936 boost::container::small_vector<Number, 200> dof_values(
dofs_per_cell);
943 this->finite_element_output.shape_function_to_row_table,
945 quadrature_points_fastest,
951template <
int dim,
int spacedim>
952template <
typename Number>
969 this->finite_element_output.shape_hessians,
975template <
int dim,
int spacedim>
976template <
typename Number>
988 boost::container::small_vector<Number, 200> dof_values(
dofs_per_cell);
998template <
int dim,
int spacedim>
999template <
typename Number>
1004 const bool quadrature_points_fastest)
const
1016 this->finite_element_output.shape_hessians,
1018 this->finite_element_output.shape_function_to_row_table,
1020 quadrature_points_fastest);
1025template <
int dim,
int spacedim>
1026template <
typename Number>
1032 const bool quadrature_points_fastest)
const
1039 boost::container::small_vector<Number, 200> dof_values(indices.
size());
1046 this->finite_element_output.shape_function_to_row_table,
1048 quadrature_points_fastest,
1054template <
int dim,
int spacedim>
1055template <
typename Number>
1059 std::vector<Number> &laplacians)
const
1072 this->finite_element_output.shape_hessians,
1078template <
int dim,
int spacedim>
1079template <
typename Number>
1084 std::vector<Number> &laplacians)
const
1091 boost::container::small_vector<Number, 200> dof_values(
dofs_per_cell);
1101template <
int dim,
int spacedim>
1102template <
typename Number>
1118 this->finite_element_output.shape_hessians,
1120 this->finite_element_output.shape_function_to_row_table,
1126template <
int dim,
int spacedim>
1127template <
typename Number>
1141 boost::container::small_vector<Number, 200> dof_values(indices.
size());
1148 this->finite_element_output.shape_function_to_row_table,
1156template <
int dim,
int spacedim>
1157template <
typename Number>
1162 std::vector<std::vector<Number>> &laplacians,
1163 const bool quadrature_points_fastest)
const
1170 boost::container::small_vector<Number, 200> dof_values(indices.
size());
1177 this->finite_element_output.shape_function_to_row_table,
1179 quadrature_points_fastest,
1185template <
int dim,
int spacedim>
1186template <
typename Number>
1203 this->finite_element_output.shape_3rd_derivatives,
1209template <
int dim,
int spacedim>
1210template <
typename Number>
1222 boost::container::small_vector<Number, 200> dof_values(
dofs_per_cell);
1231template <
int dim,
int spacedim>
1232template <
typename Number>
1237 const bool quadrature_points_fastest)
const
1249 this->finite_element_output.shape_3rd_derivatives,
1251 this->finite_element_output.shape_function_to_row_table,
1253 quadrature_points_fastest);
1258template <
int dim,
int spacedim>
1259template <
typename Number>
1265 const bool quadrature_points_fastest)
const
1272 boost::container::small_vector<Number, 200> dof_values(indices.
size());
1279 this->finite_element_output.shape_function_to_row_table,
1281 quadrature_points_fastest,
1287template <
int dim,
int spacedim>
1296template <
int dim,
int spacedim>
1297const std::vector<Tensor<1, spacedim>> &
1302 "update_normal_vectors")));
1309template <
int dim,
int spacedim>
1330template <
int dim,
int spacedim>
1349template <
int dim,
int spacedim>
1367template <
int dim,
int spacedim>
1374 if (&cell->get_triangulation() !=
1378 ->get_triangulation())
1387 cell->get_triangulation().signals.any_change.connect(
1390 cell->get_triangulation().signals.mesh_movement.connect(
1400 cell->get_triangulation().signals.post_refinement.connect(
1403 cell->get_triangulation().signals.mesh_movement.connect(
1410template <
int dim,
int spacedim>
1431 (cell->is_translation_of(
1442 ->direction_flag() != cell->direction_flag())
1451template <
int dim,
int spacedim>
1460template <
int dim,
int spacedim>
1465template <
int dim,
int spacedim>
1471#include "fe_values_base.inst"
ArrayView< std::remove_reference_t< typename std::iterator_traits< Iterator >::reference >, MemorySpaceType > make_array_view(const Iterator begin, const Iterator end)
typename LevelSelector::cell_iterator level_cell_iterator
types::global_dof_index n_dofs_for_dof_handler() const
CellIteratorWrapper()=default
bool is_initialized() const
void get_interpolated_dof_values(const ReadVector< Number > &in, Vector< Number > &out) const
CellSimilarity::Similarity cell_similarity
::internal::FEValuesViews::Cache< dim, spacedim > fe_values_views_cache
Triangulation< dim, spacedim >::cell_iterator get_cell() const
boost::signals2::connection tria_listener_mesh_transform
void get_function_values(const ReadVector< Number > &fe_function, std::vector< Number > &values) const
FEValuesBase(const unsigned int n_q_points, const unsigned int dofs_per_cell, const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe)
internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > mapping_output
virtual ~FEValuesBase() override
const SmartPointer< const Mapping< dim, spacedim >, FEValuesBase< dim, spacedim > > mapping
const unsigned int dofs_per_cell
void check_cell_similarity(const typename Triangulation< dim, spacedim >::cell_iterator &cell)
CellIteratorWrapper present_cell
const SmartPointer< const FiniteElement< dim, spacedim >, FEValuesBase< dim, spacedim > > fe
void get_function_hessians(const ReadVector< Number > &fe_function, std::vector< Tensor< 2, spacedim, Number > > &hessians) const
void always_allow_check_for_cell_similarity(const bool allow)
void get_function_laplacians(const ReadVector< Number > &fe_function, std::vector< Number > &laplacians) const
const unsigned int n_quadrature_points
CellSimilarity::Similarity get_cell_similarity() const
const std::vector< Tensor< 1, spacedim > > & get_normal_vectors() const
std::size_t memory_consumption() const
void get_function_third_derivatives(const ReadVector< Number > &fe_function, std::vector< Tensor< 3, spacedim, Number > > &third_derivatives) const
UpdateFlags compute_update_flags(const UpdateFlags update_flags) const
std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > mapping_data
void get_function_gradients(const ReadVector< Number > &fe_function, std::vector< Tensor< 1, spacedim, Number > > &gradients) const
boost::signals2::connection tria_listener_refinement
void invalidate_present_cell()
bool check_for_cell_similarity_allowed
::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > finite_element_output
std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > fe_data
void maybe_invalidate_previous_present_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell)
const unsigned int max_n_quadrature_points
unsigned int n_dofs_per_cell() const
unsigned int n_components() const
const ComponentMask & get_nonzero_components(const unsigned int i) const
bool is_primitive() const
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const =0
std::pair< unsigned int, unsigned int > system_to_component_index(const unsigned int index) const
unsigned int n_nonzero_components(const unsigned int i) const
Abstract base class for mapping classes.
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const =0
virtual size_type size() const =0
virtual void extract_subvector_to(const ArrayView< const types::global_dof_index > &indices, ArrayView< Number > &elements) const =0
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcAccessToUninitializedField(std::string arg1)
static ::ExceptionBase & ExcNotMultiple(int arg1, int arg2)
#define Assert(cond, exc)
static ::ExceptionBase & ExcNotReinited()
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
typename ActiveSelector::cell_iterator cell_iterator
TriaIterator< CellAccessor< dim, spacedim > > cell_iterator
@ update_hessians
Second derivatives of shape functions.
@ update_values
Shape function values.
@ update_normal_vectors
Normal vectors.
@ update_3rd_derivatives
Third derivatives of shape functions.
@ update_gradients
Shape function gradients.
std::enable_if_t< std::is_fundamental_v< T >, std::size_t > memory_consumption(const T &t)
void do_function_laplacians(const ArrayView< Number2 > &dof_values, const ::Table< 2, Tensor< 2, spacedim > > &shape_hessians, std::vector< Number > &laplacians)
void do_function_derivatives(const ArrayView< Number > &dof_values, const ::Table< 2, Tensor< order, spacedim > > &shape_derivatives, std::vector< Tensor< order, spacedim, Number > > &derivatives)
std::vector< unsigned int > make_shape_function_to_row_table(const FiniteElement< dim, spacedim > &fe)
void do_function_values(const ArrayView< Number2 > &dof_values, const ::Table< 2, double > &shape_values, std::vector< Number > &values)
static const unsigned int invalid_unsigned_int
const types::global_dof_index invalid_dof_index
static constexpr DEAL_II_HOST_DEVICE_ALWAYS_INLINE const T & value(const T &t)
DEAL_II_HOST constexpr Number trace(const SymmetricTensor< 2, dim2, Number > &)