16#ifndef dealii_matrix_free_tools_h
17#define dealii_matrix_free_tools_h
41 template <
int dim,
typename AdditionalData>
44 AdditionalData & additional_data);
59 typename VectorizedArrayType,
64 VectorType & diagonal_global,
70 VectorizedArrayType> &)> &local_vmult,
71 const unsigned int dof_no = 0,
72 const unsigned int quad_no = 0,
73 const unsigned int first_selected_component = 0);
78 template <
typename CLASS,
84 typename VectorizedArrayType,
89 VectorType & diagonal_global,
95 VectorizedArrayType> &) const,
96 const CLASS * owning_class,
97 const unsigned
int dof_no = 0,
98 const unsigned
int quad_no = 0,
99 const unsigned
int first_selected_component = 0);
115 typename VectorizedArrayType,
127 VectorizedArrayType> &)> &local_vmult,
128 const unsigned int dof_no = 0,
129 const unsigned int quad_no = 0,
130 const unsigned int first_selected_component = 0);
135 template <
typename CLASS,
141 typename VectorizedArrayType,
153 VectorizedArrayType> &) const,
154 const CLASS * owning_class,
155 const unsigned
int dof_no = 0,
156 const unsigned
int quad_no = 0,
157 const unsigned
int first_selected_component = 0);
209 std::vector<unsigned int> valid_fe_indices;
211 const auto &fe_collection =
212 matrix_free.get_dof_handler(additional_data.dof_index)
213 .get_fe_collection();
215 for (
unsigned int i = 0; i < fe_collection.size(); ++i)
216 if (fe_collection[i].n_dofs_per_cell() > 0)
217 valid_fe_indices.push_back(i);
231 template <
typename VectorTypeOut,
typename VectorTypeIn>
236 const VectorTypeIn &,
237 const std::pair<unsigned int, unsigned int> &)> &cell_operation,
239 const VectorTypeIn & src,
240 const bool zero_dst_vector =
false)
const
242 const auto ebd_cell_operation = [&](
const auto &
matrix_free,
246 const auto category =
matrix_free.get_cell_range_category(range);
254 matrix_free->template cell_loop<VectorTypeOut, VectorTypeIn>(
255 ebd_cell_operation, dst, src, zero_dst_vector);
265 template <
typename VectorTypeOut,
typename VectorTypeIn>
270 const VectorTypeIn &,
271 const std::pair<unsigned int, unsigned int> &)> &cell_operation,
275 const VectorTypeIn &,
276 const std::pair<unsigned int, unsigned int> &)> &face_operation,
280 const VectorTypeIn &,
281 const std::pair<unsigned int, unsigned int> &,
282 const bool)> &boundary_operation,
284 const VectorTypeIn & src,
285 const bool zero_dst_vector =
false)
const
287 const auto ebd_cell_operation = [&](
const auto &
matrix_free,
291 const auto category =
matrix_free.get_cell_range_category(range);
299 const auto ebd_internal_or_boundary_face_operation =
304 const auto category =
matrix_free.get_face_range_category(range);
306 const unsigned int type =
320 matrix_free->template loop<VectorTypeOut, VectorTypeIn>(
322 ebd_internal_or_boundary_face_operation,
323 ebd_internal_or_boundary_face_operation,
346 template <
int dim,
typename AdditionalData>
349 AdditionalData & additional_data)
352 const unsigned int level = additional_data.mg_level;
357 additional_data.cell_vectorization_category.assign(
367 std::sort(bids.begin(), bids.end());
370 unsigned int n_bids = bids.size() + 1;
372 for (
unsigned int i = 0; i < GeometryInfo<dim>::faces_per_cell;
373 i++, offset = offset * n_bids)
377 const auto to_category = [&](
const auto &cell) {
378 unsigned int c_num = 0;
379 for (
unsigned int i = 0; i < GeometryInfo<dim>::faces_per_cell; ++i)
381 const auto face = cell->face(i);
382 if (face->at_boundary() && !cell->has_periodic_neighbor(i))
384 factors[i] * (1 + std::distance(bids.begin(),
385 std::find(bids.begin(),
387 face->boundary_id())));
396 if (cell->is_locally_owned())
398 .cell_vectorization_category[cell->active_cell_index()] =
406 if (cell->is_locally_owned_on_level())
407 additional_data.cell_vectorization_category[cell->index()] =
413 additional_data.hold_all_faces_to_owned_cells =
true;
414 additional_data.cell_vectorization_categories_strict =
true;
415 additional_data.mapping_update_flags_faces_by_cells =
416 additional_data.mapping_update_flags_inner_faces |
417 additional_data.mapping_update_flags_boundary_faces;
422 template <
typename Number>
425 std::vector<unsigned int> row_lid_to_gid;
426 std::vector<unsigned int> row;
427 std::vector<unsigned int> col;
428 std::vector<Number> val;
430 std::vector<unsigned int> inverse_lookup_rows;
431 std::vector<std::pair<unsigned int, unsigned int>> inverse_lookup_origins;
439 typename VectorizedArrayType>
440 class ComputeDiagonalHelper
443 static const unsigned int n_lanes = VectorizedArrayType::size();
445 ComputeDiagonalHelper()
447 , dofs_per_component(0)
450 ComputeDiagonalHelper(
const ComputeDiagonalHelper &)
452 , dofs_per_component(0)
461 VectorizedArrayType> &phi)
465 if (dofs_per_component != phi.dofs_per_component)
467 locally_relevant_constraints_hn_map.clear();
468 dofs_per_component = phi.dofs_per_component;
474 reinit(
const unsigned int cell)
476 this->phi->reinit(cell);
479 const auto & dof_info = phi->get_dof_info();
480 const unsigned int n_fe_components = dof_info.start_components.back();
481 const auto & matrix_free = phi->get_matrix_free();
487 const unsigned int first_selected_component =
488 n_fe_components == 1 ? 0 : phi->get_first_selected_component();
490 const unsigned int n_lanes_filled =
498 std::vector<std::tuple<unsigned int, unsigned int, Number>>
499 locally_relevant_constraints, locally_relevant_constraints_tmp;
500 locally_relevant_constraints.reserve(phi->dofs_per_cell);
501 std::vector<unsigned int> constraint_position;
502 std::vector<unsigned char> is_constrained_hn;
504 for (
unsigned int v = 0; v < n_lanes_filled; ++v)
506 const unsigned int *dof_indices;
507 unsigned int index_indicators, next_index_indicators;
509 const unsigned int start =
510 (cell * n_lanes + v) * n_fe_components + first_selected_component;
512 dof_info.dof_indices.data() + dof_info.row_starts[start].first;
513 index_indicators = dof_info.row_starts[start].second;
514 next_index_indicators = dof_info.row_starts[start + 1].second;
518 locally_relevant_constraints.clear();
520 if (n_components == 1 || n_fe_components == 1)
522 unsigned int ind_local = 0;
523 for (; index_indicators != next_index_indicators;
524 ++index_indicators, ++ind_local)
526 const std::pair<unsigned short, unsigned short> indicator =
527 dof_info.constraint_indicator[index_indicators];
529 for (
unsigned int j = 0; j < indicator.first;
531 locally_relevant_constraints.emplace_back(ind_local,
535 dof_indices += indicator.first;
537 const Number *data_val =
539 const Number *end_pool =
542 for (; data_val != end_pool; ++data_val, ++dof_indices)
543 locally_relevant_constraints.emplace_back(ind_local,
550 for (; ind_local < dofs_per_component;
551 ++dof_indices, ++ind_local)
552 locally_relevant_constraints.emplace_back(ind_local,
563 for (
unsigned int comp = 0; comp < n_components; ++comp)
565 unsigned int ind_local = 0;
569 for (; index_indicators != next_index_indicators;
570 ++index_indicators, ++ind_local)
572 const std::pair<unsigned short, unsigned short>
574 dof_info.constraint_indicator[index_indicators];
577 for (
unsigned int j = 0; j < indicator.first;
579 locally_relevant_constraints.emplace_back(
580 comp * dofs_per_component + ind_local,
583 dof_indices += indicator.first;
585 const Number *data_val =
587 const Number *end_pool =
590 for (; data_val != end_pool; ++data_val, ++dof_indices)
591 locally_relevant_constraints.emplace_back(
592 comp * dofs_per_component + ind_local,
600 for (; ind_local < dofs_per_component;
601 ++dof_indices, ++ind_local)
602 locally_relevant_constraints.emplace_back(
603 comp * dofs_per_component + ind_local,
607 if (comp + 1 < n_components)
608 next_index_indicators =
609 dof_info.row_starts[start + comp + 2].second;
616 for (
unsigned int i = 1; i < locally_relevant_constraints.size();
618 Assert(std::get<0>(locally_relevant_constraints[i]) >=
619 std::get<0>(locally_relevant_constraints[i - 1]),
623 if (dof_info.hanging_node_constraint_masks.size() > 0 &&
624 dof_info.hanging_node_constraint_masks_comp.size() > 0 &&
625 dof_info.hanging_node_constraint_masks_comp
626 [phi->get_active_fe_index()][first_selected_component])
629 dof_info.hanging_node_constraint_masks[cell * n_lanes + v];
632 if (mask != ::internal::MatrixFreeFunctions::
637 if (locally_relevant_constraints_hn_map.find(mask) ==
638 locally_relevant_constraints_hn_map.end())
639 fill_constraint_type_into_map(mask);
641 const auto &locally_relevant_constraints_hn =
642 locally_relevant_constraints_hn_map[
mask];
644 locally_relevant_constraints_tmp.clear();
645 if (locally_relevant_constraints_tmp.capacity() <
646 locally_relevant_constraints.size())
647 locally_relevant_constraints_tmp.reserve(
648 locally_relevant_constraints.size() +
649 locally_relevant_constraints_hn.size());
654 constraint_position.assign(phi->dofs_per_cell,
656 for (
auto &a : locally_relevant_constraints)
657 if (constraint_position[
std::get<0>(a)] ==
659 constraint_position[
std::get<0>(a)] =
660 std::distance(locally_relevant_constraints.data(),
662 is_constrained_hn.assign(phi->dofs_per_cell,
false);
663 for (
auto &hn : locally_relevant_constraints_hn)
664 is_constrained_hn[
std::get<0>(hn)] = 1;
667 for (
const auto &a : locally_relevant_constraints)
668 if (is_constrained_hn[
std::get<0>(a)] == 0)
669 locally_relevant_constraints_tmp.push_back(a);
673 for (
const auto &hn : locally_relevant_constraints_hn)
674 if (constraint_position[
std::get<1>(hn)] !=
678 locally_relevant_constraints.size());
679 auto other = locally_relevant_constraints.begin() +
680 constraint_position[std::get<1>(hn)];
683 for (; other != locally_relevant_constraints.end() &&
684 std::get<0>(*other) == std::get<1>(hn);
686 locally_relevant_constraints_tmp.emplace_back(
689 std::get<2>(hn) * std::get<2>(*other));
692 std::swap(locally_relevant_constraints,
693 locally_relevant_constraints_tmp);
698 std::sort(locally_relevant_constraints.begin(),
699 locally_relevant_constraints.end(),
700 [](
const auto &a,
const auto &b) {
701 if (std::get<1>(a) < std::get<1>(b))
703 return (std::get<1>(a) == std::get<1>(b)) &&
704 (std::get<0>(a) < std::get<0>(b));
708 auto &c_pool = c_pools[v];
710 c_pool.row_lid_to_gid.clear();
712 c_pool.row.push_back(0);
716 if (locally_relevant_constraints.size() > 0)
717 c_pool.row_lid_to_gid.emplace_back(
718 std::get<1>(locally_relevant_constraints.front()));
719 for (
const auto &j : locally_relevant_constraints)
721 if (c_pool.row_lid_to_gid.back() != std::get<1>(j))
723 c_pool.row_lid_to_gid.push_back(std::get<1>(j));
724 c_pool.row.push_back(c_pool.val.size());
727 c_pool.col.emplace_back(std::get<0>(j));
728 c_pool.val.emplace_back(std::get<2>(j));
731 if (c_pool.val.size() > 0)
732 c_pool.row.push_back(c_pool.val.size());
734 c_pool.inverse_lookup_rows.clear();
735 c_pool.inverse_lookup_rows.resize(1 + phi->dofs_per_cell);
736 for (
const unsigned int i : c_pool.col)
737 c_pool.inverse_lookup_rows[1 + i]++;
739 std::partial_sum(c_pool.inverse_lookup_rows.begin(),
740 c_pool.inverse_lookup_rows.end(),
741 c_pool.inverse_lookup_rows.begin());
745 c_pool.inverse_lookup_origins.resize(c_pool.col.size());
746 std::vector<unsigned int> inverse_lookup_count(
748 for (
unsigned int row = 0; row < c_pool.row.size() - 1; ++row)
749 for (
unsigned int col = c_pool.row[row];
750 col < c_pool.row[row + 1];
753 const unsigned int index = c_pool.col[col];
754 c_pool.inverse_lookup_origins
755 [c_pool.inverse_lookup_rows[
index] +
756 inverse_lookup_count[
index]] = std::make_pair(row, col);
757 ++inverse_lookup_count[
index];
786 for (
unsigned int v = 0; v < n_lanes_filled; ++v)
787 diagonals_local_constrained[v].assign(
788 c_pools[v].row_lid_to_gid.size() *
789 (n_fe_components == 1 ? n_components : 1),
794 fill_constraint_type_into_map(
795 const ::internal::MatrixFreeFunctions::compressed_constraint_kind
798 auto &constraints_hn = locally_relevant_constraints_hn_map[
mask];
803 const unsigned int degree =
804 phi->get_shape_info().data.front().fe_degree;
809 values_dofs.resize(dofs_per_component);
812 VectorizedArrayType::size()>
814 constraint_mask.fill(::internal::MatrixFreeFunctions::
816 constraint_mask[0] =
mask;
818 for (
unsigned int i = 0; i < dofs_per_component; ++i)
820 for (
unsigned int j = 0; j < dofs_per_component; ++j)
821 values_dofs[j] = VectorizedArrayType();
822 values_dofs[i] = Number(1);
827 VectorizedArrayType>::apply(1,
829 phi->get_shape_info(),
834 const Number tolerance =
836 std::numeric_limits<Number>::epsilon() * 16);
837 for (
unsigned int j = 0; j < dofs_per_component; ++j)
838 if (
std::abs(values_dofs[j][0]) > tolerance &&
840 std::abs(values_dofs[j][0] - Number(1)) > tolerance))
841 constraints_hn.emplace_back(j, i, values_dofs[j][0]);
845 const unsigned int n_hn_constraints = constraints_hn.size();
846 constraints_hn.resize(n_hn_constraints * n_components);
848 for (
unsigned int c = 1; c < n_components; ++c)
849 for (
unsigned int i = 0; i < n_hn_constraints; ++i)
850 constraints_hn[c * n_hn_constraints + i] = std::make_tuple(
851 std::get<0>(constraints_hn[i]) + c * dofs_per_component,
852 std::get<1>(constraints_hn[i]) + c * dofs_per_component,
853 std::get<2>(constraints_hn[i]));
857 prepare_basis_vector(
const unsigned int i)
864 for (
unsigned int j = 0; j < phi->dofs_per_cell; ++j)
865 phi->begin_dof_values()[j] = VectorizedArrayType();
866 phi->begin_dof_values()[i] = Number(1);
875 const unsigned int n_fe_components =
876 phi->get_dof_info().start_components.back();
877 const unsigned int comp =
878 n_fe_components == 1 ? i / dofs_per_component : 0;
879 const unsigned int i_comp =
880 n_fe_components == 1 ? (i % dofs_per_component) : i;
884 for (
unsigned int v = 0;
885 v < phi->get_matrix_free().n_active_entries_per_cell_batch(
886 phi->get_current_cell_index());
889 const auto &c_pool = c_pools[v];
891 for (
unsigned int jj = c_pool.inverse_lookup_rows[i_comp];
892 jj < c_pool.inverse_lookup_rows[i_comp + 1];
895 const unsigned int j = c_pool.inverse_lookup_origins[jj].first;
898 for (
unsigned int k = c_pool.row[j]; k < c_pool.row[j + 1]; ++k)
899 temp += c_pool.val[k] *
900 phi->begin_dof_values()[comp * dofs_per_component +
904 diagonals_local_constrained
905 [v][j + comp * c_pools[v].row_lid_to_gid.size()] +=
906 temp * c_pool.val[c_pool.inverse_lookup_origins[jj].second];
911 template <
typename VectorType>
913 distribute_local_to_global(
914 std::array<VectorType *, n_components> &diagonal_global)
917 const unsigned int n_fe_components =
918 phi->get_dof_info().start_components.back();
920 for (
unsigned int v = 0;
921 v < phi->get_matrix_free().n_active_entries_per_cell_batch(
922 phi->get_current_cell_index());
927 for (
unsigned int j = 0; j < c_pools[v].row.size() - 1; ++j)
928 for (
unsigned int comp = 0;
929 comp < (n_fe_components == 1 ?
930 static_cast<unsigned int>(n_components) :
934 *diagonal_global[n_fe_components == 1 ? comp : 0],
935 c_pools[v].row_lid_to_gid[j],
936 diagonals_local_constrained
937 [v][j + comp * c_pools[v].row_lid_to_gid.size()]);
946 VectorizedArrayType> *phi;
948 unsigned int dofs_per_component;
952 std::array<internal::LocalCSR<Number>, n_lanes> c_pools;
957 std::array<std::vector<Number>, n_lanes> diagonals_local_constrained;
961 std::vector<std::tuple<unsigned int, unsigned int, Number>>>
962 locally_relevant_constraints_hn_map;
975 typename VectorizedArrayType,
980 VectorType & diagonal_global,
986 VectorizedArrayType> &)> &local_vmult,
987 const unsigned int dof_no,
988 const unsigned int quad_no,
989 const unsigned int first_selected_component)
993 std::array<typename ::internal::BlockVectorSelector<
997 diagonal_global_components;
999 for (
unsigned int d = 0;
d < n_components; ++
d)
1000 diagonal_global_components[d] = ::internal::
1002 get_vector_component(diagonal_global,
d + first_selected_component);
1004 const auto &dof_info = matrix_free.
get_dof_info(dof_no);
1006 if (dof_info.start_components.back() == 1)
1007 for (
unsigned int comp = 0; comp < n_components; ++comp)
1009 Assert(diagonal_global_components[comp] !=
nullptr,
1010 ExcMessage(
"The finite element underlying this FEEvaluation "
1011 "object is scalar, but you requested " +
1012 std::to_string(n_components) +
1013 " components via the template argument in "
1014 "FEEvaluation. In that case, you must pass an "
1015 "std::vector<VectorType> or a BlockVector to " +
1016 "read_dof_values and distribute_local_to_global."));
1018 *diagonal_global_components[comp], matrix_free, dof_info);
1023 *diagonal_global_components[0], matrix_free, dof_info);
1026 using Helper = internal::ComputeDiagonalHelper<dim,
1031 VectorizedArrayType>;
1034 matrix_free.template cell_loop<VectorType, int>(
1038 const std::pair<unsigned int, unsigned int> &range)
mutable {
1039 Helper &helper = scratch_data.
get();
1046 VectorizedArrayType>
1047 phi(matrix_free, range, dof_no, quad_no, first_selected_component);
1048 helper.initialize(phi);
1050 for (
unsigned int cell = range.first; cell < range.second; ++cell)
1052 helper.reinit(cell);
1054 for (
unsigned int i = 0; i < phi.dofs_per_cell; ++i)
1056 helper.prepare_basis_vector(i);
1061 helper.distribute_local_to_global(diagonal_global_components);
1069 template <
typename CLASS,
1075 typename VectorizedArrayType,
1076 typename VectorType>
1080 VectorType & diagonal_global,
1086 VectorizedArrayType> &) const,
1087 const CLASS * owning_class,
1088 const unsigned
int dof_no,
1089 const unsigned
int quad_no,
1090 const unsigned
int first_selected_component)
1097 VectorizedArrayType,
1101 [&](
auto &feeval) { (owning_class->*cell_operation)(feeval); },
1104 first_selected_component);
1113 template <
typename MatrixType,
1115 std::enable_if_t<std::is_same<
1116 typename std::remove_const<
typename std::remove_reference<
1117 typename MatrixType::value_type>::type>::type,
1118 typename std::remove_const<
typename std::remove_reference<
1119 Number>::type>::type>::value> * =
nullptr>
1121 create_new_affine_constraints_if_needed(
1134 template <
typename MatrixType,
1136 std::enable_if_t<!std::is_same<
1137 typename std::remove_const<
typename std::remove_reference<
1138 typename MatrixType::value_type>::type>::type,
1139 typename std::remove_const<
typename std::remove_reference<
1140 Number>::type>::type>::value> * =
nullptr>
1142 create_new_affine_constraints_if_needed(
1149 std::make_unique<AffineConstraints<typename MatrixType::value_type>>();
1150 new_constraints->
copy_from(constraints);
1152 return *new_constraints;
1161 typename VectorizedArrayType,
1162 typename MatrixType>
1167 MatrixType & matrix,
1173 VectorizedArrayType> &)> &local_vmult,
1174 const unsigned int dof_no,
1175 const unsigned int quad_no,
1176 const unsigned int first_selected_component)
1178 std::unique_ptr<AffineConstraints<typename MatrixType::value_type>>
1179 constraints_for_matrix;
1181 internal::create_new_affine_constraints_if_needed(matrix,
1183 constraints_for_matrix);
1185 matrix_free.template cell_loop<MatrixType, MatrixType>(
1186 [&](
const auto &,
auto &dst,
const auto &,
const auto range) {
1192 VectorizedArrayType>
1194 matrix_free, range, dof_no, quad_no, first_selected_component);
1196 const unsigned int dofs_per_cell = integrator.dofs_per_cell;
1198 std::vector<types::global_dof_index> dof_indices(dofs_per_cell);
1199 std::vector<types::global_dof_index> dof_indices_mf(dofs_per_cell);
1201 std::array<FullMatrix<typename MatrixType::value_type>,
1202 VectorizedArrayType::size()>
1205 std::fill_n(matrices.begin(),
1206 VectorizedArrayType::size(),
1210 const auto lexicographic_numbering =
1214 first_selected_component,
1215 integrator.get_active_fe_index(),
1216 integrator.get_active_quadrature_index())
1217 .lexicographic_numbering;
1219 for (
auto cell = range.first; cell < range.second; ++cell)
1221 integrator.reinit(cell);
1223 const unsigned int n_filled_lanes =
1226 for (
unsigned int v = 0; v < n_filled_lanes; ++v)
1229 for (
unsigned int j = 0; j < dofs_per_cell; ++j)
1231 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
1232 integrator.begin_dof_values()[i] =
1233 static_cast<Number
>(i == j);
1235 local_vmult(integrator);
1237 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
1238 for (
unsigned int v = 0; v < n_filled_lanes; ++v)
1239 matrices[v](i, j) = integrator.begin_dof_values()[i][v];
1242 for (
unsigned int v = 0; v < n_filled_lanes; ++v)
1248 cell_v->get_mg_dof_indices(dof_indices);
1250 cell_v->get_dof_indices(dof_indices);
1252 for (
unsigned int j = 0; j < dof_indices.size(); ++j)
1253 dof_indices_mf[j] = dof_indices[lexicographic_numbering[j]];
1267 template <
typename CLASS,
1273 typename VectorizedArrayType,
1274 typename MatrixType>
1279 MatrixType & matrix,
1285 VectorizedArrayType> &) const,
1286 const CLASS * owning_class,
1287 const unsigned
int dof_no,
1288 const unsigned
int quad_no,
1289 const unsigned
int first_selected_component)
1296 VectorizedArrayType,
1301 [&](
auto &feeval) { (owning_class->*cell_operation)(feeval); },
1304 first_selected_component);
void distribute_local_to_global(const InVector &local_vector, const std::vector< size_type > &local_dof_indices, OutVector &global_vector) const
void copy_from(const AffineConstraints< other_number > &other)
const internal::MatrixFreeFunctions::ShapeInfo< VectorizedArrayType > & get_shape_info(const unsigned int dof_handler_index_component=0, const unsigned int quad_index=0, const unsigned int fe_base_element=0, const unsigned int hp_active_fe_index=0, const unsigned int hp_active_quad_index=0) const
unsigned int n_active_entries_per_cell_batch(const unsigned int cell_batch_index) const
unsigned int get_mg_level() const
const internal::MatrixFreeFunctions::DoFInfo & get_dof_info(const unsigned int dof_handler_index_component=0) const
const Number * constraint_pool_begin(const unsigned int pool_index) const
DoFHandler< dim >::cell_iterator get_cell_iterator(const unsigned int cell_batch_index, const unsigned int lane_index, const unsigned int dof_handler_index=0) const
const Number * constraint_pool_end(const unsigned int pool_index) const
A class that provides a separate storage location on each thread that accesses the object.
cell_iterator begin(const unsigned int level=0) const
unsigned int n_active_cells() const
cell_iterator end() const
virtual std::vector< types::boundary_id > get_boundary_ids() const
active_cell_iterator begin_active(const unsigned int level=0) const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
@ matrix
Contents is actually a matrix.
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
constexpr T pow(const T base, const int iexp)
std::uint8_t compressed_constraint_kind
constexpr compressed_constraint_kind unconstrained_compressed_constraint_kind
void vector_access_add(VectorType &vec, const unsigned int entry, const typename VectorType::value_type &val)
void check_vector_compatibility(const VectorType &vec, const MatrixFree< dim, Number, VectorizedArrayType > &matrix_free, const internal::MatrixFreeFunctions::DoFInfo &dof_info)
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
static const unsigned int invalid_unsigned_int
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
const ::Triangulation< dim, spacedim > & tria