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);
164 template <
int dim,
typename AdditionalData>
167 AdditionalData & additional_data)
170 const unsigned int level = additional_data.mg_level;
175 additional_data.cell_vectorization_category.assign(
185 std::sort(bids.begin(), bids.end());
188 unsigned int n_bids = bids.size() + 1;
190 for (
unsigned int i = 0; i < GeometryInfo<dim>::faces_per_cell;
191 i++, offset = offset * n_bids)
195 const auto to_category = [&](
const auto &cell) {
196 unsigned int c_num = 0;
197 for (
unsigned int i = 0; i < GeometryInfo<dim>::faces_per_cell; ++i)
199 auto &face = *cell->face(i);
200 if (face.at_boundary() && !cell->has_periodic_neighbor(i))
202 factors[i] * (1 + std::distance(bids.begin(),
203 std::find(bids.begin(),
205 face.boundary_id())));
214 if (cell->is_locally_owned())
216 .cell_vectorization_category[cell->active_cell_index()] =
224 if (cell->is_locally_owned_on_level())
225 additional_data.cell_vectorization_category[cell->index()] =
231 additional_data.hold_all_faces_to_owned_cells =
true;
232 additional_data.cell_vectorization_categories_strict =
true;
233 additional_data.mapping_update_flags_faces_by_cells =
234 additional_data.mapping_update_flags_inner_faces |
235 additional_data.mapping_update_flags_boundary_faces;
240 template <
typename Number>
247 std::vector<unsigned int> row_lid_to_gid;
248 std::vector<unsigned int> row;
249 std::vector<unsigned int> col;
250 std::vector<Number> val;
258 typename VectorizedArrayType>
259 class ComputeDiagonalHelper
262 static const unsigned int n_lanes = VectorizedArrayType::size();
269 VectorizedArrayType> &phi)
274 reinit(
const unsigned int cell)
276 this->phi.reinit(cell);
278 const auto & dof_info = phi.get_dof_info();
279 const unsigned int n_fe_components = dof_info.start_components.back();
280 const unsigned int dofs_per_component = phi.dofs_per_component;
281 const auto & matrix_free = phi.get_matrix_free();
287 const unsigned int first_selected_component =
288 n_fe_components == 1 ? 0 : phi.get_first_selected_component();
290 const unsigned int n_lanes_filled =
293 std::array<const unsigned int *, n_lanes> dof_indices{};
295 for (
unsigned int v = 0; v < n_lanes_filled; ++v)
297 dof_info.dof_indices.data() +
299 .row_starts[(cell * n_lanes + v) * n_fe_components +
300 first_selected_component]
306 c_pools = std::array<internal::LocalCSR<Number>, n_lanes>();
308 for (
unsigned int v = 0; v < n_lanes_filled; ++v)
310 unsigned int index_indicators, next_index_indicators;
314 .row_starts[(cell * n_lanes + v) * n_fe_components +
315 first_selected_component]
317 next_index_indicators =
319 .row_starts[(cell * n_lanes + v) * n_fe_components +
320 first_selected_component + 1]
325 std::vector<std::tuple<unsigned int, unsigned int, Number>>
326 locally_relevant_constrains;
330 if (n_components == 1 || n_fe_components == 1)
332 unsigned int ind_local = 0;
333 for (; index_indicators != next_index_indicators;
334 ++index_indicators, ++ind_local)
336 const std::pair<unsigned short, unsigned short> indicator =
337 dof_info.constraint_indicator[index_indicators];
339 for (
unsigned int j = 0; j < indicator.first;
341 locally_relevant_constrains.emplace_back(
342 ind_local, dof_indices[v][j], 1.0);
344 dof_indices[v] += indicator.first;
346 const Number *data_val =
348 const Number *end_pool =
351 for (; data_val != end_pool; ++data_val, ++dof_indices[v])
352 locally_relevant_constrains.emplace_back(ind_local,
359 for (; ind_local < dofs_per_component;
360 ++dof_indices[v], ++ind_local)
361 locally_relevant_constrains.emplace_back(ind_local,
372 for (
unsigned int comp = 0; comp < n_components; ++comp)
374 unsigned int ind_local = 0;
378 for (; index_indicators != next_index_indicators;
379 ++index_indicators, ++ind_local)
381 const std::pair<unsigned short, unsigned short>
383 dof_info.constraint_indicator[index_indicators];
386 for (
unsigned int j = 0; j < indicator.first;
388 locally_relevant_constrains.emplace_back(
389 comp * dofs_per_component + ind_local,
392 dof_indices[v] += indicator.first;
394 const Number *data_val =
396 const Number *end_pool =
399 for (; data_val != end_pool;
400 ++data_val, ++dof_indices[v])
401 locally_relevant_constrains.emplace_back(
402 comp * dofs_per_component + ind_local,
410 for (; ind_local < dofs_per_component;
411 ++dof_indices[v], ++ind_local)
412 locally_relevant_constrains.emplace_back(
413 comp * dofs_per_component + ind_local,
417 if (comp + 1 < n_components)
419 next_index_indicators =
421 .row_starts[(cell * n_lanes + v) * n_fe_components +
422 first_selected_component + comp + 2]
430 std::sort(locally_relevant_constrains.begin(),
431 locally_relevant_constrains.end(),
432 [](
const auto &a,
const auto &b) {
433 if (std::get<0>(a) < std::get<0>(b))
435 return (std::get<0>(a) == std::get<0>(b)) &&
436 (std::get<1>(a) < std::get<1>(b));
440 locally_relevant_constrains.erase(
441 unique(locally_relevant_constrains.begin(),
442 locally_relevant_constrains.end(),
443 [](
const auto &a,
const auto &b) {
444 return (std::get<1>(a) == std::get<1>(b)) &&
445 (std::get<0>(a) == std::get<0>(b));
447 locally_relevant_constrains.end());
450 if (dof_info.hanging_node_constraint_masks.size() > 0 &&
451 dof_info.hanging_node_constraint_masks_comp.size() > 0 &&
453 .hanging_node_constraint_masks_comp[phi.get_active_fe_index()]
454 [first_selected_component])
457 dof_info.hanging_node_constraint_masks[cell * n_lanes + v];
460 if (mask != ::internal::MatrixFreeFunctions::
461 unconstrained_compressed_constraint_kind)
465 if (locally_relevant_constrains_hn_map.find(mask) ==
466 locally_relevant_constrains_hn_map.end())
475 VectorizedArrayType::size()>
477 constraint_mask.fill(
478 ::internal::MatrixFreeFunctions::
479 unconstrained_compressed_constraint_kind);
481 constraint_mask[0] =
mask;
484 std::tuple<unsigned int, unsigned int, Number>>
485 locally_relevant_constrains_hn;
487 for (
unsigned int i = 0; i < dofs_per_component; ++i)
489 for (
unsigned int j = 0; j < dofs_per_component;
491 values_dofs[j][0] =
static_cast<Number
>(i == j);
496 VectorizedArrayType>::apply(1,
500 phi.get_shape_info(),
505 for (
unsigned int j = 0; j < dofs_per_component;
507 if (1e-10 <
std::abs(values_dofs[j][0]) &&
509 1e-10 <
std::abs(values_dofs[j][0] - 1.0)))
510 locally_relevant_constrains_hn.emplace_back(
511 j, i, values_dofs[j][0]);
515 std::sort(locally_relevant_constrains_hn.begin(),
516 locally_relevant_constrains_hn.end(),
517 [](
const auto &a,
const auto &b) {
518 if (std::get<0>(a) < std::get<0>(b))
520 return (std::get<0>(a) == std::get<0>(b)) &&
521 (std::get<1>(a) < std::get<1>(b));
525 const unsigned int n_hn_constraints =
526 locally_relevant_constrains_hn.size();
527 locally_relevant_constrains_hn.resize(n_hn_constraints *
530 for (
unsigned int c = 0; c < n_components; ++c)
531 for (
unsigned int i = 0; i < n_hn_constraints; ++i)
532 locally_relevant_constrains_hn[c *
535 std::tuple<unsigned int, unsigned int, Number>{
536 std::get<0>(locally_relevant_constrains_hn[i]) +
537 c * dofs_per_component,
538 std::get<1>(locally_relevant_constrains_hn[i]) +
539 c * dofs_per_component,
540 std::get<2>(locally_relevant_constrains_hn[i])};
543 locally_relevant_constrains_hn_map[
mask] =
544 locally_relevant_constrains_hn;
547 const auto &locally_relevant_constrains_hn =
548 locally_relevant_constrains_hn_map[
mask];
551 std::vector<std::tuple<unsigned int, unsigned int, Number>>
552 locally_relevant_constrains_temp;
554 for (
unsigned int i = 0;
555 i < dofs_per_component * n_components;
558 const auto lower_bound_fu = [](
const auto &a,
560 return std::get<0>(a) <
b;
563 const auto upper_bound_fu = [](
const auto &a,
565 return a < std::get<0>(b);
568 const auto i_begin = std::lower_bound(
569 locally_relevant_constrains_hn.begin(),
570 locally_relevant_constrains_hn.end(),
573 const auto i_end = std::upper_bound(
574 locally_relevant_constrains_hn.begin(),
575 locally_relevant_constrains_hn.end(),
579 if (i_begin == i_end)
583 const auto j_begin = std::lower_bound(
584 locally_relevant_constrains.begin(),
585 locally_relevant_constrains.end(),
588 const auto j_end = std::upper_bound(
589 locally_relevant_constrains.begin(),
590 locally_relevant_constrains.end(),
594 for (
auto v = j_begin; v != j_end; ++v)
595 locally_relevant_constrains_temp.emplace_back(*v);
600 for (
auto v0 = i_begin;
v0 != i_end; ++
v0)
602 const auto j_begin = std::lower_bound(
603 locally_relevant_constrains.begin(),
604 locally_relevant_constrains.end(),
607 const auto j_end = std::upper_bound(
608 locally_relevant_constrains.begin(),
609 locally_relevant_constrains.end(),
613 for (
auto v1 = j_begin;
v1 != j_end; ++
v1)
614 locally_relevant_constrains_temp.emplace_back(
617 std::get<2>(*
v0) * std::get<2>(*
v1));
622 locally_relevant_constrains =
623 locally_relevant_constrains_temp;
628 std::sort(locally_relevant_constrains.begin(),
629 locally_relevant_constrains.end(),
630 [](
const auto &a,
const auto &b) {
631 if (std::get<1>(a) < std::get<1>(b))
633 return (std::get<1>(a) == std::get<1>(b)) &&
634 (std::get<0>(a) < std::get<0>(b));
638 auto &c_pool = c_pools[v];
640 if (locally_relevant_constrains.size() > 0)
641 c_pool.row_lid_to_gid.emplace_back(
642 std::get<1>(locally_relevant_constrains.front()));
643 for (
const auto &j : locally_relevant_constrains)
645 if (c_pool.row_lid_to_gid.back() != std::get<1>(j))
647 c_pool.row_lid_to_gid.push_back(std::get<1>(j));
648 c_pool.row.push_back(c_pool.val.size());
651 c_pool.col.emplace_back(std::get<0>(j));
652 c_pool.val.emplace_back(std::get<2>(j));
655 if (c_pool.val.size() > 0)
656 c_pool.row.push_back(c_pool.val.size());
683 for (
unsigned int v = 0; v < n_lanes_filled; ++v)
684 diagonals_local_constrained[v].assign(
685 c_pools[v].row_lid_to_gid.size() *
686 (n_fe_components == 1 ? n_components : 1),
691 prepare_basis_vector(
const unsigned int i)
698 for (
unsigned int j = 0; j < phi.dofs_per_cell; ++j)
699 phi.begin_dof_values()[j] =
static_cast<Number
>(i == j);
708 const unsigned int n_fe_components =
709 phi.get_dof_info().start_components.back();
710 const unsigned int comp =
711 n_fe_components == 1 ? i / phi.dofs_per_component : 0;
712 const unsigned int i_comp =
713 n_fe_components == 1 ? (i % phi.dofs_per_component) : i;
717 for (
unsigned int v = 0;
718 v < phi.get_matrix_free().n_active_entries_per_cell_batch(
719 phi.get_current_cell_index());
722 const auto &c_pool = c_pools[v];
724 for (
unsigned int j = 0; j < c_pool.row.size() - 1; ++j)
728 const auto scale_iterator =
729 std::lower_bound(c_pool.col.begin() + c_pool.row[j],
730 c_pool.col.begin() + c_pool.row[j + 1],
734 if (scale_iterator == c_pool.col.begin() + c_pool.row[j + 1])
738 if (*scale_iterator != i_comp)
743 for (
unsigned int k = c_pool.row[j]; k < c_pool.row[j + 1]; ++k)
744 temp += c_pool.val[k] *
745 phi.begin_dof_values()[comp * phi.dofs_per_component +
749 diagonals_local_constrained
750 [v][j + comp * c_pools[v].row_lid_to_gid.size()] +=
752 c_pool.val[std::distance(c_pool.col.begin(), scale_iterator)];
757 template <
typename VectorType>
759 distribute_local_to_global(
760 std::array<VectorType *, n_components> &diagonal_global)
763 const unsigned int n_fe_components =
764 phi.get_dof_info().start_components.back();
766 for (
unsigned int v = 0;
767 v < phi.get_matrix_free().n_active_entries_per_cell_batch(
768 phi.get_current_cell_index());
773 for (
unsigned int j = 0; j < c_pools[v].row.size() - 1; ++j)
774 for (
unsigned int comp = 0;
775 comp < (n_fe_components == 1 ?
776 static_cast<unsigned int>(n_components) :
780 *diagonal_global[n_fe_components == 1 ? comp : 0],
781 c_pools[v].row_lid_to_gid[j],
782 diagonals_local_constrained
783 [v][j + comp * c_pools[v].row_lid_to_gid.size()]);
792 VectorizedArrayType> φ
796 std::array<internal::LocalCSR<Number>, n_lanes> c_pools;
801 std::array<std::vector<Number>, n_lanes> diagonals_local_constrained;
805 std::vector<std::tuple<unsigned int, unsigned int, Number>>>
806 locally_relevant_constrains_hn_map;
816 typename VectorizedArrayType,
821 VectorType & diagonal_global,
827 VectorizedArrayType> &)> &local_vmult,
828 const unsigned int dof_no,
829 const unsigned int quad_no,
830 const unsigned int first_selected_component)
834 std::array<typename ::internal::BlockVectorSelector<
838 diagonal_global_components;
840 for (
unsigned int d = 0;
d < n_components; ++
d)
841 diagonal_global_components[d] = ::internal::
843 get_vector_component(diagonal_global, d + first_selected_component);
845 const auto &dof_info = matrix_free.
get_dof_info(dof_no);
847 if (dof_info.start_components.back() == 1)
848 for (
unsigned int comp = 0; comp < n_components; ++comp)
850 Assert(diagonal_global_components[comp] !=
nullptr,
851 ExcMessage(
"The finite element underlying this FEEvaluation "
852 "object is scalar, but you requested " +
853 std::to_string(n_components) +
854 " components via the template argument in "
855 "FEEvaluation. In that case, you must pass an "
856 "std::vector<VectorType> or a BlockVector to " +
857 "read_dof_values and distribute_local_to_global."));
859 *diagonal_global_components[comp], matrix_free, dof_info);
864 *diagonal_global_components[0], matrix_free, dof_info);
867 matrix_free.template cell_loop<VectorType, int>(
871 const std::pair<unsigned int, unsigned int> &range)
mutable {
878 phi(matrix_free, range, dof_no, quad_no, first_selected_component);
880 internal::ComputeDiagonalHelper<dim,
888 for (
unsigned int cell = range.first; cell < range.second; ++cell)
892 for (
unsigned int i = 0; i < phi.dofs_per_cell; ++i)
894 helper.prepare_basis_vector(i);
899 helper.distribute_local_to_global(diagonal_global_components);
907 template <
typename CLASS,
913 typename VectorizedArrayType,
918 VectorType & diagonal_global,
924 VectorizedArrayType> &)
const,
925 const CLASS * owning_class,
926 const unsigned int dof_no,
927 const unsigned int quad_no,
928 const unsigned int first_selected_component)
939 [&](
auto &feeval) { (owning_class->*cell_operation)(feeval); },
942 first_selected_component);
951 template <
typename MatrixType,
953 typename std::enable_if<std::is_same<
954 typename std::remove_const<
typename std::remove_reference<
955 typename MatrixType::value_type>::type>::type,
956 typename std::remove_const<
typename std::remove_reference<
957 Number>::type>::type>::value>::type * =
nullptr>
959 create_new_affine_constraints_if_needed(
972 template <
typename MatrixType,
974 typename std::enable_if<!std::is_same<
975 typename std::remove_const<
typename std::remove_reference<
976 typename MatrixType::value_type>::type>::type,
977 typename std::remove_const<
typename std::remove_reference<
978 Number>::type>::type>::value>::type * =
nullptr>
980 create_new_affine_constraints_if_needed(
987 std::make_unique<AffineConstraints<typename MatrixType::value_type>>();
990 return *new_constraints;
999 typename VectorizedArrayType,
1000 typename MatrixType>
1005 MatrixType & matrix,
1011 VectorizedArrayType> &)> &local_vmult,
1012 const unsigned int dof_no,
1013 const unsigned int quad_no,
1014 const unsigned int first_selected_component)
1016 std::unique_ptr<AffineConstraints<typename MatrixType::value_type>>
1017 constraints_for_matrix;
1019 internal::create_new_affine_constraints_if_needed(matrix,
1021 constraints_for_matrix);
1023 matrix_free.template cell_loop<MatrixType, MatrixType>(
1024 [&](
const auto &,
auto &dst,
const auto &,
const auto range) {
1030 VectorizedArrayType>
1032 matrix_free, range, dof_no, quad_no, first_selected_component);
1034 unsigned int const dofs_per_cell = integrator.dofs_per_cell;
1036 std::vector<types::global_dof_index> dof_indices(dofs_per_cell);
1037 std::vector<types::global_dof_index> dof_indices_mf(dofs_per_cell);
1039 std::array<FullMatrix<typename MatrixType::value_type>,
1040 VectorizedArrayType::size()>
1043 std::fill_n(matrices.begin(),
1044 VectorizedArrayType::size(),
1048 const auto lexicographic_numbering =
1052 first_selected_component,
1053 integrator.get_active_fe_index(),
1054 integrator.get_active_quadrature_index())
1055 .lexicographic_numbering;
1057 for (
auto cell = range.first; cell < range.second; ++cell)
1059 integrator.reinit(cell);
1061 const unsigned int n_filled_lanes =
1064 for (
unsigned int v = 0; v < n_filled_lanes; ++v)
1067 for (
unsigned int j = 0; j < dofs_per_cell; ++j)
1069 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
1070 integrator.begin_dof_values()[i] =
1071 static_cast<Number
>(i == j);
1073 local_vmult(integrator);
1075 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
1076 for (
unsigned int v = 0; v < n_filled_lanes; ++v)
1077 matrices[v](i, j) = integrator.begin_dof_values()[i][v];
1080 for (
unsigned int v = 0; v < n_filled_lanes; ++v)
1086 cell_v->get_mg_dof_indices(dof_indices);
1088 cell_v->get_dof_indices(dof_indices);
1090 for (
unsigned int j = 0; j < dof_indices.size(); ++j)
1091 dof_indices_mf[j] = dof_indices[lexicographic_numbering[j]];
1105 template <
typename CLASS,
1111 typename VectorizedArrayType,
1112 typename MatrixType>
1117 MatrixType & matrix,
1123 VectorizedArrayType> &)
const,
1124 const CLASS * owning_class,
1125 const unsigned int dof_no,
1126 const unsigned int quad_no,
1127 const unsigned int first_selected_component)
1134 VectorizedArrayType,
1139 [&](
auto &feeval) { (owning_class->*cell_operation)(feeval); },
1142 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
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 AssertIndexRange(index, range)
static ::ExceptionBase & ExcMessage(std::string arg1)
@ matrix
Contents is actually a matrix.
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
std::uint8_t 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 > abs(const ::VectorizedArray< Number, width > &)
const ::Triangulation< dim, spacedim > & tria