12#ifndef ARCCORE_ALINA_COARSENING_H
13#define ARCCORE_ALINA_COARSENING_H
34#include "arccore/alina/BuiltinBackend.h"
35#include "arccore/alina/AlinaUtils.h"
36#include "arccore/alina/QRFactorizationImpl.h"
37#include "arccore/alina/Adapters.h"
38#include "arccore/alina/ValueTypeInterface.h"
45namespace Arcane::Alina::detail
54template <
class Matrix>
55std::shared_ptr<Matrix> galerkin(
const Matrix& A,
const Matrix& P,
const Matrix& R)
57 return product(R, *product(A, P));
65template <
class Matrix>
66std::shared_ptr<Matrix> scaled_galerkin(
const Matrix& A,
const Matrix& P,
const Matrix& R,
float s)
68 auto a = galerkin(A, P, R);
81namespace Arcane::Alina
87struct nullspace_params
97 std::vector<double>
B;
111 rows = p.get(
"rows", rows);
113 precondition(
cols > 0,
"Error in nullspace parameters: "
114 "B is set, but cols is not");
116 precondition(rows > 0,
"Error in nullspace parameters: "
117 "B is set, but rows is not");
119 B.assign(b, b + rows *
cols);
122 precondition(
cols == 0,
"Error in nullspace parameters: "
123 "cols > 0, but B is empty");
126 p.check_params( {
"cols",
"rows",
"B" });
129 void get(PropertyTree&,
const std::string&)
const {}
168 : ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p,
eps_strong)
170 p.check_params( {
"eps_strong",
"block_size" });
173 void get(PropertyTree& p,
const std::string& path)
const
175 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path,
eps_strong);
179 static const ptrdiff_t undefined = -1;
180 static const ptrdiff_t removed = -2;
199 std::vector<ptrdiff_t>
id;
207 template <
class Matrix>
211 ,
id(backend::nbRow(A))
213 typedef typename backend::value_type<Matrix>::type value_type;
214 typedef typename math::scalar_of<value_type>::type scalar_type;
218 const size_t n = backend::nbRow(A);
221 auto dia = diagonal(A);
223 for (ptrdiff_t i = begin; i < (begin + size); ++i) {
224 value_type eps_dia_i = eps_squared * (*dia)[i];
226 for (ptrdiff_t j = A.ptr[i], e = A.ptr[i + 1]; j < e; ++j) {
227 ptrdiff_t c = A.col[j];
228 value_type v = A.val[j];
239 for (
size_t i = 0; i < n; ++i) {
240 ptrdiff_t j = A.ptr[i], e = A.ptr[i + 1];
241 max_neib = std::max<size_t>(max_neib, e - j);
243 ptrdiff_t state = removed;
253 std::vector<ptrdiff_t> neib;
254 neib.reserve(max_neib);
257 for (
size_t i = 0; i < n; ++i) {
258 if (
id[i] != undefined)
263 ptrdiff_t cur_id =
static_cast<ptrdiff_t
>(
count++);
268 for (ptrdiff_t j = A.ptr[i], e = A.ptr[i + 1]; j < e; ++j) {
269 ptrdiff_t c = A.col[j];
279 for (ptrdiff_t c : neib) {
280 for (ptrdiff_t j = A.ptr[c], e = A.ptr[c + 1]; j < e; ++j) {
281 ptrdiff_t cc = A.col[j];
293 std::vector<ptrdiff_t> cnt(
count, 0);
294 for (ptrdiff_t i :
id)
297 std::partial_sum(cnt.begin(), cnt.end(), cnt.begin());
299 if (
static_cast<ptrdiff_t
>(
count) > cnt.back()) {
302 for (
size_t i = 0; i < n; ++i)
304 id[i] = cnt[
id[i]] - 1;
316 const std::vector<ptrdiff_t>& key;
319 skip_negative(
const std::vector<ptrdiff_t>& key,
int block_size)
321 , block_size(block_size)
324 bool operator()(ptrdiff_t i, ptrdiff_t j)
const
327 return static_cast<size_t>(key[i]) / block_size <
328 static_cast<size_t>(key[j]) / block_size;
343template <
class Matrix>
344std::shared_ptr<Matrix>
345tentative_prolongation(
size_t n,
347 const std::vector<ptrdiff_t> aggr,
351 typedef typename backend::value_type<Matrix>::type value_type;
352 typedef typename backend::col_type<Matrix>::type
col_type;
354 auto P = std::make_shared<Matrix>();
356 ARCCORE_ALINA_TIC(
"tentative");
357 if (nullspace.
cols > 0) {
358 ptrdiff_t nba = naggr / block_size;
362 std::vector<ptrdiff_t> order(n);
363 for (
size_t i = 0; i < n; ++i)
366 std::vector<ptrdiff_t> aggr_ptr(nba + 1, 0);
367 for (ptrdiff_t i = 0; i < static_cast<ptrdiff_t>(n); ++i) {
368 ptrdiff_t a = aggr[order[i]];
371 ++aggr_ptr[a / block_size + 1];
373 std::partial_sum(aggr_ptr.begin(), aggr_ptr.end(), aggr_ptr.begin());
378 P->set_size(n, nullspace.
cols * nba);
382 for (ptrdiff_t i = begin; i < (begin + size); ++i) {
383 P->ptr[i + 1] = aggr[i] < 0 ? 0 : nullspace.
cols;
392 std::vector<double> Bnew;
393 Bnew.resize(nba * nullspace.
cols * nullspace.
cols);
396 Alina::detail::QRFactorization<double> qr;
397 std::vector<double> Bpart;
399 for (ptrdiff_t i = begin; i < (begin + size); ++i) {
400 auto aggr_beg = aggr_ptr[i];
401 auto aggr_end = aggr_ptr[i + 1];
402 auto d = aggr_end - aggr_beg;
404 Bpart.resize(d * nullspace.
cols);
406 for (ptrdiff_t j = aggr_beg, jj = 0; j < aggr_end; ++j, ++jj) {
407 ptrdiff_t ib = nullspace.
cols * order[j];
408 for (
int k = 0; k < nullspace.
cols; ++k)
409 Bpart[jj + d * k] = nullspace.
B[ib + k];
412 qr.factorize(d, nullspace.
cols, &Bpart[0], Alina::detail::col_major);
414 for (
int ii = 0, kk = 0; ii < nullspace.
cols; ++ii)
415 for (
int jj = 0; jj < nullspace.
cols; ++jj, ++kk)
416 Bnew[i * nullspace.
cols * nullspace.
cols + kk] = qr.R(ii, jj);
418 for (ptrdiff_t j = aggr_beg, ii = 0; j < aggr_end; ++j, ++ii) {
419 col_type* c = &P->col[P->ptr[order[j]]];
420 value_type* v = &P->val[P->ptr[order[j]]];
422 for (
int jj = 0; jj < nullspace.
cols; ++jj) {
423 c[jj] = i * nullspace.
cols + jj;
426 v[jj] = qr.Q(ii, jj) * math::identity<value_type>();
432 std::swap(nullspace.
B, Bnew);
435 P->set_size(n, naggr);
438 for (ptrdiff_t i = begin; i < (begin + size); ++i) {
439 P->ptr[i + 1] = (aggr[i] >= 0);
443 P->set_nonzeros(P->scan_row_sizes());
446 for (ptrdiff_t i = begin; i < (begin + size); ++i) {
448 P->col[P->ptr[i]] = aggr[i];
449 P->val[P->ptr[i]] = math::identity<value_type>();
454 ARCCORE_ALINA_TOC(
"tentative");
488 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p,
block_size)
490 p.check_params( {
"eps_strong",
"block_size" });
495 plain_aggregates::params::get(p, path);
496 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, block_size);
500 static const ptrdiff_t undefined = -1;
501 static const ptrdiff_t removed = -2;
510 std::vector<ptrdiff_t>
id;
513 template <
class Matrix>
518 plain_aggregates aggr(A, prm);
520 remove_small_aggregates(A.nbRow(), 1, min_aggregate, aggr);
523 strong_connection.swap(aggr.strong_connection);
527 strong_connection.resize(backend::nonzeros(A));
528 id.resize(backend::nbRow(A));
530 auto ap = pointwise_matrix(A, prm.block_size);
533 plain_aggregates pw_aggr(Ap, prm);
535 remove_small_aggregates(Ap.nbRow(), prm.block_size, min_aggregate, pw_aggr);
537 count = pw_aggr.count * prm.block_size;
539 arccoreParallelFor(0, Ap.nbRow(), ForLoopRunInfo{}, [&](
Int32 begin,
Int32 size) {
540 std::vector<ptrdiff_t> j(prm.block_size);
541 std::vector<ptrdiff_t> e(prm.block_size);
543 for (ptrdiff_t ip = begin; ip < (begin + size); ++ip) {
544 ptrdiff_t ia = ip * prm.block_size;
546 for (unsigned k = 0; k < prm.block_size; ++k, ++ia) {
547 id[ia] = prm.block_size * pw_aggr.id[ip] + k;
550 e[k] = A.ptr[ia + 1];
553 for (ptrdiff_t jp = Ap.ptr[ip], ep = Ap.ptr[ip + 1]; jp < ep; ++jp) {
554 ptrdiff_t cp = Ap.col[jp];
555 bool sp = (cp == ip) || pw_aggr.strong_connection[jp];
557 ptrdiff_t col_end = (cp + 1) * prm.block_size;
559 for (unsigned k = 0; k < prm.block_size; ++k) {
560 ptrdiff_t beg = j[k];
561 ptrdiff_t end = e[k];
563 while (beg < end && A.col[beg] < col_end) {
564 strong_connection[beg] = sp && A.col[beg] != (ia + k);
576 static void remove_small_aggregates(
size_t n,
unsigned block_size,
unsigned min_aggregate,
577 plain_aggregates& aggr)
579 if (min_aggregate <= 1)
583 std::vector<ptrdiff_t> count(aggr.count, 0);
585 for (
size_t i = 0; i < n; ++i) {
586 ptrdiff_t
id = aggr.id[i];
594 for (
size_t i = 0; i < aggr.count; ++i) {
595 if (block_size * count[i] < min_aggregate) {
606 for (
size_t i = 0; i < n; ++i) {
607 ptrdiff_t
id = aggr.id[i];
609 aggr.id[i] = count[id];
638template <
class Backend>
639struct AggregationCoarsening
669 :
over_interp(math::static_rows<typename Backend::value_type>::value == 1 ? 1.5f : 2.0f)
673 : ARCCORE_ALINA_PARAMS_IMPORT_CHILD(p, aggr)
674 , ARCCORE_ALINA_PARAMS_IMPORT_CHILD(p, nullspace)
675 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, over_interp)
677 p.check_params( {
"aggr",
"nullspace",
"over_interp" });
680 void get(PropertyTree& p,
const std::string& path)
const
682 ARCCORE_ALINA_PARAMS_EXPORT_CHILD(p, path, aggr);
683 ARCCORE_ALINA_PARAMS_EXPORT_CHILD(p, path, nullspace);
684 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, over_interp);
688 explicit AggregationCoarsening(
const params& prm = params())
699 template <
class Matrix>
700 std::tuple<std::shared_ptr<Matrix>, std::shared_ptr<Matrix>>
703 const size_t n = backend::nbRow(A);
705 ARCCORE_ALINA_TIC(
"aggregates");
706 Aggregates aggr(A, prm.aggr, prm.nullspace.cols);
707 ARCCORE_ALINA_TOC(
"aggregates");
709 ARCCORE_ALINA_TIC(
"interpolation");
710 auto P = tentative_prolongation<Matrix>(
711 n, aggr.
count, aggr.
id, prm.nullspace, prm.aggr.block_size);
712 ARCCORE_ALINA_TOC(
"interpolation");
714 return std::make_tuple(P, transpose(*P));
725 template <
class Matrix>
726 std::shared_ptr<Matrix>
729 return detail::scaled_galerkin(A, P, R, 1 / prm.over_interp);
741template <
template <
class>
class Coarsening>
744 template <
class Backend>
747 typedef math::element_of<typename Backend::value_type>::type Scalar;
749 typedef Coarsening<BaseBackend> Base;
751 typedef typename Base::params params;
754 type(
const params& prm = params())
757 template <
class Matrix>
758 typename std::enable_if<backend::coarsening_is_supported<BaseBackend, Coarsening>::value &&
760 std::tuple<std::shared_ptr<Matrix>, std::shared_ptr<Matrix>>>::type
761 transfer_operators(
const Matrix& B)
763 typedef typename backend::value_type<Matrix>::type Block;
764 auto T = base.transfer_operators(*adapter::unblock_matrix(B));
766 auto& P = *std::get<0>(T);
767 auto& R = *std::get<1>(T);
772 return std::make_tuple(
773 std::make_shared<Matrix>(adapter::block_matrix<Block>(P)),
774 std::make_shared<Matrix>(adapter::block_matrix<Block>(R)));
777 template <
class Matrix>
778 typename std::enable_if<backend::coarsening_is_supported<BaseBackend, Coarsening>::value &&
780 std::tuple<std::shared_ptr<Matrix>, std::shared_ptr<Matrix>>>::type
781 transfer_operators(
const Matrix& A)
783 return base.transfer_operators(A);
786 template <
class Matrix>
787 typename std::enable_if<!backend::coarsening_is_supported<BaseBackend, Coarsening>::value,
788 std::tuple<std::shared_ptr<Matrix>, std::shared_ptr<Matrix>>>::type
789 transfer_operators(
const Matrix&)
791 throw std::logic_error(
"The coarsening is not supported by the backend");
794 template <
class Matrix>
795 std::shared_ptr<Matrix>
798 return base.coarse_operator(A, P, R);
811template <
class Vector>
812int rigid_body_modes(
int ndim,
const Vector& coo, std::vector<double>& B,
bool transpose =
false)
814 precondition(ndim == 2 || ndim == 3,
"Only 2D or 3D problems are supported");
815 precondition(coo.size() % ndim == 0,
"Coordinate vector size should be divisible by ndim");
817 size_t n = coo.size();
818 int nmodes = (ndim == 2 ? 3 : 6);
819 B.resize(n * nmodes, 0.0);
821 const int stride1 = transpose ? 1 : nmodes;
822 const int stride2 = transpose ? n : 1;
824 double sn = 1 / sqrt(n);
827 for (
size_t i = 0; i < n; ++i) {
828 size_t nod = i / ndim;
829 size_t dim = i % ndim;
831 double x = coo[nod * 2 + 0];
832 double y = coo[nod * 2 + 1];
835 B[i * stride1 + dim * stride2] = sn;
840 B[i * stride1 + 2 * stride2] = -y;
843 B[i * stride1 + 2 * stride2] = x;
848 else if (ndim == 3) {
849 for (
size_t i = 0; i < n; ++i) {
850 size_t nod = i / ndim;
851 size_t dim = i % ndim;
853 double x = coo[nod * 3 + 0];
854 double y = coo[nod * 3 + 1];
855 double z = coo[nod * 3 + 2];
858 B[i * stride1 + dim * stride2] = sn;
863 B[i * stride1 + 3 * stride2] = y;
864 B[i * stride1 + 5 * stride2] = z;
867 B[i * stride1 + 3 * stride2] = -x;
868 B[i * stride1 + 4 * stride2] = -z;
871 B[i * stride1 + 4 * stride2] = y;
872 B[i * stride1 + 5 * stride2] = -x;
879 std::array<double, 6>
dot;
880 for (
int i = ndim; i < nmodes; ++i) {
881 std::fill(
dot.begin(),
dot.end(), 0.0);
882 for (
size_t j = 0; j < n; ++j) {
883 for (
int k = 0; k < i; ++k)
884 dot[k] += B[j * stride1 + k * stride2] * B[j * stride1 + i * stride2];
887 for (
size_t j = 0; j < n; ++j) {
888 for (
int k = 0; k < i; ++k)
889 B[j * stride1 + i * stride2] -= dot[k] * B[j * stride1 + k * stride2];
890 s += B[j * stride1 + i * stride2] * B[j * stride1 + i * stride2];
893 for (
size_t j = 0; j < n; ++j)
894 B[j * stride1 + i * stride2] /= s;
908template <
class Backend>
909struct RugeStubenCoarsening
944 : ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p,
eps_strong)
945 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p,
do_trunc)
946 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p,
eps_trunc)
948 p.check_params( {
"eps_strong",
"do_trunc",
"eps_trunc" });
951 void get(
PropertyTree& p,
const std::string& path)
const
953 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, eps_strong);
954 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, do_trunc);
955 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, eps_trunc);
959 explicit RugeStubenCoarsening(
const params& prm = params())
963 template <
class Matrix>
964 std::tuple<std::shared_ptr<Matrix>, std::shared_ptr<Matrix>>
965 transfer_operators(
const Matrix& A)
const
967 typedef typename backend::value_type<Matrix>::type Val;
968 typedef typename backend::col_type<Matrix>::type Col;
969 typedef typename backend::ptr_type<Matrix>::type Ptr;
970 typedef typename math::scalar_of<Val>::type Scalar;
972 const size_t n = backend::nbRow(A);
974 static const Scalar eps = Alina::detail::eps<Scalar>(1);
976 static const Val zero = math::zero<Val>();
978 std::vector<char> cf(n,
'U');
979 CSRMatrix<char, Col, Ptr> S;
981 ARCCORE_ALINA_TIC(
"C/F split");
982 connect(A, prm.eps_strong, S, cf);
984 ARCCORE_ALINA_TOC(
"C/F split");
986 ARCCORE_ALINA_TIC(
"interpolation");
988 std::vector<ptrdiff_t> cidx(n);
989 for (
size_t i = 0; i < n; ++i)
991 cidx[i] =
static_cast<ptrdiff_t
>(nc++);
994 throw error::empty_level();
996 auto P = std::make_shared<Matrix>();
997 P->set_size(n, nc,
true);
999 std::vector<Val> Amin, Amax;
1006 arccoreParallelFor(0, n, ForLoopRunInfo{}, [&](Int32 begin, Int32 size) {
1007 for (ptrdiff_t i = begin; i < (begin + size); ++i) {
1014 Val amin = zero, amax = zero;
1016 for (ptrdiff_t j = A.ptr[i], e = A.ptr[i + 1]; j < e; ++j) {
1017 if (!S.val[j] || cf[A.col[j]] !=
'C')
1020 amin = std::min(amin, A.val[j]);
1021 amax = std::max(amax, A.val[j]);
1024 Amin[i] = (amin *= prm.eps_trunc);
1025 Amax[i] = (amax *= prm.eps_trunc);
1027 for (ptrdiff_t j = A.ptr[i], e = A.ptr[i + 1]; j < e; ++j) {
1028 if (!S.val[j] || cf[A.col[j]] !=
'C')
1031 if (A.val[j] < amin || amax < A.val[j])
1036 for (ptrdiff_t j = A.ptr[i], e = A.ptr[i + 1]; j < e; ++j)
1037 if (S.val[j] && cf[A.col[j]] ==
'C')
1043 P->set_nonzeros(P->scan_row_sizes());
1045 arccoreParallelFor(0, n, ForLoopRunInfo{}, [&](Int32 begin, Int32 size) {
1046 for (ptrdiff_t i = begin; i < (begin + size); ++i) {
1047 ptrdiff_t row_head = P->ptr[i];
1050 P->col[row_head] = cidx[i];
1051 P->val[row_head] = math::identity<Val>();
1056 Val a_num = zero, a_den = zero;
1057 Val b_num = zero, b_den = zero;
1058 Val d_neg = zero, d_pos = zero;
1060 for (ptrdiff_t j = A.ptr[i], e = A.ptr[i + 1]; j < e; ++j) {
1061 ptrdiff_t c = A.col[j];
1071 if (S.val[j] && cf[c] ==
'C') {
1073 if (prm.do_trunc && Amin[i] < v)
1079 if (S.val[j] && cf[c] ==
'C') {
1081 if (prm.do_trunc && v < Amax[i])
1091 if (math::norm(
static_cast<Val
>(a_den - d_neg)) > eps)
1092 cf_neg = math::norm(a_den) / math::norm(
static_cast<Val
>(a_den - d_neg));
1094 if (math::norm(
static_cast<Val
>(b_den - d_pos)) > eps)
1095 cf_pos = math::norm(b_den) / math::norm(
static_cast<Val
>(b_den - d_pos));
1098 if (zero < b_num && math::norm(b_den) < eps)
1101 Scalar alpha = math::norm(a_den) > eps ? -cf_neg * math::norm(a_num) / (math::norm(dia) * math::norm(a_den)) : 0;
1102 Scalar beta = math::norm(b_den) > eps ? -cf_pos * math::norm(b_num) / (math::norm(dia) * math::norm(b_den)) : 0;
1104 for (ptrdiff_t j = A.ptr[i], e = A.ptr[i + 1]; j < e; ++j) {
1105 ptrdiff_t c = A.col[j];
1108 if (!S.val[j] || cf[c] !=
'C')
1110 if (prm.do_trunc && Amin[i] <= v && v <= Amax[i])
1113 P->col[row_head] = cidx[c];
1114 P->val[row_head] = (v < zero ? alpha : beta) * v;
1119 ARCCORE_ALINA_TOC(
"interpolation");
1121 return std::make_tuple(P,
transpose(*P));
1124 template <
class Matrix>
1125 std::shared_ptr<Matrix>
1126 coarse_operator(
const Matrix& A,
const Matrix& P,
const Matrix& R)
const
1128 return detail::galerkin(A, P, R);
1140 template <
typename Val,
typename Col,
typename Ptr>
1141 static void connect(CSRMatrix<Val, Col, Ptr>
const& A,
float eps_strong,
1142 CSRMatrix<char, Col, Ptr>& S,
1143 std::vector<char>& cf)
1145 typedef typename math::scalar_of<Val>::type Scalar;
1147 const size_t n = backend::nbRow(A);
1148 const size_t nnz = backend::nonzeros(A);
1149 const Scalar eps = Alina::detail::eps<Scalar>(1);
1153 S.ptr.resize(n + 1);
1157 arccoreParallelFor(0, n, ForLoopRunInfo{}, [&](Int32 begin, Int32 size) {
1158 for (ptrdiff_t i = begin; i < (begin + size); ++i) {
1162 Val a_min = math::zero<Val>();
1164 for (
auto a = backend::row_begin(A, i); a; ++a)
1166 a_min = std::min(a_min, a.value());
1168 if (math::norm(a_min) < eps) {
1173 a_min *= eps_strong;
1175 for (Ptr j = A.ptr[i], e = A.ptr[i + 1]; j < e; ++j)
1176 S.val[j] = (A.col[j] != i && A.val[j] < a_min);
1181 for (
size_t i = 0; i < nnz; ++i)
1183 ++(S.ptr[A.col[i] + 1]);
1186 S.col.resize(S.ptr[n]);
1188 for (
size_t i = 0; i < n; ++i)
1189 for (Ptr j = A.ptr[i], e = A.ptr[i + 1]; j < e; ++j)
1191 S.col[S.ptr[A.col[j]]++] = i;
1193 std::rotate(S.ptr.data(), S.ptr.data() + n, S.ptr.data() + n + 1);
1198 template <
typename Val,
typename Col,
typename Ptr>
1199 static void cfsplit(CSRMatrix<Val, Col, Ptr>
const& A,
1200 CSRMatrix<char, Col, Ptr>
const& S,
1201 std::vector<char>& cf)
1203 const size_t n = A.nbRow();
1205 std::vector<Col> lambda(n);
1208 for (
size_t i = 0; i < n; ++i) {
1210 for (Ptr j = S.ptr[i], e = S.ptr[i + 1]; j < e; ++j)
1211 temp += (cf[S.col[j]] ==
'U' ? 1 : 2);
1220 std::vector<Ptr> ptr(n + 1, 0);
1221 std::vector<Ptr> cnt(n, 0);
1222 std::vector<Ptr> i2n(n);
1223 std::vector<Ptr> n2i(n);
1225 for (
size_t i = 0; i < n; ++i)
1226 ++ptr[lambda[i] + 1];
1228 std::partial_sum(ptr.begin(), ptr.end(), ptr.begin());
1230 for (
size_t i = 0; i < n; ++i) {
1231 Col lam = lambda[i];
1232 Ptr idx = ptr[lam] + cnt[lam]++;
1241 for (
size_t top = n; top-- > 0;) {
1243 Col lam = lambda[i];
1246 std::replace(cf.begin(), cf.end(),
'U',
'C');
1260 for (Ptr j = S.ptr[i], e = S.ptr[i + 1]; j < e; ++j) {
1269 for (Ptr aj = A.ptr[c], ae = A.ptr[c + 1]; aj < ae; ++aj) {
1274 Col lam_a = lambda[ac];
1276 if (cf[ac] !=
'U' ||
static_cast<size_t>(lam_a) + 1 >= n)
1279 Ptr old_pos = n2i[ac];
1280 Ptr new_pos = ptr[lam_a] + cnt[lam_a] - 1;
1282 n2i[i2n[old_pos]] = new_pos;
1283 n2i[i2n[new_pos]] = old_pos;
1285 std::swap(i2n[old_pos], i2n[new_pos]);
1289 ptr[lam_a + 1] = ptr[lam_a] + cnt[lam_a];
1291 lambda[ac] = lam_a + 1;
1296 for (Ptr j = A.ptr[i], e = A.ptr[i + 1]; j < e; j++) {
1301 Col lam = lambda[c];
1303 if (cf[c] !=
'U' || lam == 0)
1306 Ptr old_pos = n2i[c];
1307 Ptr new_pos = ptr[lam];
1309 n2i[i2n[old_pos]] = new_pos;
1310 n2i[i2n[new_pos]] = old_pos;
1312 std::swap(i2n[old_pos], i2n[new_pos]);
1317 lambda[c] = lam - 1;
1331template <
class Backend>
1332struct SmoothedAggregationCoarserning
1377 bool estimate_spectral_radius =
false;
1381 int power_iters = 0;
1386 : ARCCORE_ALINA_PARAMS_IMPORT_CHILD(p,
aggr)
1387 , ARCCORE_ALINA_PARAMS_IMPORT_CHILD(p,
nullspace)
1388 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p,
relax)
1389 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, estimate_spectral_radius)
1390 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, power_iters)
1392 p.check_params({
"aggr",
"nullspace",
"relax",
"estimate_spectral_radius",
"power_iters" });
1395 void get(
PropertyTree& p,
const std::string& path)
const
1397 ARCCORE_ALINA_PARAMS_EXPORT_CHILD(p, path, aggr);
1398 ARCCORE_ALINA_PARAMS_EXPORT_CHILD(p, path, nullspace);
1399 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, relax);
1400 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, estimate_spectral_radius);
1401 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, power_iters);
1405 SmoothedAggregationCoarserning(
const params& prm = params())
1409 template <
class Matrix>
1410 std::tuple<std::shared_ptr<Matrix>, std::shared_ptr<Matrix>>
1411 transfer_operators(
const Matrix& A)
1413 typedef typename backend::value_type<Matrix>::type value_type;
1414 typedef typename math::scalar_of<value_type>::type scalar_type;
1416 const size_t n = backend::nbRow(A);
1418 ARCCORE_ALINA_TIC(
"aggregates");
1419 Aggregates aggr(A, prm.aggr, prm.nullspace.cols);
1420 prm.aggr.eps_strong *= 0.5;
1421 ARCCORE_ALINA_TOC(
"aggregates");
1423 auto P_tent = tentative_prolongation<Matrix>(n, aggr.count, aggr.id, prm.nullspace, prm.aggr.
block_size);
1425 auto P = std::make_shared<Matrix>();
1426 P->set_size(backend::nbRow(*P_tent), backend::nbColumn(*P_tent),
true);
1428 scalar_type omega = prm.relax;
1429 if (prm.estimate_spectral_radius) {
1430 omega *=
static_cast<scalar_type
>(4.0 / 3) / spectral_radius<true>(A, prm.power_iters);
1433 omega *=
static_cast<scalar_type
>(2.0 / 3);
1436 ARCCORE_ALINA_TIC(
"smoothing");
1437 arccoreParallelFor(0, n, ForLoopRunInfo{}, [&](Int32 begin, Int32 size) {
1438 std::vector<ptrdiff_t> marker(P->ncols, -1);
1441 for (ptrdiff_t i = begin; i < (begin + size); ++i) {
1442 for (ptrdiff_t ja = A.ptr[i], ea = A.ptr[i + 1]; ja < ea; ++ja) {
1443 ptrdiff_t ca = A.col[ja];
1446 if (ca != i && !aggr.strong_connection[ja])
1449 for (ptrdiff_t jp = P_tent->ptr[ca], ep = P_tent->ptr[ca + 1]; jp < ep; ++jp) {
1450 ptrdiff_t cp = P_tent->col[jp];
1452 if (marker[cp] != i) {
1461 P->scan_row_sizes();
1464 arccoreParallelFor(0, n, ForLoopRunInfo{}, [&](Int32 begin, Int32 size) {
1465 std::vector<ptrdiff_t> marker(P->ncols, -1);
1468 for (ptrdiff_t i = begin; i < (begin + size); ++i) {
1472 value_type dia = math::zero<value_type>();
1473 for (ptrdiff_t j = A.ptr[i], e = A.ptr[i + 1]; j < e; ++j) {
1474 if (A.col[j] == i || !aggr.strong_connection[j])
1477 if (!math::is_zero(dia))
1478 dia = -omega * math::inverse(dia);
1480 ptrdiff_t row_beg = P->ptr[i];
1481 ptrdiff_t row_end = row_beg;
1482 for (ptrdiff_t ja = A.ptr[i], ea = A.ptr[i + 1]; ja < ea; ++ja) {
1483 ptrdiff_t ca = A.col[ja];
1486 if (ca != i && !aggr.strong_connection[ja])
1489 value_type va = (ca == i)
1490 ?
static_cast<value_type
>(
static_cast<scalar_type
>(1 - omega) * math::identity<value_type>())
1493 for (ptrdiff_t jp = P_tent->ptr[ca], ep = P_tent->ptr[ca + 1]; jp < ep; ++jp) {
1494 ptrdiff_t cp = P_tent->col[jp];
1495 value_type vp = P_tent->val[jp];
1497 if (marker[cp] < row_beg) {
1498 marker[cp] = row_end;
1499 P->col[row_end] = cp;
1500 P->val[row_end] = va * vp;
1504 P->val[marker[cp]] += va * vp;
1510 ARCCORE_ALINA_TOC(
"smoothing");
1512 return std::make_tuple(P,
transpose(*P));
1515 template <
class Matrix>
1516 std::shared_ptr<Matrix>
1517 coarse_operator(
const Matrix& A,
const Matrix& P,
const Matrix& R)
const
1519 return detail::galerkin(A, P, R);
1533template <
class Backend>
1534struct SmoothedAggregationEnergyMinCoarsening
1550 : ARCCORE_ALINA_PARAMS_IMPORT_CHILD(p, aggr)
1551 , ARCCORE_ALINA_PARAMS_IMPORT_CHILD(p, nullspace)
1553 p.check_params( {
"aggr",
"nullspace" });
1556 void get(PropertyTree& p,
const std::string& path)
const
1558 ARCCORE_ALINA_PARAMS_EXPORT_CHILD(p, path, aggr);
1559 ARCCORE_ALINA_PARAMS_EXPORT_CHILD(p, path, nullspace);
1563 SmoothedAggregationEnergyMinCoarsening(
const params& prm = params())
1567 template <
class Matrix>
1568 std::tuple<std::shared_ptr<Matrix>, std::shared_ptr<Matrix>>
1569 transfer_operators(
const Matrix& A)
1571 typedef typename backend::value_type<Matrix>::type Val;
1572 typedef typename backend::col_type<Matrix>::type Col;
1573 typedef typename backend::ptr_type<Matrix>::type Ptr;
1574 typedef ptrdiff_t Idx;
1576 ARCCORE_ALINA_TIC(
"aggregates");
1577 Aggregates aggr(A, prm.aggr, prm.nullspace.cols);
1578 prm.aggr.eps_strong *= 0.5;
1579 ARCCORE_ALINA_TOC(
"aggregates");
1581 ARCCORE_ALINA_TIC(
"interpolation");
1582 auto P_tent = tentative_prolongation<Matrix>(backend::nbRow(A), aggr.count, aggr.id, prm.nullspace, prm.aggr.
block_size);
1585 CSRMatrix<Val, Col, Ptr> Af;
1586 Af.set_size(backend::nbRow(A), backend::nbColumn(A));
1589 std::vector<Val> dia(Af.nbRow());
1591 arccoreParallelFor(0, Af.nbRow(), ForLoopRunInfo{}, [&](Int32 begin, Int32 size) {
1592 for (Idx i = begin; i < (begin + size); ++i) {
1593 Idx row_begin = A.ptr[i];
1594 Idx row_end = A.ptr[i + 1];
1595 Idx row_width = row_end - row_begin;
1597 Val D = math::zero<Val>();
1598 for (Idx j = row_begin; j < row_end; ++j) {
1604 else if (!aggr.strong_connection[j]) {
1611 Af.ptr[i + 1] = row_width;
1615 Af.set_nonzeros(Af.scan_row_sizes());
1617 arccoreParallelFor(0, Af.nbRow(), ForLoopRunInfo{}, [&](Int32 begin, Int32 size) {
1618 for (Idx i = begin; i < (begin + size); ++i) {
1620 Idx row_begin = A.ptr[i];
1621 Idx row_end = A.ptr[i + 1];
1622 Idx row_head = Af.ptr[i];
1624 for (Idx j = row_begin; j < row_end; ++j) {
1628 Af.col[row_head] = i;
1629 Af.val[row_head] = dia[i];
1632 else if (aggr.strong_connection[j]) {
1633 Af.col[row_head] = c;
1634 Af.val[row_head] = A.val[j];
1641 std::vector<Val> omega;
1643 auto P = interpolation(Af, dia, *P_tent, omega);
1644 auto R = restriction(Af, dia, *P_tent, omega);
1645 ARCCORE_ALINA_TOC(
"interpolation");
1647 return std::make_tuple(P, R);
1650 template <
class Matrix>
1651 std::shared_ptr<Matrix>
1652 coarse_operator(
const Matrix& A,
const Matrix& P,
const Matrix& R)
const
1654 return detail::galerkin(A, P, R);
1659 template <
class AMatrix,
typename Val,
typename Col,
typename Ptr>
1660 static std::shared_ptr<CSRMatrix<Val, Col, Ptr>>
1661 interpolation(
const AMatrix& A,
const std::vector<Val>& Adia,
1662 const CSRMatrix<Val, Col, Ptr>& P_tent,
1663 std::vector<Val>& omega)
1665 const size_t n = backend::nbRow(P_tent);
1666 const size_t nc = backend::nbColumn(P_tent);
1668 auto AP = product(A, P_tent,
true);
1670 omega.resize(nc, math::zero<Val>());
1671 std::vector<Val> denum(nc, math::zero<Val>());
1677 arccoreParallelFor(0, n, ForLoopRunInfo{}, [&](Int32 begin, Int32 size) {
1678 std::vector<ptrdiff_t> marker(nc, -1);
1683 std::vector<Col> adap_col(128);
1684 std::vector<Val> adap_val(128);
1686 for (ptrdiff_t ia = begin; ia < (begin + size); ++ia) {
1691 for (
auto a = A.row_begin(ia); a; ++a) {
1693 Val va = math::inverse(Adia[ca]) * a.value();
1695 for (
auto p = AP->row_begin(ca); p; ++p) {
1697 Val v = va * p.value();
1699 if (marker[c] < 0) {
1700 marker[c] = adap_col.size();
1701 adap_col.push_back(c);
1702 adap_val.push_back(v);
1705 adap_val[marker[c]] += v;
1710 detail::sort_row(&adap_col[0], &adap_val[0], adap_col.size());
1715 Ptr ja = AP->ptr[ia], ea = AP->ptr[ia + 1],
1716 jb = 0, eb = adap_col.size();
1717 ja < ea && jb < eb;) {
1718 Col ca = AP->col[ja];
1719 Col cb = adap_col[jb];
1726 Val v = AP->val[ja] * adap_val[jb];
1728 std::scoped_lock lock(mutex1);
1737 for (
size_t j = 0, e = adap_col.size(); j < e; ++j) {
1738 Col c = adap_col[j];
1739 Val v = adap_val[j];
1741 std::scoped_lock lock(mutex2);
1749 for (
size_t i = 0, m = omega.size(); i < m; ++i)
1750 omega[i] = math::inverse(denum[i]) * omega[i];
1759 arccoreParallelFor(0, n, ForLoopRunInfo{}, [&](Int32 begin, Int32 size) {
1760 for (ptrdiff_t i = begin; i < (begin + size); ++i) {
1761 Val dia = math::inverse(Adia[i]);
1763 for (Ptr ja = AP->ptr[i], ea = AP->ptr[i + 1],
1764 jp = P_tent.ptr[i], ep = P_tent.ptr[i + 1];
1766 Col ca = AP->col[ja];
1767 Val va = -dia * AP->val[ja] * omega[ca];
1769 for (; jp < ep; ++jp) {
1770 Col cp = P_tent.col[jp];
1775 va += P_tent.val[jp];
1788 template <
typename AMatrix,
typename Val,
typename Col,
typename Ptr>
1789 static std::shared_ptr<CSRMatrix<Val, Col, Ptr>>
1790 restriction(
const AMatrix& A,
const std::vector<Val>& Adia,
1791 const CSRMatrix<Val, Col, Ptr>& P_tent,
1792 const std::vector<Val>& omega)
1794 const size_t nc = backend::nbColumn(P_tent);
1799 auto RA = product(*R_tent, A,
true);
1808 arccoreParallelFor(0, nc, ForLoopRunInfo{}, [&](Int32 begin, Int32 size) {
1809 for (ptrdiff_t i = begin; i < (begin + size); ++i) {
1812 for (Ptr ja = RA->ptr[i], ea = RA->ptr[i + 1],
1813 jr = R_tent->ptr[i], er = R_tent->ptr[i + 1];
1815 Col ca = RA->col[ja];
1816 Val va = -w * math::inverse(Adia[ca]) * RA->val[ja];
1818 for (; jr < er; ++jr) {
1819 Col cr = R_tent->col[jr];
1824 va += R_tent->val[jr];
1846namespace Arcane::Alina::backend
1849template <
class Backend>
1851 typename std::enable_if<!std::is_arithmetic<typename backend::value_type<Backend>::type>::value>
::type> : std::false_type
std::vector< char > strong_connection
Strong connectivity matrix.
std::vector< ptrdiff_t > id
Aggerate id that each fine-level variable belongs to.
size_t count
Number of aggregates.
pointwise_aggregates(const Matrix &A, const params &prm, unsigned min_aggregate)
Constructs aggregates for a given matrix.
Loop execution information.
Matrix class, to be used by user.
Vector class, to be used by user.
__host__ __device__ Real dot(Real2 u, Real2 v)
Dot product of u by v in .
__host__ __device__ double sqrt(double v)
Square root of v.
__host__ __device__ Real3x3 transpose(const Real3x3 &t)
Transpose the matrix.
void arccoreParallelFor(const ComplexForLoopRanges< RankValue > &loop_ranges, const ForLoopRunInfo &run_info, const LambdaType &lambda_function, const ReducerArgs &... reducer_args)
Applies the lambda function lambda_function concurrently over the iteration interval given by loop_ra...
std::int32_t Int32
Signed integer type of 32 bits.
float over_interp
Over-interpolation factor .
nullspace_params nullspace
Near nullspace parameters.
Aggregates::params aggr
Aggregation parameters.
std::tuple< std::shared_ptr< Matrix >, std::shared_ptr< Matrix > > transfer_operators(const Matrix &A)
Creates transfer operators for the given system matrix.
std::shared_ptr< Matrix > coarse_operator(const Matrix &A, const Matrix &P, const Matrix &R) const
Creates system matrix for the coarser level.
Apply a scalar coarsening on a block matrix.
bool do_trunc
Truncate prolongation operator?
float eps_trunc
Truncation parameter .
float eps_strong
Parameter defining strong couplings.
Classic Ruge-Stuben coarsening with direct interpolation.
float relax
Relaxation factor.
Aggregates::params aggr
Aggregation parameters.
nullspace_params nullspace
Near nullspace parameters.
Aggregates::params aggr
Aggregation parameters.
nullspace_params nullspace
Near nullspace parameters.
Is the coarsening supported by the backend?
Number of rows for statically sized matrix types.
std::vector< double > B
Near nullspace vectors.
int cols
Number of vectors in near nullspace.
float eps_strong
Parameter defining strong couplings.
std::vector< ptrdiff_t > id
Aggerate id that each fine-level variable belongs to.
std::vector< char > strong_connection
Strong connectivity matrix.
plain_aggregates(const Matrix &A, const params &prm)
Constructs aggregates for a given matrix.
size_t count
Number of aggregates.
Int32 block_size
Block size for the system matrix.