15#ifndef dealii_tensor_product_matrix_h
16#define dealii_tensor_product_matrix_h
114template <
int dim,
typename Number,
int n_rows_1d = -1>
139 template <
typename T>
160 template <
typename T>
263 template <
typename Number>
272 std::pair<std::bitset<width>,
282 const auto &M_0 = left.second.first;
283 const auto &K_0 = left.second.second;
284 const auto &M_1 = right.second.first;
285 const auto &K_1 = right.second.second;
287 std::bitset<width>
mask;
289 for (
unsigned int v = 0; v <
width; ++v)
290 mask[v] = left.first[v] && right.first[v];
295 if (comparator(M_0, M_1))
297 else if (comparator(M_1, M_0))
299 else if (comparator(K_0, K_1))
336template <
int dim,
typename Number,
int n_rows_1d = -1>
342 std::bitset<::internal::VectorizedArrayTrait<Number>::width()>,
386 template <
typename T>
388 insert(
const unsigned int index,
const T &Ms,
const T &Ks);
516 template <
typename Number>
518 spectral_assembly(
const Number *mass_matrix,
519 const Number *derivative_matrix,
520 const unsigned int n_rows,
521 const unsigned int n_cols,
527 std::vector<bool> constrained_dofs(n_rows,
false);
529 for (
unsigned int i = 0; i < n_rows; ++i)
531 if (mass_matrix[i + i * n_rows] == 0.0)
533 Assert(derivative_matrix[i + i * n_rows] == 0.0,
536 for (
unsigned int j = 0; j < n_rows; ++j)
538 Assert(derivative_matrix[i + j * n_rows] == 0,
540 Assert(derivative_matrix[j + i * n_rows] == 0,
544 constrained_dofs[i] =
true;
548 const auto transpose_fill_nm = [&constrained_dofs](Number *out,
550 const unsigned int n,
551 const unsigned int m) {
552 for (
unsigned int mm = 0, c = 0; mm < m; ++mm)
553 for (
unsigned int nn = 0; nn < n; ++nn, ++c)
555 (mm == nn && constrained_dofs[mm]) ? Number(1.0) : in[c];
558 std::vector<::Vector<Number>> eigenvecs(n_rows);
562 transpose_fill_nm(&(mass_copy(0, 0)), mass_matrix, n_rows, n_cols);
563 transpose_fill_nm(&(deriv_copy(0, 0)), derivative_matrix, n_rows, n_cols);
568 for (
unsigned int i = 0, c = 0; i < n_rows; ++i)
569 for (
unsigned int j = 0; j < n_cols; ++j, ++c)
570 if (constrained_dofs[i] ==
false)
573 for (
unsigned int i = 0; i < n_rows; ++i, ++
eigenvalues)
579 template <std::
size_t dim,
typename Number>
586 const unsigned int n_rows_1d = mass_matrix[0].n_cols();
588 for (
unsigned int dir = 0; dir < dim; ++dir)
593 derivative_matrix[dir].n_rows());
595 derivative_matrix[dir].n_cols());
598 mass_matrix[dir].n_rows());
599 eigenvalues[dir].resize(mass_matrix[dir].n_cols());
600 internal::TensorProductMatrixSymmetricSum::spectral_assembly<Number>(
601 &(mass_matrix[dir](0, 0)),
602 &(derivative_matrix[dir](0, 0)),
603 mass_matrix[dir].n_rows(),
604 mass_matrix[dir].n_cols(),
612 template <std::
size_t dim,
typename Number, std::
size_t n_lanes>
623 const unsigned int n_rows_1d = mass_matrix[0].n_cols();
624 constexpr unsigned int macro_size =
626 const std::size_t nm_flat_size_max = n_rows_1d * n_rows_1d * macro_size;
627 const std::size_t n_flat_size_max = n_rows_1d * macro_size;
629 std::vector<Number> mass_matrix_flat;
630 std::vector<Number> deriv_matrix_flat;
631 std::vector<Number> eigenvalues_flat;
632 std::vector<Number> eigenvectors_flat;
633 mass_matrix_flat.resize(nm_flat_size_max);
634 deriv_matrix_flat.resize(nm_flat_size_max);
635 eigenvalues_flat.resize(n_flat_size_max);
636 eigenvectors_flat.resize(nm_flat_size_max);
637 std::array<unsigned int, macro_size> offsets_nm;
638 std::array<unsigned int, macro_size> offsets_n;
639 for (
unsigned int dir = 0; dir < dim; ++dir)
644 derivative_matrix[dir].n_rows());
646 derivative_matrix[dir].n_cols());
648 const unsigned int n_rows = mass_matrix[dir].n_rows();
649 const unsigned int n_cols = mass_matrix[dir].n_cols();
650 const unsigned int nm = n_rows * n_cols;
651 for (
unsigned int vv = 0; vv < macro_size; ++vv)
652 offsets_nm[vv] = nm * vv;
657 &(mass_matrix[dir](0, 0)),
659 mass_matrix_flat.data());
663 &(derivative_matrix[dir](0, 0)),
665 deriv_matrix_flat.data());
667 const Number *mass_cbegin = mass_matrix_flat.data();
668 const Number *deriv_cbegin = deriv_matrix_flat.data();
669 Number *eigenvec_begin = eigenvectors_flat.data();
670 Number *eigenval_begin = eigenvalues_flat.data();
671 for (
unsigned int lane = 0; lane < macro_size; ++lane)
672 internal::TensorProductMatrixSymmetricSum::spectral_assembly<
673 Number>(mass_cbegin + nm * lane,
674 deriv_cbegin + nm * lane,
677 eigenval_begin + n_rows * lane,
678 eigenvec_begin + nm * lane);
682 for (
unsigned int vv = 0; vv < macro_size; ++vv)
683 offsets_n[vv] = n_rows * vv;
686 eigenvalues_flat.data(),
691 eigenvectors_flat.data(),
699 template <std::
size_t dim,
typename Number>
700 inline std::array<Table<2, Number>, dim>
708 template <std::
size_t dim,
typename Number>
709 inline std::array<Table<2, Number>, dim>
712 std::array<Table<2, Number>, dim> mass_copy;
714 std::transform(mass_matrix.cbegin(),
726 template <std::
size_t dim,
typename Number>
727 inline std::array<Table<2, Number>, dim>
730 std::array<Table<2, Number>, dim> matrices;
732 std::fill(matrices.begin(), matrices.end(), matrix);
739 template <
int n_rows_1d_templated, std::
size_t dim,
typename Number>
744 const unsigned int n_rows_1d_non_templated,
745 const std::array<const Number *, dim> &mass_matrix,
746 const std::array<const Number *, dim> &derivative_matrix)
748 const unsigned int n_rows_1d = n_rows_1d_templated == 0 ?
749 n_rows_1d_non_templated :
754 Number *t = tmp.
begin();
761 eval({}, {}, {}, n_rows_1d, n_rows_1d);
765 const Number *A = derivative_matrix[0];
766 eval.template apply<0, false, false>(A, src, dst);
771 const Number *A0 = derivative_matrix[0];
772 const Number *M0 = mass_matrix[0];
773 const Number *A1 = derivative_matrix[1];
774 const Number *M1 = mass_matrix[1];
775 eval.template apply<0, false, false>(M0, src, t);
776 eval.template apply<1, false, false>(A1, t, dst);
777 eval.template apply<0, false, false>(A0, src, t);
778 eval.template apply<1, false, true>(M1, t, dst);
783 const Number *A0 = derivative_matrix[0];
784 const Number *M0 = mass_matrix[0];
785 const Number *A1 = derivative_matrix[1];
786 const Number *M1 = mass_matrix[1];
787 const Number *A2 = derivative_matrix[2];
788 const Number *M2 = mass_matrix[2];
789 eval.template apply<0, false, false>(M0, src, t + n);
790 eval.template apply<1, false, false>(M1, t + n, t);
791 eval.template apply<2, false, false>(A2, t, dst);
792 eval.template apply<1, false, false>(A1, t + n, t);
793 eval.template apply<0, false, false>(A0, src, t + n);
794 eval.template apply<1, false, true>(M1, t + n, t);
795 eval.template apply<2, false, true>(M2, t, dst);
804 template <
int n_rows_1d_templated, std::
size_t dim,
typename Number>
806 apply_inverse(Number *dst,
808 const unsigned int n_rows_1d_non_templated,
810 const std::array<const Number *, dim> &
eigenvalues,
811 const Number *inverted_eigenvalues =
nullptr)
813 const unsigned int n_rows_1d = n_rows_1d_templated == 0 ?
814 n_rows_1d_non_templated :
822 eval({}, {}, {}, n_rows_1d, n_rows_1d);
832 eval.template apply<0, true, false>(S, src, dst);
834 for (
unsigned int i = 0; i < n_rows_1d; ++i)
835 if (inverted_eigenvalues)
836 dst[i] *= inverted_eigenvalues[i];
840 eval.template apply<0, false, false>(S, dst, dst);
847 eval.template apply<0, true, false>(S0, src, dst);
848 eval.template apply<1, true, false>(S1, dst, dst);
850 for (
unsigned int i1 = 0, c = 0; i1 < n_rows_1d; ++i1)
851 for (
unsigned int i0 = 0; i0 < n_rows_1d; ++i0, ++c)
852 if (inverted_eigenvalues)
853 dst[c] *= inverted_eigenvalues[c];
857 eval.template apply<1, false, false>(S1, dst, dst);
858 eval.template apply<0, false, false>(S0, dst, dst);
866 eval.template apply<0, true, false>(S0, src, dst);
867 eval.template apply<1, true, false>(S1, dst, dst);
868 eval.template apply<2, true, false>(S2, dst, dst);
870 for (
unsigned int i2 = 0, c = 0; i2 < n_rows_1d; ++i2)
871 for (
unsigned int i1 = 0; i1 < n_rows_1d; ++i1)
872 for (
unsigned int i0 = 0; i0 < n_rows_1d; ++i0, ++c)
873 if (inverted_eigenvalues)
874 dst[c] *= inverted_eigenvalues[c];
879 eval.template apply<2, false, false>(S2, dst, dst);
880 eval.template apply<1, false, false>(S1, dst, dst);
881 eval.template apply<0, false, false>(S0, dst, dst);
890 template <
int n_rows_1d_templated, std::
size_t dim,
typename Number>
892 select_vmult(Number *dst,
895 const unsigned int n_rows_1d,
896 const std::array<const Number *, dim> &mass_matrix,
897 const std::array<const Number *, dim> &derivative_matrix);
901 template <
int n_rows_1d_templated, std::
size_t dim,
typename Number>
903 select_apply_inverse(Number *dst,
905 const unsigned int n_rows_1d,
907 const std::array<const Number *, dim> &
eigenvalues,
908 const Number *inverted_eigenvalues =
nullptr);
913template <
int dim,
typename Number,
int n_rows_1d>
918 for (
unsigned int d = 1;
d < dim; ++
d)
919 m *= mass_matrix[d].n_rows();
925template <
int dim,
typename Number,
int n_rows_1d>
930 for (
unsigned int d = 1;
d < dim; ++
d)
931 n *= mass_matrix[d].n_cols();
937template <
int dim,
typename Number,
int n_rows_1d>
943 std::lock_guard<std::mutex> lock(this->mutex);
944 this->vmult(dst_view, src_view, this->tmp_array);
949template <
int dim,
typename Number,
int n_rows_1d>
959 Number *dst = dst_view.
begin();
960 const Number *src = src_view.
begin();
962 std::array<const Number *, dim>
mass_matrix, derivative_matrix;
964 for (
unsigned int d = 0;
d < dim; ++
d)
967 derivative_matrix[
d] = &this->derivative_matrix[
d](0, 0);
970 const unsigned int n_rows_1d_non_templated = this->mass_matrix[0].n_rows();
973 internal::TensorProductMatrixSymmetricSum::vmult<
974 n_rows_1d == -1 ? 0 : n_rows_1d>(dst,
977 n_rows_1d_non_templated,
981 internal::TensorProductMatrixSymmetricSum::select_vmult<1>(
985 n_rows_1d_non_templated,
992template <
int dim,
typename Number,
int n_rows_1d>
1001 Number *dst = dst_view.
begin();
1002 const Number *src = src_view.
begin();
1006 for (
unsigned int d = 0;
d < dim; ++
d)
1012 const unsigned int n_rows_1d_non_templated = this->mass_matrix[0].n_rows();
1014 if (n_rows_1d != -1)
1015 internal::TensorProductMatrixSymmetricSum::apply_inverse<
1016 n_rows_1d == -1 ? 0 : n_rows_1d>(
1019 internal::TensorProductMatrixSymmetricSum::select_apply_inverse<1>(
1025template <
int dim,
typename Number,
int n_rows_1d>
1039template <
int dim,
typename Number,
int n_rows_1d>
1040template <
typename T>
1043 const T &derivative_matrix)
1045 reinit(mass_matrix, derivative_matrix);
1050template <
int dim,
typename Number,
int n_rows_1d>
1051template <
typename T>
1054 const T &mass_matrix,
1055 const T &derivative_matrix)
1058 internal::TensorProductMatrixSymmetricSum::convert<dim>(mass_matrix);
1059 this->derivative_matrix =
1060 internal::TensorProductMatrixSymmetricSum::convert<dim>(derivative_matrix);
1062 internal::TensorProductMatrixSymmetricSum::setup(this->mass_matrix,
1063 this->derivative_matrix,
1070template <
int dim,
typename Number,
int n_rows_1d>
1073 const bool precompute_inverse_diagonal)
1074 : compress_matrices(compress_matrices)
1075 , precompute_inverse_diagonal(precompute_inverse_diagonal)
1080template <
int dim,
typename Number,
int n_rows_1d>
1083 const AdditionalData &additional_data)
1084 : compress_matrices(additional_data.compress_matrices)
1085 , precompute_inverse_diagonal(additional_data.precompute_inverse_diagonal)
1090template <
int dim,
typename Number,
int n_rows_1d>
1093 const unsigned int size)
1095 if (compress_matrices ==
false)
1096 mass_and_derivative_matrices.resize(size * dim);
1103template <
int dim,
typename Number,
int n_rows_1d>
1104template <
typename T>
1107 const unsigned int index,
1112 internal::TensorProductMatrixSymmetricSum::convert<dim>(Ms_in);
1114 internal::TensorProductMatrixSymmetricSum::convert<dim>(Ks_in);
1116 for (
unsigned int d = 0;
d < dim; ++
d)
1118 if (compress_matrices ==
false)
1120 const MatrixPairType
matrix(Ms[d], Ks[d]);
1121 mass_and_derivative_matrices[
index * dim +
d] =
matrix;
1125 using VectorizedArrayTrait =
1128 std::bitset<VectorizedArrayTrait::width()>
mask;
1130 for (
unsigned int v = 0; v < VectorizedArrayTrait::width(); ++v)
1132 typename VectorizedArrayTrait::value_type a = 0.0;
1134 for (
unsigned int i = 0; i < Ms[
d].size(0); ++i)
1135 for (
unsigned int j = 0; j < Ms[
d].size(1); ++j)
1137 a +=
std::abs(VectorizedArrayTrait::get(Ms[d][i][j], v));
1138 a +=
std::abs(VectorizedArrayTrait::get(Ks[d][i][j], v));
1141 mask[v] = (a != 0.0);
1144 const MatrixPairTypeWithMask
matrix{
mask, {Ms[
d], Ks[
d]}};
1146 const auto ptr = cache.find(matrix);
1148 if (ptr != cache.end())
1150 const auto ptr_index = ptr->second;
1151 indices[
index * dim +
d] = ptr_index;
1154 for (
unsigned int v = 0; v < VectorizedArrayTrait::width();
1156 if ((mask[v] ==
true) && (ptr->first.first[v] ==
false))
1166 auto mask_new = ptr->first.first;
1167 auto Ms_new = ptr->first.second.first;
1168 auto Ks_new = ptr->first.second.second;
1170 for (
unsigned int v = 0; v < VectorizedArrayTrait::width();
1172 if (mask_new[v] ==
false && mask[v] ==
true)
1176 for (
unsigned int i = 0; i < Ms_new.size(0); ++i)
1177 for (
unsigned int j = 0; j < Ms_new.size(1); ++j)
1179 VectorizedArrayTrait::get(Ms_new[i][j], v) =
1180 VectorizedArrayTrait::get(Ms[d][i][j], v);
1181 VectorizedArrayTrait::get(Ks_new[i][j], v) =
1182 VectorizedArrayTrait::get(Ks[d][i][j], v);
1188 const MatrixPairTypeWithMask entry_new{mask_new,
1191 const auto ptr_ = cache.find(entry_new);
1194 cache[entry_new] = ptr_index;
1199 const auto size = cache.size();
1200 indices[
index * dim +
d] = size;
1209template <
int dim,
typename Number,
int n_rows_1d>
1213 const auto store = [&](
const unsigned int index,
1214 const MatrixPairType &M_and_K) {
1218 std::array<Table<2, Number>, 1> derivative_matrix;
1219 derivative_matrix[0] = M_and_K.second;
1224 internal::TensorProductMatrixSymmetricSum::setup(mass_matrix,
1229 for (
unsigned int i = 0, m = matrix_ptr[index], v = vector_ptr[index];
1233 for (
unsigned int j = 0; j <
mass_matrix[0].n_cols(); ++j, ++m)
1236 this->derivative_matrices[m] = derivative_matrix[0][i][j];
1244 if (compress_matrices ==
false)
1251 this->vector_ptr.resize(mass_and_derivative_matrices.size() + 1);
1252 this->matrix_ptr.resize(mass_and_derivative_matrices.size() + 1);
1254 for (
unsigned int i = 0; i < mass_and_derivative_matrices.size(); ++i)
1256 const auto &M = mass_and_derivative_matrices[i].first;
1258 this->vector_ptr[i + 1] = M.n_rows();
1259 this->matrix_ptr[i + 1] = M.n_rows() * M.n_cols();
1262 for (
unsigned int i = 0; i < mass_and_derivative_matrices.size(); ++i)
1264 this->vector_ptr[i + 1] += this->vector_ptr[i];
1265 this->matrix_ptr[i + 1] += this->matrix_ptr[i];
1268 this->mass_matrices.resize_fast(matrix_ptr.back());
1269 this->derivative_matrices.resize_fast(matrix_ptr.back());
1270 this->eigenvectors.resize_fast(matrix_ptr.back());
1271 this->eigenvalues.resize_fast(vector_ptr.back());
1273 for (
unsigned int i = 0; i < mass_and_derivative_matrices.size(); ++i)
1274 store(i, mass_and_derivative_matrices[i]);
1276 mass_and_derivative_matrices.clear();
1278 else if (cache.size() == indices.size())
1282 this->vector_ptr.resize(cache.size() + 1);
1283 this->matrix_ptr.resize(cache.size() + 1);
1285 std::map<unsigned int, MatrixPairType> inverted_cache;
1287 for (
const auto &i : cache)
1288 inverted_cache[i.second] = i.first.second;
1290 for (
unsigned int i = 0; i < indices.size(); ++i)
1292 const auto &M = inverted_cache[indices[i]].first;
1294 this->vector_ptr[i + 1] = M.n_rows();
1295 this->matrix_ptr[i + 1] = M.n_rows() * M.n_cols();
1298 for (
unsigned int i = 0; i < cache.size(); ++i)
1300 this->vector_ptr[i + 1] += this->vector_ptr[i];
1301 this->matrix_ptr[i + 1] += this->matrix_ptr[i];
1304 this->mass_matrices.resize_fast(matrix_ptr.back());
1305 this->derivative_matrices.resize_fast(matrix_ptr.back());
1306 this->eigenvectors.resize_fast(matrix_ptr.back());
1307 this->eigenvalues.resize_fast(vector_ptr.back());
1309 for (
unsigned int i = 0; i < indices.size(); ++i)
1310 store(i, inverted_cache[indices[i]]);
1319 this->vector_ptr.resize(cache.size() + 1);
1320 this->matrix_ptr.resize(cache.size() + 1);
1322 for (
const auto &i : cache)
1324 const auto &M = i.first.second.first;
1326 this->vector_ptr[i.second + 1] = M.n_rows();
1327 this->matrix_ptr[i.second + 1] = M.n_rows() * M.n_cols();
1330 for (
unsigned int i = 0; i < cache.size(); ++i)
1332 this->vector_ptr[i + 1] += this->vector_ptr[i];
1333 this->matrix_ptr[i + 1] += this->matrix_ptr[i];
1336 this->mass_matrices.resize_fast(matrix_ptr.back());
1337 this->derivative_matrices.resize_fast(matrix_ptr.back());
1338 this->eigenvectors.resize_fast(matrix_ptr.back());
1339 this->eigenvalues.resize_fast(vector_ptr.back());
1341 for (
const auto &i : cache)
1342 store(i.second, i.first.second);
1347 if (precompute_inverse_diagonal)
1352 for (
unsigned int i = 0; i < this->eigenvalues.size(); ++i)
1353 this->eigenvalues[i] = Number(1.0) / this->eigenvalues[i];
1354 std::swap(this->inverted_eigenvalues,
eigenvalues);
1364 std::vector<unsigned int> indices_ev;
1366 if (indices.size() > 0)
1369 const unsigned int n_cells = indices.size() / dim;
1370 std::map<std::array<unsigned int, dim>,
unsigned int> cache_ev;
1371 std::vector<unsigned int> cache_ev_idx(n_cells);
1373 for (
unsigned int i = 0, c = 0; i <
n_cells; ++i)
1375 std::array<unsigned int, dim> id;
1377 for (
unsigned int d = 0;
d < dim; ++
d, ++c)
1380 const auto id_ptr = cache_ev.find(
id);
1382 if (id_ptr == cache_ev.end())
1384 const auto size = cache_ev.size();
1385 cache_ev_idx[i] = size;
1386 cache_ev[id] = size;
1390 cache_ev_idx[i] = id_ptr->second;
1395 std::vector<unsigned int> new_indices;
1396 new_indices.reserve(indices.size() / dim * (dim + 1));
1398 for (
unsigned int i = 0, c = 0; i <
n_cells; ++i)
1400 for (
unsigned int d = 0;
d < dim; ++
d, ++c)
1401 new_indices.push_back(indices[c]);
1402 new_indices.push_back(cache_ev_idx[i]);
1406 indices_ev.resize(cache_ev.size() * dim);
1407 for (
const auto &entry : cache_ev)
1408 for (
unsigned int d = 0;
d < dim; ++
d)
1409 indices_ev[entry.second * dim + d] = entry.first[d];
1411 std::swap(this->indices, new_indices);
1415 const unsigned int n_diag =
1416 ((indices_ev.size() > 0) ? indices_ev.size() :
1417 (matrix_ptr.size() - 1)) /
1420 std::vector<unsigned int> new_vector_ptr(n_diag + 1, 0);
1421 std::vector<unsigned int> new_vector_n_rows_1d(n_diag, 0);
1423 for (
unsigned int i = 0; i < n_diag; ++i)
1425 const unsigned int c = (indices_ev.size() > 0) ?
1426 indices_ev[dim * i + 0] :
1429 const unsigned int n_rows = vector_ptr[c + 1] - vector_ptr[c];
1431 new_vector_n_rows_1d[i] = n_rows;
1435 for (
unsigned int i = 0; i < n_diag; ++i)
1436 new_vector_ptr[i + 1] += new_vector_ptr[i];
1438 this->inverted_eigenvalues.resize(new_vector_ptr.back());
1441 for (
unsigned int i = 0; i < n_diag; ++i)
1443 std::array<Number *, dim> evs;
1445 for (
unsigned int d = 0;
d < dim; ++
d)
1448 ->
eigenvalues[this->vector_ptr[(indices_ev.size() > 0) ?
1449 indices_ev[dim * i +
d] :
1452 const unsigned int mm = new_vector_n_rows_1d[i];
1455 for (
unsigned int i1 = 0, c = 0; i1 < mm; ++i1)
1456 for (
unsigned int i0 = 0; i0 < mm; ++i0, ++c)
1457 this->inverted_eigenvalues[new_vector_ptr[i] + c] =
1458 Number(1.0) / (evs[1][i1] + evs[0][i0]);
1462 for (
unsigned int i2 = 0, c = 0; i2 < mm; ++i2)
1463 for (
unsigned int i1 = 0; i1 < mm; ++i1)
1464 for (
unsigned int i0 = 0; i0 < mm; ++i0, ++c)
1465 this->inverted_eigenvalues[new_vector_ptr[i] + c] =
1466 Number(1.0) / (evs[2][i2] + evs[1][i1] + evs[0][i0]);
1471 std::swap(this->vector_ptr, new_vector_ptr);
1472 std::swap(this->vector_n_rows_1d, new_vector_n_rows_1d);
1475 this->eigenvalues.clear();
1481template <
int dim,
typename Number,
int n_rows_1d>
1488 Number *dst = dst_in.
begin();
1489 const Number *src = src_in.
begin();
1491 if (this->eigenvalues.empty() ==
false)
1495 unsigned int n_rows_1d_non_templated = 0;
1497 for (
unsigned int d = 0;
d < dim; ++
d)
1499 const unsigned int translated_index =
1500 (indices.size() > 0) ? indices[dim *
index +
d] : (dim *
index +
d);
1503 this->eigenvectors.data() + matrix_ptr[translated_index];
1505 this->eigenvalues.data() + vector_ptr[translated_index];
1506 n_rows_1d_non_templated =
1507 vector_ptr[translated_index + 1] - vector_ptr[translated_index];
1510 if (n_rows_1d != -1)
1511 internal::TensorProductMatrixSymmetricSum::apply_inverse<
1512 n_rows_1d == -1 ? 0 : n_rows_1d>(
1515 internal::TensorProductMatrixSymmetricSum::select_apply_inverse<1>(
1521 const Number *inverted_eigenvalues =
nullptr;
1522 unsigned int n_rows_1d_non_templated = 0;
1524 for (
unsigned int d = 0;
d < dim; ++
d)
1526 const unsigned int translated_index =
1527 (indices.size() > 0) ?
1528 indices[((dim == 1) ? 1 : (dim + 1)) * index + d] :
1532 this->eigenvectors.data() + matrix_ptr[translated_index];
1536 const unsigned int translated_index =
1537 ((indices.size() > 0) && (dim != 1)) ?
1538 indices[(dim + 1) * index + dim] :
1541 inverted_eigenvalues =
1542 this->inverted_eigenvalues.data() + vector_ptr[translated_index];
1543 n_rows_1d_non_templated =
1545 (vector_ptr[translated_index + 1] - vector_ptr[translated_index]) :
1546 vector_n_rows_1d[translated_index];
1549 if (n_rows_1d != -1)
1550 internal::TensorProductMatrixSymmetricSum::apply_inverse<
1551 n_rows_1d == -1 ? 0 : n_rows_1d>(dst,
1553 n_rows_1d_non_templated,
1556 inverted_eigenvalues);
1558 internal::TensorProductMatrixSymmetricSum::select_apply_inverse<1>(
1561 n_rows_1d_non_templated,
1564 inverted_eigenvalues);
1570template <
int dim,
typename Number,
int n_rows_1d>
1586template <
int dim,
typename Number,
int n_rows_1d>
1591 if (matrix_ptr.empty())
1594 return matrix_ptr.size() - 1;
void resize_fast(const size_type new_size)
std::complex< typename numbers::NumberTraits< number >::real_type > eigenvalue(const size_type i) const
void compute_generalized_eigenvalues_symmetric(LAPACKFullMatrix< number > &B, const number lower_bound, const number upper_bound, const number abs_accuracy, Vector< number > &eigenvalues, std::vector< Vector< number > > &eigenvectors, const types::blas_int itype=1)
std::size_t storage_size() const
AlignedVector< Number > mass_matrices
const bool precompute_inverse_diagonal
std::size_t memory_consumption() const
std::vector< unsigned int > matrix_ptr
void apply_inverse(const unsigned int index, const ArrayView< Number > &dst_in, const ArrayView< const Number > &src_in) const
AlignedVector< Number > eigenvectors
std::vector< unsigned int > vector_ptr
const bool compress_matrices
std::vector< MatrixPairType > mass_and_derivative_matrices
std::pair< Table< 2, Number >, Table< 2, Number > > MatrixPairType
void reserve(const unsigned int size)
void insert(const unsigned int index, const T &Ms, const T &Ks)
AlignedVector< Number > inverted_eigenvalues
std::pair< std::bitset<::internal::VectorizedArrayTrait< Number >::width()>, MatrixPairType > MatrixPairTypeWithMask
std::vector< unsigned int > indices
std::vector< unsigned int > vector_n_rows_1d
std::map< MatrixPairTypeWithMask, unsigned int, internal::TensorProductMatrixSymmetricSum::MatrixPairComparator< Number > > cache
AlignedVector< Number > eigenvalues
AlignedVector< Number > derivative_matrices
TensorProductMatrixSymmetricSumCollection(const AdditionalData &additional_data=AdditionalData())
void vmult(const ArrayView< Number > &dst, const ArrayView< const Number > &src, AlignedVector< Number > &tmp) const
std::array< Table< 2, Number >, dim > eigenvectors
std::array< Table< 2, Number >, dim > derivative_matrix
void reinit(const T &mass_matrix, const T &derivative_matrix)
void apply_inverse(const ArrayView< Number > &dst, const ArrayView< const Number > &src) const
void vmult(const ArrayView< Number > &dst, const ArrayView< const Number > &src) const
AlignedVector< Number > tmp_array
static constexpr int n_rows_1d_static
std::size_t memory_consumption() const
TensorProductMatrixSymmetricSum()=default
std::array< Table< 2, Number >, dim > mass_matrix
std::array< AlignedVector< Number >, dim > eigenvalues
TensorProductMatrixSymmetricSum(const T &mass_matrix, const T &derivative_matrix)
static constexpr std::size_t size()
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NOT_IMPLEMENTED()
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcInternalError()
#define AssertThrow(cond, exc)
@ matrix
Contents is actually a matrix.
void mass_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const double factor=1.)
std::enable_if_t< std::is_fundamental_v< T >, std::size_t > memory_consumption(const T &t)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
constexpr T fixed_power(const T t)
constexpr T pow(const T base, const int iexp)
unsigned int n_cells(const internal::TriangulationImplementation::NumberCache< 1 > &c)
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
constexpr unsigned int invalid_unsigned_int
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
bool precompute_inverse_diagonal
AdditionalData(const bool compress_matrices=true, const bool precompute_inverse_diagonal=true)
::internal::VectorizedArrayTrait< Number > VectorizedArrayTrait
std::pair< std::bitset< width >, std::pair< Table< 2, Number >, Table< 2, Number > > > MatrixPairType
typename VectorizedArrayTrait::value_type ScalarNumber
bool operator()(const MatrixPairType &left, const MatrixPairType &right) const
static constexpr std::size_t width
static constexpr std::size_t width()
std::array< Number, 1 > eigenvalues(const SymmetricTensor< 2, 1, Number > &T)
std::array< std::pair< Number, Tensor< 1, dim, Number > >, dim > eigenvectors(const SymmetricTensor< 2, dim, Number > &T, const SymmetricTensorEigenvectorMethod method=SymmetricTensorEigenvectorMethod::ql_implicit_shifts)
void vectorized_load_and_transpose(const unsigned int n_entries, const Number *in, const unsigned int *offsets, VectorizedArray< Number, width > *out)
void vectorized_transpose_and_store(const bool add_into, const unsigned int n_entries, const VectorizedArray< Number, width > *in, const unsigned int *offsets, Number *out)