44class DistributedCPRPreconditioner
46 static_assert(std::is_same<
47 typename PPrecond::BackendType,
48 typename SPrecond::BackendType>::value,
49 "Backends for pressure and flow preconditioners should coinside!");
53 using BackendType = PPrecond::BackendType;
55 typedef BackendType::value_type value_type;
56 typedef math::scalar_of<value_type>::type scalar_type;
57 typedef BackendType::matrix bmatrix;
58 typedef BackendType::vector vector;
59 typedef BackendType::params backend_params;
62 typedef typename BuiltinBackend<value_type>::matrix build_matrix;
66 typedef typename PPrecond::params pprecond_params;
67 typedef typename SPrecond::params sprecond_params;
69 pprecond_params pprecond;
70 sprecond_params sprecond;
77 : ARCCORE_ALINA_PARAMS_IMPORT_CHILD(p, pprecond)
78 , ARCCORE_ALINA_PARAMS_IMPORT_CHILD(p, sprecond)
79 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, block_size)
81 p.check_params({
"pprecond",
"sprecond",
"block_size",
"active_rows" });
84 void get(
PropertyTree& p,
const std::string& path =
"")
const
86 ARCCORE_ALINA_PARAMS_EXPORT_CHILD(p, path, pprecond);
87 ARCCORE_ALINA_PARAMS_EXPORT_CHILD(p, path, sprecond);
88 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, block_size);
92 template <
class Matrix>
96 const backend_params& bprm = backend_params())
99 , n(backend::nbRow(K))
101 init(std::make_shared<matrix>(comm, K, backend::nbRow(K)), bprm);
105 std::shared_ptr<matrix> K,
107 const backend_params& bprm = backend_params())
115 template <
class Vec1,
class Vec2>
116 void apply(
const Vec1& rhs, Vec2&& x)
const
118 const auto one = math::identity<scalar_type>();
119 const auto zero = math::zero<scalar_type>();
122 backend::residual(rhs, S->system_matrix(), x, *rs);
124 backend::spmv(one, *Fpp, *rs, zero, *rp);
127 backend::spmv(one, *Scatter, *xp, one, x);
130 std::shared_ptr<matrix> system_matrix_ptr()
const
132 return S->system_matrix_ptr();
135 const matrix& system_matrix()
const
137 return S->system_matrix();
146 mpi_communicator comm;
149 std::shared_ptr<PPrecond> P;
150 std::shared_ptr<SPrecond> S;
152 std::shared_ptr<bmatrix> Fpp, Scatter;
153 std::shared_ptr<vector> rp, xp, rs;
155 void init(std::shared_ptr<matrix> K,
const backend_params& bprm)
157 typedef typename backend::row_iterator<build_matrix>::type row_iterator;
158 const int B = prm.block_size;
160 auto _K_loc = K->local();
161 auto _K_rem = K->remote();
163 auto& K_loc = *_K_loc;
164 auto& K_rem = *_K_rem;
168 auto fpp = std::make_shared<build_matrix>();
169 fpp->set_size(np, n);
170 fpp->set_nonzeros(n);
173 auto App_loc = std::make_shared<build_matrix>();
174 auto App_rem = std::make_shared<build_matrix>();
176 App_loc->set_size(np, np,
true);
177 App_rem->set_size(np, 0,
true);
182 std::vector<row_iterator> k;
184 multi_array<value_type, 2> v(B, B);
186 for (ptrdiff_t ip = begin; ip < (begin + size); ++ip) {
187 ptrdiff_t ik = ip * B;
189 ptrdiff_t cur_col = 0;
193 for (
int i = 0; i < B; ++i) {
194 k.push_back(backend::row_begin(K_loc, ik + i));
197 ptrdiff_t col = k.back().col() / B;
203 cur_col = std::min(cur_col, col);
207 fpp->col[ik + i] = ik + i;
209 fpp->ptr[ip + 1] = ik + B;
212 ++App_loc->ptr[ip + 1];
214 ptrdiff_t end = (cur_col + 1) * B;
220 for (
int i = 0; i < B; ++i)
221 for (
int j = 0; j < B; ++j)
224 for (
int i = 0; i < B; ++i)
225 for (; k[i] && k[i].col() < end; ++k[i])
226 v(k[i].col() % B, i) = k[i].value();
228 invert(v, &fpp->val[ik]);
233 for (
int i = 0; i < B; ++i)
234 while (k[i] && k[i].col() < end)
240 for (
int i = 0; i < B; ++i) {
242 ptrdiff_t col = k[i].col() / B;
248 cur_col = std::min(cur_col, col);
256 for (
int i = 0; i < B; ++i) {
257 k.push_back(backend::row_begin(K_rem, ik + i));
260 ptrdiff_t col = k.back().col() / B;
266 cur_col = std::min(cur_col, col);
272 ++App_rem->ptr[ip + 1];
274 ptrdiff_t end = (cur_col + 1) * B;
276 for (
int i = 0; i < B; ++i)
277 while (k[i] && k[i].col() < end)
282 for (
int i = 0; i < B; ++i) {
284 ptrdiff_t col = k[i].col() / B;
290 cur_col = std::min(cur_col, col);
298 App_loc->set_nonzeros(App_loc->scan_row_sizes());
299 App_rem->set_nonzeros(App_rem->scan_row_sizes());
301 auto scatter = std::make_shared<build_matrix>();
302 scatter->set_size(n, np);
303 scatter->set_nonzeros(np);
307 std::vector<row_iterator> k;
310 for (ptrdiff_t ip = begin; ip < (begin + size); ++ip) {
311 ptrdiff_t ik = ip * B;
315 value_type* d = &fpp->val[ik];
318 ptrdiff_t head = App_loc->ptr[ip];
320 for (
int i = 0; i < B; ++i) {
321 k.push_back(backend::row_begin(K_loc, ik + i));
324 ptrdiff_t col = k.back().col() / B;
330 cur_col = std::min(cur_col, col);
336 ptrdiff_t end = (cur_col + 1) * B;
339 for (
int i = 0; i < B; ++i) {
340 for (; k[i] && k[i].col() < end; ++k[i]) {
341 if (k[i].col() % B == 0) {
342 app += d[i] * k[i].value();
347 App_loc->col[head] = cur_col;
348 App_loc->val[head] = app;
353 for (
int i = 0; i < B; ++i) {
355 ptrdiff_t col = k[i].col() / B;
361 cur_col = std::min(cur_col, col);
368 head = App_rem->ptr[ip];
370 for (
int i = 0; i < B; ++i) {
371 k.push_back(backend::row_begin(K_rem, ik + i));
374 ptrdiff_t col = k.back().col() / B;
380 cur_col = std::min(cur_col, col);
386 ptrdiff_t end = (cur_col + 1) * B;
389 for (
int i = 0; i < B; ++i) {
390 for (; k[i] && k[i].col() < end; ++k[i]) {
391 if (k[i].col() % B == 0) {
392 app += d[i] * k[i].value();
397 App_rem->col[head] = cur_col;
398 App_rem->val[head] = app;
403 for (
int i = 0; i < B; ++i) {
405 ptrdiff_t col = k[i].col() / B;
411 cur_col = std::min(cur_col, col);
417 scatter->col[ip] = ip;
418 scatter->val[ip] = math::identity<value_type>();
421 for (
int i = 0; i < B; ++i) {
424 scatter->ptr[ik + i + 1] = nnz;
429 auto App = std::make_shared<matrix>(comm, App_loc, App_rem);
431 P = std::make_shared<PPrecond>(comm, App, prm.pprecond, bprm);
432 S = std::make_shared<SPrecond>(comm, K, prm.sprecond, bprm);
434 Fpp = BackendType::copy_matrix(fpp, bprm);
435 Scatter = BackendType::copy_matrix(scatter, bprm);
437 rp = BackendType::create_vector(np, bprm);
438 xp = BackendType::create_vector(np, bprm);
439 rs = BackendType::create_vector(n, bprm);
444 void invert(multi_array<value_type, 2>& A, value_type* y)
446 const int B = prm.block_size;
449 for (
int k = 0; k < B; ++k) {
450 value_type d = A(k, k);
451 assert(!math::is_zero(d));
452 for (
int i = k + 1; i < B; ++i) {
454 for (
int j = k + 1; j < B; ++j)
455 A(i, j) -= A(i, k) * A(k, j);
461 for (
int i = 0; i < B; ++i) {
462 value_type b =
static_cast<value_type
>(i == 0);
463 for (
int j = 0; j < i; ++j)
469 for (
int i = B; i-- > 0;) {
470 for (
int j = i + 1; j < B; ++j)
471 y[i] -= A(i, j) * y[j];
476 template <
class P,
class S>
477 friend std::ostream& operator<<(std::ostream& os,
const DistributedCPRPreconditioner<P, S>& cpr)
479 os <<
"CPR (two-stage preconditioner)\n"
480 "### Pressure preconditioner:\n"
482 "### Global preconditioner:\n"
483 << *cpr.S << std::endl;