39 template <
int dim,
typename Number,
bool is_face_>
52 std::function<std::vector<std::unique_ptr<FEEvalType>>(
53 const std::pair<unsigned int, unsigned int> &)>
55 std::function<void(std::vector<std::unique_ptr<FEEvalType>> &,
58 std::function<void(std::vector<std::unique_ptr<FEEvalType>> &)>
68 template <
int dim,
typename AdditionalData>
71 AdditionalData &additional_data);
91 typename VectorizedArrayType,
96 VectorType &diagonal_global,
102 VectorizedArrayType> &)>
104 const unsigned int dof_no = 0,
105 const unsigned int quad_no = 0,
106 const unsigned int first_selected_component = 0,
107 const unsigned int first_vector_component = 0);
118 typename QuadOperation>
123 const QuadOperation &quad_operation,
126 const unsigned int dof_no = 0,
127 const unsigned int quad_no = 0,
128 const unsigned int first_selected_component = 0,
129 const unsigned int first_vector_component = 0);
134 template <
typename CLASS,
140 typename VectorizedArrayType,
145 VectorType &diagonal_global,
151 VectorizedArrayType> &) const,
152 const CLASS *owning_class,
153 const unsigned
int dof_no = 0,
154 const unsigned
int quad_no = 0,
155 const unsigned
int first_selected_component = 0,
156 const unsigned
int first_vector_component = 0);
175 typename VectorizedArrayType,
180 VectorType &diagonal_global,
186 VectorizedArrayType> &)>
193 VectorizedArrayType> &,
199 VectorizedArrayType> &)>
206 VectorizedArrayType> &)>
208 const unsigned int dof_no = 0,
209 const unsigned int quad_no = 0,
210 const unsigned int first_selected_component = 0,
211 const unsigned int first_vector_component = 0);
218 template <
typename CLASS,
224 typename VectorizedArrayType,
229 VectorType &diagonal_global,
235 VectorizedArrayType> &) const,
241 VectorizedArrayType> &,
247 VectorizedArrayType> &)
254 VectorizedArrayType> &)
256 const CLASS *owning_class,
257 const unsigned
int dof_no = 0,
258 const unsigned
int quad_no = 0,
259 const unsigned
int first_selected_component = 0,
260 const unsigned
int first_vector_component = 0);
277 typename VectorizedArrayType,
289 VectorizedArrayType> &)>
291 const unsigned int dof_no = 0,
292 const unsigned int quad_no = 0,
293 const unsigned int first_selected_component = 0);
300 template <
typename CLASS,
306 typename VectorizedArrayType,
318 VectorizedArrayType> &) const,
319 const CLASS *owning_class,
320 const unsigned
int dof_no = 0,
321 const unsigned
int quad_no = 0,
322 const unsigned
int first_selected_component = 0);
333 typename VectorizedArrayType,
335 typename VectorType2>
345 VectorType &diagonal_global,
346 std::vector<VectorType2 *> &diagonal_global_components);
354 typename VectorizedArrayType,
385 typename VectorizedArrayType,
397 VectorizedArrayType> &)>
404 VectorizedArrayType> &,
410 VectorizedArrayType> &)>
417 VectorizedArrayType> &)>
419 const unsigned int dof_no = 0,
420 const unsigned int quad_no = 0,
421 const unsigned int first_selected_component = 0);
428 template <
typename CLASS,
434 typename VectorizedArrayType,
446 VectorizedArrayType> &) const,
452 VectorizedArrayType> &,
458 VectorizedArrayType> &)
465 VectorizedArrayType> &)
467 const CLASS *owning_class,
468 const unsigned
int dof_no = 0,
469 const unsigned
int quad_no = 0,
470 const unsigned
int first_selected_component = 0);
522 std::vector<unsigned int> valid_fe_indices;
524 const auto &fe_collection =
526 .get_fe_collection();
528 for (
unsigned int i = 0; i < fe_collection.size(); ++i)
529 if (fe_collection[i].n_dofs_per_cell() > 0)
530 valid_fe_indices.push_back(i);
544 template <
typename VectorTypeOut,
typename VectorTypeIn>
549 const VectorTypeIn &,
550 const std::pair<unsigned int, unsigned int> &)> &cell_operation,
552 const VectorTypeIn &src,
553 const bool zero_dst_vector =
false)
const
555 const auto ebd_cell_operation = [&](
const auto &
matrix_free,
559 const auto category =
matrix_free.get_cell_range_category(range);
567 matrix_free->template cell_loop<VectorTypeOut, VectorTypeIn>(
568 ebd_cell_operation, dst, src, zero_dst_vector);
578 template <
typename VectorTypeOut,
typename VectorTypeIn>
583 const VectorTypeIn &,
584 const std::pair<unsigned int, unsigned int> &)> &cell_operation,
588 const VectorTypeIn &,
589 const std::pair<unsigned int, unsigned int> &)> &face_operation,
593 const VectorTypeIn &,
594 const std::pair<unsigned int, unsigned int> &,
595 const bool)> &boundary_operation,
597 const VectorTypeIn &src,
598 const bool zero_dst_vector =
false)
const
600 const auto ebd_cell_operation = [&](
const auto &
matrix_free,
604 const auto category =
matrix_free.get_cell_range_category(range);
612 const auto ebd_internal_or_boundary_face_operation =
617 const auto category =
matrix_free.get_face_range_category(range);
619 const unsigned int type =
633 matrix_free->template loop<VectorTypeOut, VectorTypeIn>(
635 ebd_internal_or_boundary_face_operation,
636 ebd_internal_or_boundary_face_operation,
659 template <
int dim,
typename AdditionalData>
662 AdditionalData &additional_data)
665 const unsigned int level = additional_data.mg_level;
670 additional_data.cell_vectorization_category.assign(
673 additional_data.cell_vectorization_category.assign(tria.
n_active_cells(),
680 std::sort(bids.begin(), bids.end());
683 unsigned int n_bids = bids.size() + 1;
685 for (
unsigned int i = 0; i < GeometryInfo<dim>::faces_per_cell;
686 i++, offset = offset * n_bids)
690 const auto to_category = [&](
const auto &cell) {
691 unsigned int c_num = 0;
692 for (
unsigned int i = 0; i < GeometryInfo<dim>::faces_per_cell; ++i)
694 const auto face = cell->face(i);
695 if (face->at_boundary() && !cell->has_periodic_neighbor(i))
697 factors[i] * (1 + std::distance(bids.begin(),
698 std::find(bids.begin(),
700 face->boundary_id())));
709 if (cell->is_locally_owned())
711 .cell_vectorization_category[cell->active_cell_index()] =
719 if (cell->is_locally_owned_on_level())
720 additional_data.cell_vectorization_category[cell->index()] =
726 additional_data.hold_all_faces_to_owned_cells =
true;
727 additional_data.cell_vectorization_categories_strict =
true;
728 additional_data.mapping_update_flags_faces_by_cells =
729 additional_data.mapping_update_flags_inner_faces |
730 additional_data.mapping_update_flags_boundary_faces;
735 template <
typename Number>
738 std::vector<unsigned int> row_lid_to_gid;
739 std::vector<unsigned int> row;
740 std::vector<unsigned int> col;
741 std::vector<Number> val;
743 std::vector<unsigned int> inverse_lookup_rows;
744 std::vector<std::pair<unsigned int, unsigned int>> inverse_lookup_origins;
747 template <
int dim,
typename VectorizedArrayType,
bool is_face>
748 class ComputeDiagonalHelper
751 using FEEvaluationType =
754 using Number =
typename VectorizedArrayType::value_type;
755 static const unsigned int n_lanes = VectorizedArrayType::size();
757 ComputeDiagonalHelper()
759 , matrix_free(nullptr)
760 , dofs_per_component(0)
764 ComputeDiagonalHelper(
const ComputeDiagonalHelper &)
766 , matrix_free(nullptr)
767 , dofs_per_component(0)
773 FEEvaluationType &phi,
775 const unsigned int n_components)
779 if (dofs_per_component !=
780 phi.get_shape_info().dofs_per_component_on_cell)
782 locally_relevant_constraints_hn_map.
clear();
784 phi.get_shape_info().dofs_per_component_on_cell;
786 this->n_components = n_components;
787 this->dofs_per_cell = n_components * dofs_per_component;
789 this->matrix_free = &matrix_free;
793 reinit(
const unsigned int cell)
796 const auto &matrix_free = *this->matrix_free;
804 const unsigned int first_selected_component =
805 n_fe_components == 1 ? 0 : phi->get_first_selected_component();
807 this->n_lanes_filled =
816 std::vector<std::tuple<unsigned int, unsigned int, Number>>
817 locally_relevant_constraints, locally_relevant_constraints_tmp;
818 locally_relevant_constraints.reserve(dofs_per_cell);
819 std::vector<unsigned int> constraint_position;
820 std::vector<unsigned char> is_constrained_hn;
822 const std::array<unsigned int, n_lanes> &cells =
823 this->phi->get_cell_ids();
825 for (
unsigned int v = 0; v < n_lanes_filled; ++v)
830 const unsigned int *dof_indices;
831 unsigned int index_indicators, next_index_indicators;
833 const unsigned int start =
834 cells[v] * n_fe_components + first_selected_component;
836 dof_info.dof_indices.data() + dof_info.row_starts[start].first;
837 index_indicators = dof_info.row_starts[start].second;
838 next_index_indicators = dof_info.row_starts[start + 1].second;
842 locally_relevant_constraints.clear();
844 if (n_components == 1 || n_fe_components == 1)
846 unsigned int ind_local = 0;
847 for (; index_indicators != next_index_indicators;
848 ++index_indicators, ++ind_local)
850 const std::pair<unsigned short, unsigned short> indicator =
851 dof_info.constraint_indicator[index_indicators];
853 for (
unsigned int j = 0; j < indicator.first;
855 locally_relevant_constraints.emplace_back(ind_local,
859 dof_indices += indicator.first;
861 const Number *data_val =
863 const Number *end_pool =
866 for (; data_val != end_pool; ++data_val, ++dof_indices)
867 locally_relevant_constraints.emplace_back(ind_local,
874 for (; ind_local < dofs_per_component;
875 ++dof_indices, ++ind_local)
876 locally_relevant_constraints.emplace_back(ind_local,
887 for (
unsigned int comp = 0; comp < n_components; ++comp)
889 unsigned int ind_local = 0;
893 for (; index_indicators != next_index_indicators;
894 ++index_indicators, ++ind_local)
896 const std::pair<unsigned short, unsigned short>
898 dof_info.constraint_indicator[index_indicators];
901 for (
unsigned int j = 0; j < indicator.first;
903 locally_relevant_constraints.emplace_back(
904 comp * dofs_per_component + ind_local,
907 dof_indices += indicator.first;
909 const Number *data_val =
911 const Number *end_pool =
914 for (; data_val != end_pool; ++data_val, ++dof_indices)
915 locally_relevant_constraints.emplace_back(
916 comp * dofs_per_component + ind_local,
924 for (; ind_local < dofs_per_component;
925 ++dof_indices, ++ind_local)
926 locally_relevant_constraints.emplace_back(
927 comp * dofs_per_component + ind_local,
931 if (comp + 1 < n_components)
932 next_index_indicators =
933 dof_info.row_starts[start + comp + 2].second;
940 for (
unsigned int i = 1; i < locally_relevant_constraints.size();
942 Assert(std::get<0>(locally_relevant_constraints[i]) >=
943 std::get<0>(locally_relevant_constraints[i - 1]),
947 if (dof_info.hanging_node_constraint_masks.size() > 0 &&
948 dof_info.hanging_node_constraint_masks_comp.size() > 0 &&
949 dof_info.hanging_node_constraint_masks_comp
950 [phi->get_active_fe_index()][first_selected_component])
953 dof_info.hanging_node_constraint_masks[cells[v]];
956 if (mask != ::internal::MatrixFreeFunctions::
961 if (locally_relevant_constraints_hn_map.find(mask) ==
962 locally_relevant_constraints_hn_map.end())
963 fill_constraint_type_into_map(mask);
965 const auto &locally_relevant_constraints_hn =
966 locally_relevant_constraints_hn_map[
mask];
968 locally_relevant_constraints_tmp.clear();
969 if (locally_relevant_constraints_tmp.capacity() <
970 locally_relevant_constraints.size())
971 locally_relevant_constraints_tmp.reserve(
972 locally_relevant_constraints.size() +
973 locally_relevant_constraints_hn.size());
978 constraint_position.assign(dofs_per_cell,
980 for (
auto &a : locally_relevant_constraints)
981 if (constraint_position[
std::get<0>(a)] ==
983 constraint_position[
std::get<0>(a)] =
984 std::distance(locally_relevant_constraints.
data(),
986 is_constrained_hn.assign(dofs_per_cell,
false);
987 for (
auto &hn : locally_relevant_constraints_hn)
988 is_constrained_hn[
std::get<0>(hn)] = 1;
991 for (
const auto &a : locally_relevant_constraints)
992 if (is_constrained_hn[
std::get<0>(a)] == 0)
993 locally_relevant_constraints_tmp.push_back(a);
997 for (
const auto &hn : locally_relevant_constraints_hn)
998 if (constraint_position[
std::get<1>(hn)] !=
1002 locally_relevant_constraints.size());
1003 auto other = locally_relevant_constraints.begin() +
1004 constraint_position[std::get<1>(hn)];
1007 for (; other != locally_relevant_constraints.end() &&
1008 std::get<0>(*other) == std::get<1>(hn);
1010 locally_relevant_constraints_tmp.emplace_back(
1012 std::get<1>(*other),
1013 std::get<2>(hn) * std::get<2>(*other));
1016 std::swap(locally_relevant_constraints,
1017 locally_relevant_constraints_tmp);
1022 std::sort(locally_relevant_constraints.begin(),
1023 locally_relevant_constraints.end(),
1024 [](
const auto &a,
const auto &b) {
1025 if (std::get<1>(a) < std::get<1>(b))
1027 return (std::get<1>(a) == std::get<1>(b)) &&
1028 (std::get<0>(a) < std::get<0>(b));
1032 auto &c_pool = c_pools[v];
1034 c_pool.row_lid_to_gid.clear();
1036 c_pool.row.push_back(0);
1040 if (locally_relevant_constraints.size() > 0)
1041 c_pool.row_lid_to_gid.emplace_back(
1042 std::get<1>(locally_relevant_constraints.front()));
1043 for (
const auto &j : locally_relevant_constraints)
1045 if (c_pool.row_lid_to_gid.back() != std::get<1>(j))
1047 c_pool.row_lid_to_gid.push_back(std::get<1>(j));
1048 c_pool.row.push_back(c_pool.val.size());
1051 c_pool.col.emplace_back(std::get<0>(j));
1052 c_pool.val.emplace_back(std::get<2>(j));
1055 if (c_pool.val.size() > 0)
1056 c_pool.row.push_back(c_pool.val.size());
1058 c_pool.inverse_lookup_rows.clear();
1059 c_pool.inverse_lookup_rows.resize(1 + dofs_per_cell);
1060 for (
const unsigned int i : c_pool.col)
1061 c_pool.inverse_lookup_rows[1 + i]++;
1063 std::partial_sum(c_pool.inverse_lookup_rows.begin(),
1064 c_pool.inverse_lookup_rows.end(),
1065 c_pool.inverse_lookup_rows.begin());
1069 c_pool.inverse_lookup_origins.resize(c_pool.col.size());
1070 std::vector<unsigned int> inverse_lookup_count(dofs_per_cell);
1071 for (
unsigned int row = 0; row < c_pool.row.size() - 1; ++row)
1072 for (
unsigned int col = c_pool.row[row];
1073 col < c_pool.row[row + 1];
1076 const unsigned int index = c_pool.col[col];
1077 c_pool.inverse_lookup_origins
1078 [c_pool.inverse_lookup_rows[
index] +
1079 inverse_lookup_count[
index]] = std::make_pair(row, col);
1080 ++inverse_lookup_count[
index];
1109 for (
unsigned int v = 0; v < n_lanes_filled; ++v)
1110 diagonals_local_constrained[v].assign(
1111 c_pools[v].row_lid_to_gid.size() *
1112 (n_fe_components == 1 ? n_components : 1),
1116 bool use_fast_path =
true;
1118 for (
unsigned int v = 0; v < n_lanes_filled; ++v)
1120 auto &c_pool = c_pools[v];
1122 for (
unsigned int i = 0; i < c_pool.row.size() - 1; ++i)
1124 if ((c_pool.row[i + 1] - c_pool.row[i]) > 1)
1126 use_fast_path =
false;
1129 else if (((c_pool.row[i + 1] - c_pool.row[i]) == 1) &&
1130 (c_pool.val[c_pool.row[i]] != 1.0))
1132 use_fast_path =
false;
1137 if (use_fast_path ==
false)
1141 this->has_simple_constraints_ = use_fast_path;
1145 fill_constraint_type_into_map(
1146 const ::internal::MatrixFreeFunctions::compressed_constraint_kind
1149 auto &constraints_hn = locally_relevant_constraints_hn_map[
mask];
1154 const unsigned int degree =
1155 phi->get_shape_info().data.front().fe_degree;
1160 values_dofs.resize(dofs_per_component);
1163 VectorizedArrayType::size()>
1165 constraint_mask.fill(::internal::MatrixFreeFunctions::
1167 constraint_mask[0] =
mask;
1169 for (
unsigned int i = 0; i < dofs_per_component; ++i)
1171 for (
unsigned int j = 0; j < dofs_per_component; ++j)
1172 values_dofs[j] = VectorizedArrayType();
1173 values_dofs[i] = Number(1);
1178 VectorizedArrayType>::apply(1,
1180 phi->get_shape_info(),
1183 values_dofs.data());
1185 const Number tolerance =
1187 std::numeric_limits<Number>::epsilon() * 16);
1188 for (
unsigned int j = 0; j < dofs_per_component; ++j)
1189 if (
std::abs(values_dofs[j][0]) > tolerance &&
1191 std::abs(values_dofs[j][0] - Number(1)) > tolerance))
1192 constraints_hn.emplace_back(j, i, values_dofs[j][0]);
1196 const unsigned int n_hn_constraints = constraints_hn.size();
1197 constraints_hn.resize(n_hn_constraints * n_components);
1199 for (
unsigned int c = 1; c < n_components; ++c)
1200 for (
unsigned int i = 0; i < n_hn_constraints; ++i)
1201 constraints_hn[c * n_hn_constraints + i] = std::make_tuple(
1202 std::get<0>(constraints_hn[i]) + c * dofs_per_component,
1203 std::get<1>(constraints_hn[i]) + c * dofs_per_component,
1204 std::get<2>(constraints_hn[i]));
1208 prepare_basis_vector(
const unsigned int i)
1215 VectorizedArrayType *dof_values = phi->begin_dof_values();
1216 for (
unsigned int j = 0; j < dofs_per_cell; ++j)
1217 dof_values[j] = VectorizedArrayType();
1218 dof_values[i] = Number(1);
1224 VectorizedArrayType *dof_values = phi->begin_dof_values();
1225 for (
unsigned int j = 0; j < dofs_per_cell; ++j)
1226 dof_values[j] = VectorizedArrayType();
1235 const unsigned int n_fe_components =
1236 phi->get_dof_info().start_components.back();
1237 const unsigned int comp =
1238 n_fe_components == 1 ? i / dofs_per_component : 0;
1239 const unsigned int i_comp =
1240 n_fe_components == 1 ? (i % dofs_per_component) : i;
1244 for (
unsigned int v = 0; v < n_lanes_filled; ++v)
1246 const auto &c_pool = c_pools[v];
1248 for (
unsigned int jj = c_pool.inverse_lookup_rows[i_comp];
1249 jj < c_pool.inverse_lookup_rows[i_comp + 1];
1252 const unsigned int j = c_pool.inverse_lookup_origins[jj].first;
1255 for (
unsigned int k = c_pool.row[j]; k < c_pool.row[j + 1]; ++k)
1256 temp += c_pool.val[k] *
1257 phi->begin_dof_values()[comp * dofs_per_component +
1261 diagonals_local_constrained
1262 [v][j + comp * c_pools[v].row_lid_to_gid.size()] +=
1263 temp * c_pool.val[c_pool.inverse_lookup_origins[jj].second];
1268 template <
typename VectorType>
1270 distribute_local_to_global(std::vector<VectorType *> &diagonal_global)
1273 const unsigned int n_fe_components =
1274 phi->get_dof_info().start_components.back();
1276 if (n_fe_components == 1)
1279 for (
unsigned int v = 0; v < n_lanes_filled; ++v)
1283 for (
unsigned int j = 0; j < c_pools[v].row.size() - 1; ++j)
1284 for (
unsigned int comp = 0;
1285 comp < (n_fe_components == 1 ?
1286 static_cast<unsigned int>(n_components) :
1290 *diagonal_global[n_fe_components == 1 ? comp : 0],
1291 c_pools[v].row_lid_to_gid[j],
1292 diagonals_local_constrained
1293 [v][j + comp * c_pools[v].row_lid_to_gid.
size()]);
1297 has_simple_constraints()
const
1299 return has_simple_constraints_;
1303 FEEvaluationType *phi;
1306 unsigned int dofs_per_component;
1307 unsigned int dofs_per_cell;
1308 unsigned int n_components;
1312 unsigned int n_lanes_filled;
1314 std::array<internal::LocalCSR<Number>, n_lanes> c_pools;
1319 std::array<std::vector<Number>, n_lanes> diagonals_local_constrained;
1323 std::vector<std::tuple<unsigned int, unsigned int, Number>>>
1324 locally_relevant_constraints_hn_map;
1329 bool has_simple_constraints_;
1332 template <
bool is_face,
1335 typename VectorizedArrayType>
1339 const std::pair<unsigned int, unsigned int> &range,
1340 const unsigned int dof_no,
1341 const unsigned int quad_no,
1342 const unsigned int first_selected_component,
1343 const unsigned int fe_degree,
1344 const unsigned int n_q_points_1d,
1345 const bool is_interior_face =
true)
1347 const unsigned int static_n_q_points =
1351 unsigned int active_fe_index = 0;
1354 else if (is_interior_face)
1359 const auto init_data = ::internal::
1360 extract_initialization_data<is_face, dim, Number, VectorizedArrayType>(
1363 first_selected_component,
1371 return init_data.shape_info->dofs_per_component_on_cell == 0;
1381 typename QuadOperation>
1385 CellAction(
const QuadOperation &quad_operation,
1388 : m_quad_operation(quad_operation)
1389 , m_evaluation_flags(evaluation_flags)
1390 , m_integration_flags(integration_flags)
1393 KOKKOS_FUNCTION
void
1394 operator()(
const unsigned int cell,
1401 gpu_data, shared_data);
1402 m_quad_operation.set_matrix_free_data(*gpu_data);
1403 m_quad_operation.set_cell(cell);
1404 constexpr int dofs_per_cell =
decltype(fe_eval)::tensor_dofs_per_cell;
1405 Number
diagonal[dofs_per_cell] = {};
1406 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
1408 Kokkos::parallel_for(
1409 Kokkos::TeamThreadRange(shared_data->
team_member, dofs_per_cell),
1410 [&](
int j) { fe_eval.submit_dof_value(i == j ? 1 : 0, j); });
1412 Portable::internal::
1413 resolve_hanging_nodes<dim, fe_degree, false, Number>(
1417 Kokkos::subview(shared_data->
values, Kokkos::ALL, 0));
1419 fe_eval.evaluate(m_evaluation_flags);
1420 fe_eval.apply_for_each_quad_point(m_quad_operation);
1421 fe_eval.integrate(m_integration_flags);
1423 Portable::internal::
1424 resolve_hanging_nodes<dim, fe_degree, true, Number>(
1428 Kokkos::subview(shared_data->
values, Kokkos::ALL, 0));
1430 Kokkos::single(Kokkos::PerTeam(shared_data->
team_member),
1431 [&] { diagonal[i] = fe_eval.get_dof_value(i); });
1434 Kokkos::single(Kokkos::PerTeam(shared_data->
team_member), [&] {
1435 for (unsigned int i = 0; i < dofs_per_cell; ++i)
1436 fe_eval.submit_dof_value(diagonal[i], i);
1443 Kokkos::parallel_for(Kokkos::TeamThreadRange(shared_data->
team_member,
1446 dst[gpu_data->local_to_global(cell, i)] +=
1447 shared_data->values(i, 0);
1452 Kokkos::parallel_for(
1453 Kokkos::TeamThreadRange(shared_data->
team_member, dofs_per_cell),
1455 Kokkos::atomic_add(&dst[gpu_data->local_to_global(cell, i)],
1456 shared_data->values(i, 0));
1461 static constexpr unsigned int n_local_dofs = QuadOperation::n_local_dofs;
1464 mutable QuadOperation m_quad_operation;
1476 typename QuadOperation>
1481 const QuadOperation &quad_operation,
1484 const unsigned int dof_no,
1485 const unsigned int quad_no,
1486 const unsigned int first_selected_component,
1487 const unsigned int first_vector_component)
1497 CellAction<dim, fe_degree, n_q_points_1d, Number, QuadOperation>
1498 cell_action(quad_operation, evaluation_flags, integration_flags);
1500 matrix_free.
cell_loop(cell_action, dummy, diagonal_global);
1510 typename VectorizedArrayType,
1515 VectorType &diagonal_global,
1521 VectorizedArrayType> &)>
1523 const unsigned int dof_no,
1524 const unsigned int quad_no,
1525 const unsigned int first_selected_component,
1526 const unsigned int first_vector_component)
1533 VectorizedArrayType,
1541 first_selected_component,
1542 first_vector_component);
1545 template <
typename CLASS,
1551 typename VectorizedArrayType,
1556 VectorType &diagonal_global,
1562 VectorizedArrayType> &) const,
1563 const CLASS *owning_class,
1564 const unsigned
int dof_no,
1565 const unsigned
int quad_no,
1566 const unsigned
int first_selected_component,
1567 const unsigned
int first_vector_component)
1574 VectorizedArrayType,
1578 [&](
auto &phi) { (owning_class->*cell_operation)(phi); },
1581 first_selected_component,
1582 first_vector_component);
1590 typename VectorizedArrayType,
1595 VectorType &diagonal_global,
1601 VectorizedArrayType> &)>
1608 VectorizedArrayType> &,
1614 VectorizedArrayType> &)>
1621 VectorizedArrayType> &)>
1622 &boundary_operation,
1623 const unsigned int dof_no,
1624 const unsigned int quad_no,
1625 const unsigned int first_selected_component,
1626 const unsigned int first_vector_component)
1628 std::vector<typename ::internal::BlockVectorSelector<
1631 diagonal_global_components(n_components);
1633 for (
unsigned int d = 0;
d < n_components; ++
d)
1634 diagonal_global_components[d] = ::internal::
1636 get_vector_component(diagonal_global,
d + first_vector_component);
1638 const auto &dof_info = matrix_free.
get_dof_info(dof_no);
1640 if (dof_info.start_components.back() == 1)
1641 for (
unsigned int comp = 0; comp < n_components; ++comp)
1643 Assert(diagonal_global_components[comp] !=
nullptr,
1644 ExcMessage(
"The finite element underlying this FEEvaluation "
1645 "object is scalar, but you requested " +
1646 std::to_string(n_components) +
1647 " components via the template argument in "
1648 "FEEvaluation. In that case, you must pass an "
1649 "std::vector<VectorType> or a BlockVector to " +
1650 "read_dof_values and distribute_local_to_global."));
1652 *diagonal_global_components[comp], matrix_free, dof_info);
1657 *diagonal_global_components[0], matrix_free, dof_info);
1665 VectorizedArrayType>;
1672 VectorizedArrayType>;
1674 internal::ComputeMatrixScratchData<dim, VectorizedArrayType, false>
1677 data_cell.dof_numbers = {dof_no};
1678 data_cell.quad_numbers = {quad_no};
1679 data_cell.n_components = {n_components};
1680 data_cell.first_selected_components = {first_selected_component};
1681 data_cell.batch_type = {0};
1683 data_cell.op_create =
1684 [&](
const std::pair<unsigned int, unsigned int> &range) {
1686 std::unique_ptr<FEEvaluationData<dim, VectorizedArrayType, false>>>
1689 if (!internal::is_fe_nothing<false>(matrix_free,
1693 first_selected_component,
1696 phi.emplace_back(std::make_unique<FEEvalType>(
1697 matrix_free, range, dof_no, quad_no, first_selected_component));
1702 data_cell.op_reinit = [](
auto &phi,
const unsigned batch) {
1703 if (phi.size() == 1)
1704 static_cast<FEEvalType &
>(*phi[0]).reinit(batch);
1708 data_cell.op_compute = [&](
auto &phi) {
1709 cell_operation(
static_cast<FEEvalType &
>(*phi[0]));
1712 internal::ComputeMatrixScratchData<dim, VectorizedArrayType, true>
1715 data_face.dof_numbers = {dof_no, dof_no};
1716 data_face.quad_numbers = {quad_no, quad_no};
1717 data_face.n_components = {n_components, n_components};
1718 data_face.first_selected_components = {first_selected_component,
1719 first_selected_component};
1720 data_face.batch_type = {1, 2};
1722 data_face.op_create =
1723 [&](
const std::pair<unsigned int, unsigned int> &range) {
1725 std::unique_ptr<FEEvaluationData<dim, VectorizedArrayType, true>>>
1728 if (!internal::is_fe_nothing<true>(matrix_free,
1732 first_selected_component,
1736 !internal::is_fe_nothing<true>(matrix_free,
1740 first_selected_component,
1746 std::make_unique<FEFaceEvalType>(matrix_free,
1751 first_selected_component));
1753 std::make_unique<FEFaceEvalType>(matrix_free,
1758 first_selected_component));
1764 data_face.op_reinit = [](
auto &phi,
const unsigned batch) {
1765 if (phi.size() == 2)
1767 static_cast<FEFaceEvalType &
>(*phi[0]).
reinit(batch);
1768 static_cast<FEFaceEvalType &
>(*phi[1]).
reinit(batch);
1773 data_face.op_compute = [&](
auto &phi) {
1774 face_operation(
static_cast<FEFaceEvalType &
>(*phi[0]),
1775 static_cast<FEFaceEvalType &
>(*phi[1]));
1778 internal::ComputeMatrixScratchData<dim, VectorizedArrayType, true>
1781 data_boundary.dof_numbers = {dof_no};
1782 data_boundary.quad_numbers = {quad_no};
1783 data_boundary.n_components = {n_components};
1784 data_boundary.first_selected_components = {first_selected_component};
1785 data_boundary.batch_type = {1};
1788 .op_create = [&](
const std::pair<unsigned int, unsigned int> &range) {
1790 std::unique_ptr<FEEvaluationData<dim, VectorizedArrayType, true>>>
1793 if (!internal::is_fe_nothing<true>(matrix_free,
1797 first_selected_component,
1801 phi.emplace_back(std::make_unique<FEFaceEvalType>(
1802 matrix_free, range,
true, dof_no, quad_no, first_selected_component));
1807 data_boundary.op_reinit = [](
auto &phi,
const unsigned batch) {
1808 if (phi.size() == 1)
1809 static_cast<FEFaceEvalType &
>(*phi[0]).reinit(batch);
1812 if (boundary_operation)
1813 data_boundary.op_compute = [&](
auto &phi) {
1814 boundary_operation(
static_cast<FEFaceEvalType &
>(*phi[0]));
1822 diagonal_global_components);
1829 typename VectorizedArrayType,
1831 typename VectorType2>
1835 const internal::ComputeMatrixScratchData<dim, VectorizedArrayType, false>
1837 const internal::ComputeMatrixScratchData<dim, VectorizedArrayType, true>
1839 const internal::ComputeMatrixScratchData<dim, VectorizedArrayType, true>
1841 VectorType &diagonal_global,
1842 std::vector<VectorType2 *> &diagonal_global_components)
1849 internal::ComputeDiagonalHelper<dim, VectorizedArrayType, false>;
1852 internal::ComputeDiagonalHelper<dim, VectorizedArrayType, true>;
1856 scratch_data_internal;
1859 const auto batch_operation =
1862 const std::pair<unsigned int, unsigned int> &range) {
1863 if (!
data.op_compute)
1866 auto phi =
data.op_create(range);
1868 const unsigned int n_blocks = phi.size();
1870 auto &helpers = scratch_data.
get();
1871 helpers.resize(n_blocks);
1874 helpers[b].initialize(*phi[b], matrix_free,
data.n_components[b]);
1876 for (
unsigned int batch = range.first; batch < range.second; ++batch)
1878 data.op_reinit(phi, batch);
1881 helpers[b].
reinit(batch);
1885 Assert(std::all_of(helpers.begin(),
1887 [](
const auto &helper) {
1888 return helper.has_simple_constraints();
1895 for (
unsigned int i = 0;
1896 i < phi[
b]->get_shape_info().dofs_per_component_on_cell *
1897 data.n_components[b];
1900 for (
unsigned int bb = 0; bb <
n_blocks; ++bb)
1902 helpers[bb].prepare_basis_vector(i);
1904 helpers[bb].zero_basis_vector();
1906 data.op_compute(phi);
1907 helpers[
b].submit();
1910 helpers[
b].distribute_local_to_global(
1911 diagonal_global_components);
1916 const auto cell_operation_wrapped =
1917 [&](
const auto &,
auto &,
const auto &,
const auto range) {
1918 batch_operation(data_cell, scratch_data, range);
1921 const auto face_operation_wrapped =
1922 [&](
const auto &,
auto &,
const auto &,
const auto range) {
1923 batch_operation(data_face, scratch_data_internal, range);
1926 const auto boundary_operation_wrapped =
1927 [&](
const auto &,
auto &,
const auto &,
const auto range) {
1928 batch_operation(data_boundary, scratch_data_bc, range);
1931 if (data_face.op_compute || data_boundary.op_compute)
1932 matrix_free.template loop<VectorType, int>(cell_operation_wrapped,
1933 face_operation_wrapped,
1934 boundary_operation_wrapped,
1939 matrix_free.template cell_loop<VectorType, int>(cell_operation_wrapped,
1946 template <
typename CLASS,
1952 typename VectorizedArrayType,
1957 VectorType &diagonal_global,
1963 VectorizedArrayType> &) const,
1969 VectorizedArrayType> &,
1975 VectorizedArrayType> &)
1982 VectorizedArrayType> &)
1984 const CLASS *owning_class,
1985 const unsigned
int dof_no,
1986 const unsigned
int quad_no,
1987 const unsigned
int first_selected_component,
1988 const unsigned
int first_vector_component)
1995 VectorizedArrayType,
1999 [&](
auto &phi) { (owning_class->*cell_operation)(phi); },
2000 [&](
auto &phi_m,
auto &phi_p) {
2001 (owning_class->*face_operation)(phi_m, phi_p);
2003 [&](
auto &phi) { (owning_class->*boundary_operation)(phi); },
2006 first_selected_component,
2007 first_vector_component);
2019 std::enable_if_t<std::is_same_v<
2020 std::remove_const_t<
2021 std::remove_reference_t<typename MatrixType::value_type>>,
2022 std::remove_const_t<std::remove_reference_t<Number>>>> * =
nullptr>
2024 create_new_affine_constraints_if_needed(
2040 std::enable_if_t<!std::is_same_v<
2041 std::remove_const_t<
2042 std::remove_reference_t<typename MatrixType::value_type>>,
2043 std::remove_const_t<std::remove_reference_t<Number>>>> * =
nullptr>
2045 create_new_affine_constraints_if_needed(
2052 std::make_unique<AffineConstraints<typename MatrixType::value_type>>();
2053 new_constraints->
copy_from(constraints);
2055 return *new_constraints;
2064 typename VectorizedArrayType,
2076 VectorizedArrayType> &)>
2078 const unsigned int dof_no,
2079 const unsigned int quad_no,
2080 const unsigned int first_selected_component)
2087 VectorizedArrayType,
2096 first_selected_component);
2099 template <
typename CLASS,
2105 typename VectorizedArrayType,
2117 VectorizedArrayType> &) const,
2118 const CLASS *owning_class,
2119 const unsigned
int dof_no,
2120 const unsigned
int quad_no,
2121 const unsigned
int first_selected_component)
2128 VectorizedArrayType,
2133 [&](
auto &phi) { (owning_class->*cell_operation)(phi); },
2136 first_selected_component);
2144 typename VectorizedArrayType,
2156 VectorizedArrayType> &)>
2163 VectorizedArrayType> &,
2169 VectorizedArrayType> &)>
2176 VectorizedArrayType> &)>
2177 &boundary_operation,
2178 const unsigned int dof_no,
2179 const unsigned int quad_no,
2180 const unsigned int first_selected_component)
2187 VectorizedArrayType>;
2194 VectorizedArrayType>;
2196 internal::ComputeMatrixScratchData<dim, VectorizedArrayType, false>
2199 data_cell.dof_numbers = {dof_no};
2200 data_cell.quad_numbers = {quad_no};
2201 data_cell.n_components = {n_components};
2202 data_cell.first_selected_components = {first_selected_component};
2203 data_cell.batch_type = {0};
2205 data_cell.op_create =
2206 [&](
const std::pair<unsigned int, unsigned int> &range) {
2208 std::unique_ptr<FEEvaluationData<dim, VectorizedArrayType, false>>>
2211 if (!internal::is_fe_nothing<false>(matrix_free,
2215 first_selected_component,
2218 phi.emplace_back(std::make_unique<FEEvalType>(
2219 matrix_free, range, dof_no, quad_no, first_selected_component));
2224 data_cell.op_reinit = [](
auto &phi,
const unsigned batch) {
2225 if (phi.size() == 1)
2226 static_cast<FEEvalType &
>(*phi[0]).reinit(batch);
2230 data_cell.op_compute = [&](
auto &phi) {
2231 cell_operation(
static_cast<FEEvalType &
>(*phi[0]));
2234 internal::ComputeMatrixScratchData<dim, VectorizedArrayType, true>
2237 data_face.dof_numbers = {dof_no, dof_no};
2238 data_face.quad_numbers = {quad_no, quad_no};
2239 data_face.n_components = {n_components, n_components};
2240 data_face.first_selected_components = {first_selected_component,
2241 first_selected_component};
2242 data_face.batch_type = {1, 2};
2244 data_face.op_create =
2245 [&](
const std::pair<unsigned int, unsigned int> &range) {
2247 std::unique_ptr<FEEvaluationData<dim, VectorizedArrayType, true>>>
2250 if (!internal::is_fe_nothing<true>(matrix_free,
2254 first_selected_component,
2258 !internal::is_fe_nothing<true>(matrix_free,
2262 first_selected_component,
2268 std::make_unique<FEFaceEvalType>(matrix_free,
2273 first_selected_component));
2275 std::make_unique<FEFaceEvalType>(matrix_free,
2280 first_selected_component));
2286 data_face.op_reinit = [](
auto &phi,
const unsigned batch) {
2287 if (phi.size() == 2)
2289 static_cast<FEFaceEvalType &
>(*phi[0]).
reinit(batch);
2290 static_cast<FEFaceEvalType &
>(*phi[1]).
reinit(batch);
2295 data_face.op_compute = [&](
auto &phi) {
2296 face_operation(
static_cast<FEFaceEvalType &
>(*phi[0]),
2297 static_cast<FEFaceEvalType &
>(*phi[1]));
2300 internal::ComputeMatrixScratchData<dim, VectorizedArrayType, true>
2303 data_boundary.dof_numbers = {dof_no};
2304 data_boundary.quad_numbers = {quad_no};
2305 data_boundary.n_components = {n_components};
2306 data_boundary.first_selected_components = {first_selected_component};
2307 data_boundary.batch_type = {1};
2310 .op_create = [&](
const std::pair<unsigned int, unsigned int> &range) {
2312 std::unique_ptr<FEEvaluationData<dim, VectorizedArrayType, true>>>
2315 if (!internal::is_fe_nothing<true>(matrix_free,
2319 first_selected_component,
2323 phi.emplace_back(std::make_unique<FEFaceEvalType>(
2324 matrix_free, range,
true, dof_no, quad_no, first_selected_component));
2329 data_boundary.op_reinit = [](
auto &phi,
const unsigned batch) {
2330 if (phi.size() == 1)
2331 static_cast<FEFaceEvalType &
>(*phi[0]).reinit(batch);
2334 if (boundary_operation)
2335 data_boundary.op_compute = [&](
auto &phi) {
2336 boundary_operation(
static_cast<FEFaceEvalType &
>(*phi[0]));
2340 matrix_free, constraints_in, data_cell, data_face, data_boundary, matrix);
2347 typename VectorizedArrayType,
2353 const internal::ComputeMatrixScratchData<dim, VectorizedArrayType, false>
2355 const internal::ComputeMatrixScratchData<dim, VectorizedArrayType, true>
2357 const internal::ComputeMatrixScratchData<dim, VectorizedArrayType, true>
2361 std::unique_ptr<AffineConstraints<typename MatrixType::value_type>>
2362 constraints_for_matrix;
2364 internal::create_new_affine_constraints_if_needed(
2365 matrix, constraints_in, constraints_for_matrix);
2367 const auto batch_operation =
2368 [&matrix_free, &constraints, &
matrix](
2369 auto &
data,
const std::pair<unsigned int, unsigned int> &range) {
2370 if (!
data.op_compute)
2373 auto phi =
data.op_create(range);
2375 const unsigned int n_blocks = phi.size();
2384 n_blocks, VectorizedArrayType::size());
2387 std::array<FullMatrix<typename MatrixType::value_type>,
2388 VectorizedArrayType::size()>>
2389 matrices(n_blocks, n_blocks);
2394 .get_fe(phi[b]->get_active_fe_index());
2396 const auto component_base =
2399 const auto component_in_base =
2400 data.first_selected_components[
b] -
2405 data.dof_numbers[b],
2406 data.quad_numbers[b],
2408 phi[b]->get_active_fe_index(),
2409 phi[b]->get_active_quadrature_index());
2412 shape_info.dofs_per_component_on_cell *
data.n_components[
b];
2414 dof_indices[
b].resize(fe.n_dofs_per_cell());
2416 for (
unsigned int v = 0; v < VectorizedArrayType::size(); ++v)
2417 dof_indices_mf[b][v].resize(dofs_per_cell[b]);
2419 lexicographic_numbering[
b].insert(
2420 lexicographic_numbering[b].
begin(),
2421 shape_info.lexicographic_numbering.begin() +
2422 component_in_base * shape_info.dofs_per_component_on_cell,
2423 shape_info.lexicographic_numbering.begin() +
2424 (component_in_base +
data.n_components[
b]) *
2425 shape_info.dofs_per_component_on_cell);
2428 for (
unsigned int bj = 0; bj <
n_blocks; ++bj)
2429 for (
unsigned int bi = 0; bi <
n_blocks; ++bi)
2430 std::fill_n(matrices[bi][bj].
begin(),
2431 VectorizedArrayType::size(),
2433 dofs_per_cell[bi], dofs_per_cell[bj]));
2435 for (
auto batch = range.first; batch < range.second; ++batch)
2437 data.op_reinit(phi, batch);
2439 const unsigned int n_filled_lanes =
2444 for (
unsigned int v = 0; v < n_filled_lanes; ++v)
2448 (
data.batch_type[
b] == 0) ?
2449 (batch * VectorizedArrayType::size() + v) :
2450 ((
data.batch_type[
b] == 1) ?
2451 matrix_free.get_face_info(batch).cells_interior[v] :
2452 matrix_free.get_face_info(batch).cells_exterior[v]);
2457 data.dof_numbers[b]);
2461 cell_iterator->get_mg_dof_indices(dof_indices[b]);
2463 cell_iterator->get_dof_indices(dof_indices[b]);
2465 for (
unsigned int j = 0; j < dofs_per_cell[
b]; ++j)
2466 dof_indices_mf[b][v][j] =
2467 dof_indices[b][lexicographic_numbering[b][j]];
2470 for (
unsigned int bj = 0; bj <
n_blocks; ++bj)
2472 for (
unsigned int j = 0; j < dofs_per_cell[bj]; ++j)
2474 for (
unsigned int bi = 0; bi <
n_blocks; ++bi)
2475 for (
unsigned int i = 0; i < dofs_per_cell[bi]; ++i)
2476 phi[bi]->begin_dof_values()[i] =
2477 (bj == bi) ?
static_cast<Number
>(i == j) : 0.0;
2479 data.op_compute(phi);
2481 for (
unsigned int bi = 0; bi <
n_blocks; ++bi)
2482 for (
unsigned int i = 0; i < dofs_per_cell[bi]; ++i)
2483 for (
unsigned int v = 0; v < n_filled_lanes; ++v)
2484 matrices[bi][bj][v](i, j) =
2485 phi[bi]->begin_dof_values()[i][v];
2488 for (
unsigned int v = 0; v < n_filled_lanes; ++v)
2489 for (
unsigned int bi = 0; bi <
n_blocks; ++bi)
2497 matrices[bi][bi][v], dof_indices_mf[bi][v], matrix);
2500 matrices[bi][bj][v],
2501 dof_indices_mf[bi][v],
2502 dof_indices_mf[bj][v],
2508 const auto cell_operation_wrapped =
2509 [&](
const auto &,
auto &,
const auto &,
const auto range) {
2510 batch_operation(data_cell, range);
2513 const auto face_operation_wrapped =
2514 [&](
const auto &,
auto &,
const auto &,
const auto range) {
2515 batch_operation(data_face, range);
2518 const auto boundary_operation_wrapped =
2519 [&](
const auto &,
auto &,
const auto &,
const auto range) {
2520 batch_operation(data_boundary, range);
2523 if (data_face.op_compute || data_boundary.op_compute)
2525 matrix_free.template loop<MatrixType, MatrixType>(
2526 cell_operation_wrapped,
2527 face_operation_wrapped,
2528 boundary_operation_wrapped,
2533 matrix_free.template cell_loop<MatrixType, MatrixType>(
2534 cell_operation_wrapped, matrix, matrix);
2540 template <
typename CLASS,
2546 typename VectorizedArrayType,
2558 VectorizedArrayType> &) const,
2564 VectorizedArrayType> &,
2570 VectorizedArrayType> &)
2577 VectorizedArrayType> &)
2579 const CLASS *owning_class,
2580 const unsigned
int dof_no,
2581 const unsigned
int quad_no,
2582 const unsigned
int first_selected_component)
2589 VectorizedArrayType,
2594 [&](
auto &phi) { (owning_class->*cell_operation)(phi); },
2595 [&](
auto &phi_m,
auto &phi_p) {
2596 (owning_class->*face_operation)(phi_m, phi_p);
2598 [&](
auto &phi) { (owning_class->*boundary_operation)(phi); },
2601 first_selected_component);