Arcane  v4.1.10.0
Documentation développeur
Chargement...
Recherche...
Aucune correspondance
Adapters.h
1// -*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
2//-----------------------------------------------------------------------------
3// Copyright 2000-2026 CEA (www.cea.fr) IFPEN (www.ifpenergiesnouvelles.com)
4// See the top-level COPYRIGHT file for details.
5// SPDX-License-Identifier: Apache-2.0
6//-----------------------------------------------------------------------------
7/*---------------------------------------------------------------------------*/
8/* Adapters.h (C) 2000-2026 */
9/* */
10/*---------------------------------------------------------------------------*/
11#ifndef ARCCORE_ALINA_ADAPTERS_H
12#define ARCCORE_ALINA_ADAPTERS_H
13/*---------------------------------------------------------------------------*/
14/*---------------------------------------------------------------------------*/
15
16/*
17 * This file is based on the work on AMGCL library (version march 2026)
18 * which can be found at https://github.com/ddemidov/amgcl.
19 *
20 * Copyright (c) 2012-2022 Denis Demidov <dennis.demidov@gmail.com>
21 * SPDX-License-Identifier: MIT
22 */
23
24/*---------------------------------------------------------------------------*/
25/*---------------------------------------------------------------------------*/
26
27#include <type_traits>
28#include <vector>
29#include <tuple>
30
31#include <boost/iterator/permutation_iterator.hpp>
32
33#include "arccore/alina/AlinaUtils.h"
34#include "arccore/alina/BuiltinBackend.h"
35#include "arccore/alina/ValueTypeInterface.h"
36#include "arccore/alina/MatrixOperationsImpl.h"
37#include "arccore/alina/CuthillMcKeeReorderer.h"
38
39#include <span>
40
41/*---------------------------------------------------------------------------*/
42/*---------------------------------------------------------------------------*/
43
44namespace Arcane::Alina::backend
45{
46
47//---------------------------------------------------------------------------
48// Specialization of matrix interface
49//---------------------------------------------------------------------------
50template <typename N, typename PRng, typename CRng, typename VRng>
51struct value_type<std::tuple<N, PRng, CRng, VRng>>
52{
53 typedef std::decay_t<decltype(std::declval<VRng>()[0])> type;
54};
55
56template <typename N, typename PRng, typename CRng, typename VRng>
57struct rows_impl<std::tuple<N, PRng, CRng, VRng>>
58{
59 static size_t get(const std::tuple<N, PRng, CRng, VRng>& A)
60 {
61 return std::get<0>(A);
62 }
63};
64
65template <typename N, typename PRng, typename CRng, typename VRng>
66struct cols_impl<std::tuple<N, PRng, CRng, VRng>>
67{
68 static size_t get(const std::tuple<N, PRng, CRng, VRng>& A)
69 {
70 return std::get<0>(A);
71 }
72};
73
74template <typename N, typename PRng, typename CRng, typename VRng>
75struct nonzeros_impl<std::tuple<N, PRng, CRng, VRng>>
76{
77 static size_t get(const std::tuple<N, PRng, CRng, VRng>& A)
78 {
79 return std::get<1>(A)[std::get<0>(A)];
80 }
81};
82
83template <typename N, typename PRng, typename CRng, typename VRng>
84struct row_iterator<std::tuple<N, PRng, CRng, VRng>>
85{
86 class type
87 {
88 public:
89
90 typedef std::decay_t<decltype(std::declval<CRng>()[0])> col_type;
91 typedef std::decay_t<decltype(std::declval<VRng>()[0])> val_type;
92
93 type(const std::tuple<N, PRng, CRng, VRng>& A, size_t row)
94 : m_col(std::begin(std::get<2>(A)))
95 , m_end(std::begin(std::get<2>(A)))
96 , m_val(std::begin(std::get<3>(A)))
97 {
98 typedef std::decay_t<decltype(std::declval<PRng>()[0])> ptr_type;
99
100 ptr_type row_begin = std::get<1>(A)[row];
101 ptr_type row_end = std::get<1>(A)[row + 1];
102
103 m_col += row_begin;
104 m_end += row_end;
105 m_val += row_begin;
106 }
107
108 operator bool() const
109 {
110 return m_col != m_end;
111 }
112
113 type& operator++()
114 {
115 ++m_col;
116 ++m_val;
117 return *this;
118 }
119
120 col_type col() const
121 {
122 return *m_col;
123 }
124
125 val_type value() const
126 {
127 return *m_val;
128 }
129
130 private:
131
132 typedef decltype(std::begin(std::declval<VRng>())) val_iterator;
133 typedef decltype(std::begin(std::declval<CRng>())) col_iterator;
134
135 col_iterator m_col;
136 col_iterator m_end;
137 val_iterator m_val;
138 };
139};
140
141template <typename N, typename PRng, typename CRng, typename VRng>
142struct row_begin_impl<std::tuple<N, PRng, CRng, VRng>>
143{
144 typedef std::tuple<N, PRng, CRng, VRng> Matrix;
145 static typename row_iterator<Matrix>::type
146 get(const Matrix& matrix, size_t row)
147 {
148 return typename row_iterator<Matrix>::type(matrix, row);
149 }
150};
151
152template <typename N, typename PRng, typename CRng, typename VRng>
153struct row_nonzeros_impl<std::tuple<N, PRng, CRng, VRng>>
154{
155 typedef std::tuple<N, PRng, CRng, VRng> Matrix;
156
157 static size_t get(const Matrix& A, size_t row)
158 {
159 return std::get<1>(A)[row + 1] - std::get<1>(A)[row];
160 }
161};
162
163template <typename N, typename PRng, typename CRng, typename VRng>
164struct ptr_data_impl<std::tuple<N, PRng, CRng, VRng>>
165{
166 typedef std::tuple<N, PRng, CRng, VRng> Matrix;
167 typedef std::decay_t<decltype(std::declval<PRng>()[0])> ptr_type;
168 typedef const ptr_type* type;
169 static type get(const Matrix& A)
170 {
171 return &std::get<1>(A)[0];
172 }
173};
174
175template <typename N, typename PRng, typename CRng, typename VRng>
176struct col_data_impl<std::tuple<N, PRng, CRng, VRng>>
177{
178 typedef std::tuple<N, PRng, CRng, VRng> Matrix;
179 typedef std::decay_t<decltype(std::declval<CRng>()[0])> col_type;
180 typedef const col_type* type;
181 static type get(const Matrix& A)
182 {
183 return &std::get<2>(A)[0];
184 }
185};
186
187template <typename N, typename PRng, typename CRng, typename VRng>
188struct val_data_impl<std::tuple<N, PRng, CRng, VRng>>
189{
190 typedef std::tuple<N, PRng, CRng, VRng> Matrix;
191 typedef std::decay_t<decltype(std::declval<VRng>()[0])> val_type;
192 typedef const val_type* type;
193 static type get(const Matrix& A)
194 {
195 return &std::get<3>(A)[0];
196 }
197};
198
199/*---------------------------------------------------------------------------*/
200/*---------------------------------------------------------------------------*/
201
202} // namespace Arcane::Alina::backend
203
204/*---------------------------------------------------------------------------*/
205/*---------------------------------------------------------------------------*/
206
207namespace Arcane::Alina::adapter
208{
209
210/*---------------------------------------------------------------------------*/
211/*---------------------------------------------------------------------------*/
212
213template <class Matrix, class BlockType>
214struct block_matrix_adapter
215{
216 typedef BlockType value_type;
217 static const int BlockSize = math::static_rows<BlockType>::value;
218
219 const Matrix& A;
220
221 block_matrix_adapter(const Matrix& A)
222 : A(A)
223 {
224 precondition(
225 backend::nbRow(A) % BlockSize == 0 &&
226 backend::nbColumn(A) % BlockSize == 0,
227 "Matrix size is not divisible by block size!");
228 }
229
230 size_t rows() const
231 {
232 return backend::nbRow(A) / BlockSize;
233 }
234
235 size_t cols() const
236 {
237 return backend::nbColumn(A) / BlockSize;
238 }
239
240 size_t nonzeros() const
241 {
242 // Just an estimate:
243 return backend::nonzeros(A) / (BlockSize * BlockSize);
244 }
245
246 struct row_iterator
247 {
248 typedef typename backend::row_iterator<Matrix>::type Base;
249 typedef ptrdiff_t col_type;
250 typedef BlockType val_type;
251
252 std::array<char, sizeof(Base) * BlockSize> buf;
253 Base* base;
254
255 bool done;
256 col_type cur_col;
257 val_type cur_val;
258
259 row_iterator(const Matrix& A, col_type row)
260 : done(true)
261 {
262 base = reinterpret_cast<Base*>(buf.data());
263 for (int i = 0; i < BlockSize; ++i) {
264 new (base + i) Base(backend::row_begin(A, row * BlockSize + i));
265
266 if (base[i]) {
267 col_type col = base[i].col() / BlockSize;
268 if (done) {
269 cur_col = col;
270 done = false;
271 }
272 else {
273 cur_col = std::min<col_type>(cur_col, col);
274 }
275 }
276 }
277
278 if (done)
279 return;
280
281 // While we are gathering the current value,
282 // base iteratirs are advanced to the next block-column.
283 cur_val = math::zero<val_type>();
284 col_type end = (cur_col + 1) * BlockSize;
285 for (int i = 0; i < BlockSize; ++i) {
286 for (; base[i] && static_cast<ptrdiff_t>(base[i].col()) < end; ++base[i]) {
287 cur_val(i, base[i].col() % BlockSize) = base[i].value();
288 }
289 }
290 }
291
292 ~row_iterator()
293 {
294 for (int i = 0; i < BlockSize; ++i)
295 base[i].~Base();
296 }
297
298 operator bool() const
299 {
300 return !done;
301 }
302
303 row_iterator& operator++()
304 {
305 // Base iterators are already at the next block-column.
306 // We just need to gather the current column and value.
307 done = true;
308
309 col_type end = (cur_col + 1) * BlockSize;
310 for (int i = 0; i < BlockSize; ++i) {
311 if (base[i]) {
312 col_type col = base[i].col() / BlockSize;
313 if (done) {
314 cur_col = col;
315 done = false;
316 }
317 else {
318 cur_col = std::min<col_type>(cur_col, col);
319 }
320 }
321 }
322
323 if (done)
324 return *this;
325
326 cur_val = math::zero<val_type>();
327 end = (cur_col + 1) * BlockSize;
328 for (int i = 0; i < BlockSize; ++i) {
329 for (; base[i] && static_cast<ptrdiff_t>(base[i].col()) < end; ++base[i]) {
330 cur_val(i, base[i].col() % BlockSize) = base[i].value();
331 }
332 }
333
334 return *this;
335 }
336
337 col_type col() const
338 {
339 return cur_col;
340 }
341
342 val_type value() const
343 {
344 return cur_val;
345 }
346 };
347
348 row_iterator row_begin(size_t i) const
349 {
350 return row_iterator(A, i);
351 }
352};
353
354/*---------------------------------------------------------------------------*/
355/*---------------------------------------------------------------------------*/
356
358template <class BlockType, class Matrix>
359block_matrix_adapter<Matrix, BlockType> block_matrix(const Matrix& A)
360{
361 return block_matrix_adapter<Matrix, BlockType>(A);
362}
363
364/*---------------------------------------------------------------------------*/
365/*---------------------------------------------------------------------------*/
366
367template <class Matrix>
368std::shared_ptr<CSRMatrix<typename math::element_of<
369 typename backend::value_type<Matrix>::type>::type,
370 typename backend::col_type<Matrix>::type,
371 typename backend::ptr_type<Matrix>::type>>
372unblock_matrix(const Matrix& B)
373{
374 typedef typename backend::value_type<Matrix>::type Block;
375 typedef typename math::element_of<Block>::type Scalar;
376 typedef typename backend::col_type<Matrix>::type Col;
377 typedef typename backend::ptr_type<Matrix>::type Ptr;
378
379 const int brows = math::static_rows<Block>::value;
380 const int bcols = math::static_cols<Block>::value;
381
382 static_assert(brows > 1 || bcols > 1, "Can not unblock scalar matrix!");
383
384 auto A = std::make_shared<CSRMatrix<Scalar, Col, Ptr>>();
385
386 A->set_size(backend::nbRow(B) * brows, backend::nbColumn(B) * bcols);
387 A->ptr[0] = 0;
388
389 const ptrdiff_t nb = backend::nbRow(B);
390
391 arccoreParallelFor(0, nb, ForLoopRunInfo{}, [&](Int32 begin, Int32 size) {
392 for (ptrdiff_t ib = begin; ib < (begin + size); ++ib) {
393 auto w = backend::row_nonzeros(B, ib);
394 for (ptrdiff_t i = 0, ia = ib * brows; i < brows; ++i, ++ia) {
395 A->ptr[ia + 1] = w * bcols;
396 }
397 }
398 });
399
400 A->scan_row_sizes();
401 A->set_nonzeros();
402
403 arccoreParallelFor(0, nb, ForLoopRunInfo{}, [&](Int32 begin, Int32 size) {
404 for (ptrdiff_t ib = begin; ib < (begin + size); ++ib) {
405 for (auto b = backend::row_begin(B, ib); b; ++b) {
406 auto c = b.col();
407 auto v = b.value();
408
409 for (ptrdiff_t i = 0, ia = ib * brows; i < brows; ++i, ++ia) {
410 auto row_head = A->ptr[ia];
411 for (int j = 0; j < bcols; ++j) {
412 A->col[row_head] = c * bcols + j;
413 A->val[row_head] = v(i, j);
414 ++row_head;
415 }
416 A->ptr[ia] = row_head;
417 }
418 }
419 }
420 });
421
422 std::rotate(A->ptr.data(), A->ptr.data() + A->nbRow(), A->ptr.data() + A->nbRow() + 1);
423 A->ptr[0] = 0;
424
425 return A;
426}
427
428/*---------------------------------------------------------------------------*/
429/*---------------------------------------------------------------------------*/
430
431template <class Matrix>
432struct complex_adapter
433{
435 "value type should be complex");
436
438
439 const Matrix& A;
440
441 complex_adapter(const Matrix& A)
442 : A(A)
443 {}
444
445 size_t rows() const
446 {
447 return 2 * backend::nbRow(A);
448 }
449
450 size_t cols() const
451 {
452 return 2 * backend::nbColumn(A);
453 }
454
455 size_t nonzeros() const
456 {
457 return 4 * backend::nonzeros(A);
458 }
459
460 struct row_iterator
461 {
462 typedef typename backend::row_iterator<Matrix>::type Base;
463 typedef typename Base::col_type col_type;
464
465 row_iterator(const Base& base, bool row_real)
466 : base(base)
467 , row_real(row_real)
468 , col_real(true)
469 {}
470
471 operator bool() const
472 {
473 return static_cast<bool>(base);
474 }
475
476 row_iterator& operator++()
477 {
478 col_real = !col_real;
479 if (col_real)
480 ++base;
481
482 return *this;
483 }
484
485 col_type col() const
486 {
487 if (col_real)
488 return base.col() * 2;
489 else
490 return base.col() * 2 + 1;
491 }
492
493 value_type value() const
494 {
495 if (row_real) {
496 if (col_real)
497 return std::real(base.value());
498 else
499 return -std::imag(base.value());
500 }
501 else {
502 if (col_real)
503 return std::imag(base.value());
504 else
505 return std::real(base.value());
506 }
507 }
508
509 private:
510
511 Base base;
512 bool row_real;
513 bool col_real;
514 };
515
516 row_iterator row_begin(size_t i) const
517 {
518 return row_iterator(backend::row_begin(A, i / 2), i % 2 == 0);
519 }
520};
521
522template <class Matrix>
523complex_adapter<Matrix> complex_matrix(const Matrix& A)
524{
525 return complex_adapter<Matrix>(A);
526}
527
528template <class DataType, class Range>
529auto complex_range(Range& rng) -> Span<DataType>
530{
531 DataType* b = reinterpret_cast<DataType*>(&rng[0]);
532 size_t s = 2 * std::size(rng);
533
534 return Span<DataType>(b, s);
535}
536
537/*---------------------------------------------------------------------------*/
538/*---------------------------------------------------------------------------*/
544template <class RowBuilder>
545struct matrix_builder
546{
547 typedef typename RowBuilder::val_type value_type;
548 typedef typename RowBuilder::col_type col_type;
549
550 RowBuilder build_row;
551
552 matrix_builder(const RowBuilder& row_builder)
553 : build_row(row_builder)
554 {}
555
556 size_t rows() const { return build_row.rows(); }
557 size_t cols() const { return build_row.rows(); }
558 size_t nonzeros() const { return build_row.nonzeros(); }
559
560 struct row_iterator
561 {
562 typedef RowBuilder::col_type col_type;
563 typedef RowBuilder::val_type val_type;
564
565 typedef std::vector<col_type>::const_iterator col_iterator;
566 typedef std::vector<val_type>::const_iterator val_iterator;
567
568 row_iterator(const RowBuilder& build_row, size_t i)
569 : ptr(0)
570 {
571 build_row(i, m_col, m_val);
572 }
573
574 operator bool() const
575 {
576 return m_col.size() - ptr;
577 }
578
579 row_iterator& operator++()
580 {
581 ++ptr;
582 return *this;
583 }
584
585 col_type col() const
586 {
587 return m_col[ptr];
588 }
589
590 val_type value() const
591 {
592 return m_val[ptr];
593 }
594
595 private:
596
597 int ptr;
598 std::vector<col_type> m_col;
599 std::vector<value_type> m_val;
600 };
601
602 row_iterator row_begin(size_t i) const
603 {
604 return row_iterator(build_row, i);
605 }
606};
607
608/*---------------------------------------------------------------------------*/
609/*---------------------------------------------------------------------------*/
610
612template <class RowBuilder>
613matrix_builder<RowBuilder> make_matrix(const RowBuilder& row_builder)
614{
615 return matrix_builder<RowBuilder>(row_builder);
616}
617
618/*---------------------------------------------------------------------------*/
619/*---------------------------------------------------------------------------*/
620
621template <class Matrix>
622struct reordered_matrix
623{
624 typedef backend::value_type<Matrix>::type value_type;
625 typedef backend::row_iterator<Matrix>::type base_iterator;
626
627 const Matrix& A;
628 const ptrdiff_t* perm;
629 const ptrdiff_t* iperm;
630
631 reordered_matrix(const Matrix& A, const ptrdiff_t* perm, const ptrdiff_t* iperm)
632 : A(A)
633 , perm(perm)
634 , iperm(iperm)
635 {}
636
637 size_t rows() const
638 {
639 return backend::nbRow(A);
640 }
641
642 size_t cols() const
643 {
644 return backend::nbColumn(A);
645 }
646
647 size_t nonzeros() const
648 {
649 return backend::nonzeros(A);
650 }
651
652 struct row_iterator
653 {
654 base_iterator base;
655 const ptrdiff_t* iperm;
656
657 row_iterator(const base_iterator& base, const ptrdiff_t* iperm)
658 : base(base)
659 , iperm(iperm)
660 {}
661
662 operator bool() const
663 {
664 return base;
665 }
666
667 row_iterator& operator++()
668 {
669 ++base;
670 return *this;
671 }
672
673 ptrdiff_t col() const
674 {
675 return iperm[base.col()];
676 }
677
678 value_type value() const
679 {
680 return base.value();
681 }
682 };
683
684 row_iterator row_begin(size_t i) const
685 {
686 return row_iterator(backend::row_begin(A, perm[i]), iperm);
687 }
688};
689
690/*---------------------------------------------------------------------------*/
691/*---------------------------------------------------------------------------*/
692
693template <class Vector>
694struct reordered_vector
695{
696 typedef backend::value_type<std::decay_t<Vector>>::type raw_value_type;
697 typedef std::conditional_t<std::is_const_v<Vector>, const raw_value_type, raw_value_type> value_type;
698
699 Vector& x;
700 const ptrdiff_t* perm;
701
702 reordered_vector(Vector& x, const ptrdiff_t* perm)
703 : x(x)
704 , perm(perm)
705 {}
706
707 size_t size() const
708 {
709 return std::size(x);
710 }
711
712 value_type& operator[](size_t i) const
713 {
714 return x[perm[i]];
715 }
716
717 boost::permutation_iterator<typename std::decay_t<Vector>::iterator, const ptrdiff_t*>
718 begin()
719 {
720 return boost::make_permutation_iterator(std::begin(x), perm);
721 }
722
723 boost::permutation_iterator<typename std::decay_t<Vector>::const_iterator, const ptrdiff_t*>
724 begin() const
725 {
726 return boost::make_permutation_iterator(std::begin(x), perm);
727 }
728
729 boost::permutation_iterator<typename std::decay_t<Vector>::iterator, const ptrdiff_t*>
730 end()
731 {
732 return boost::make_permutation_iterator(std::end(x), perm + size());
733 }
734
735 boost::permutation_iterator<typename std::decay_t<Vector>::const_iterator, const ptrdiff_t*>
736 end() const
737 {
738 return boost::make_permutation_iterator(std::end(x), perm + size());
739 }
740};
741
742/*---------------------------------------------------------------------------*/
743/*---------------------------------------------------------------------------*/
744
745template <class ordering = CuthillMcKeeReorderer<false>>
746class reorder
747{
748 public:
749
750 template <class Matrix>
751 explicit reorder(const Matrix& A)
752 : n(backend::nbRow(A))
753 , perm(n)
754 , iperm(n)
755 {
756 ordering::get(A, perm);
757 arccoreParallelFor(0, n, ForLoopRunInfo{}, [&](Int32 begin, Int32 size) {
758 for (ptrdiff_t i = begin; i < (begin + size); ++i) {
759 iperm[perm[i]] = i;
760 }
761 });
762 }
763
764 template <class Matrix>
765 std::enable_if_t<!backend::is_builtin_vector<Matrix>::value, reordered_matrix<Matrix>>
766 operator()(const Matrix& A) const
767 {
768 return reordered_matrix<Matrix>(A, perm.data(), iperm.data());
769 }
770
771 template <class Vector>
772 std::enable_if_t<backend::is_builtin_vector<Vector>::value, reordered_vector<Vector>>
773 operator()(Vector& x) const
774 {
775 return reordered_vector<Vector>(x, perm.data());
776 }
777
778 template <class Vector>
779 std::enable_if_t<backend::is_builtin_vector<Vector>::value, reordered_vector<const Vector>>
780 operator()(const Vector& x) const
781 {
782 return reordered_vector<const Vector>(x, perm.data());
783 }
784
785 template <class Vector1, class Vector2>
786 void forward(const Vector1& x, Vector2& y) const
787 {
788 arccoreParallelFor(0, n, ForLoopRunInfo{}, [&](Int32 begin, Int32 size) {
789 for (ptrdiff_t i = begin; i < (begin + size); ++i) {
790 y[i] = x[perm[i]];
791 }
792 });
793 }
794
795 template <class Vector1, class Vector2>
796 void inverse(const Vector1& x, Vector2& y) const
797 {
798 arccoreParallelFor(0, n, ForLoopRunInfo{}, [&](Int32 begin, Int32 size) {
799 for (ptrdiff_t i = begin; i < (begin + size); ++i) {
800 y[perm[i]] = x[i];
801 }
802 });
803 }
804
805 private:
806
807 ptrdiff_t n;
808 numa_vector<ptrdiff_t> perm, iperm;
809};
810
811/*---------------------------------------------------------------------------*/
812/*---------------------------------------------------------------------------*/
813
814template <class Matrix, class Scale>
815struct scaled_matrix
816{
817 typedef typename backend::value_type<Matrix>::type value_type;
818 typedef typename backend::value_type<Scale>::type scale_type;
819
820 const Matrix& A;
821 const Scale& s;
822
823 scaled_matrix(const Matrix& A, const Scale& s)
824 : A(A)
825 , s(s)
826 {}
827
828 size_t rows() const { return backend::nbRow(A); }
829 size_t cols() const { return backend::nbColumn(A); }
830 size_t nonzeros() const { return backend::nonzeros(A); }
831
832 struct row_iterator : public backend::row_iterator<Matrix>::type
833 {
834 typedef typename backend::row_iterator<Matrix>::type Base;
835
836 scale_type si;
837 const Scale& s;
838
839 row_iterator(const Matrix& A, const Scale& s, size_t i)
840 : Base(A, i)
841 , si(s[i])
842 , s(s)
843 {}
844
845 value_type value() const
846 {
847 return si * static_cast<const Base*>(this)->value() * s[this->col()];
848 }
849 };
850
851 row_iterator row_begin(size_t i) const
852 {
853 return row_iterator(A, s, i);
854 }
855};
856
857/*---------------------------------------------------------------------------*/
858/*---------------------------------------------------------------------------*/
859
860template <class Backend, class Scale>
861struct scaled_problem
862{
863 typedef typename Backend::params backend_params;
864
865 const std::shared_ptr<Scale> s;
866 const backend_params& bprm;
867
868 scaled_problem(std::shared_ptr<Scale> s, const backend_params& bprm = backend_params())
869 : s(s)
870 , bprm(bprm)
871 {}
872
873 template <class Matrix>
874 scaled_matrix<Matrix, Scale> matrix(const Matrix& A) const
875 {
876 return scaled_matrix<Matrix, Scale>(A, *s);
877 }
878
879 template <class Vector>
880 std::shared_ptr<typename Backend::vector> rhs(const Vector& v) const
881 {
882 auto t = Backend::copy_vector(v, bprm);
883 (*this)(*t);
884 return t;
885 }
886
887 template <class Vector>
888 void operator()(Vector& x) const
889 {
890 typedef typename backend::value_type<Vector>::type value_type;
891 typedef typename math::scalar_of<value_type>::type scalar_type;
892
893 const auto one = math::identity<scalar_type>();
894 const auto zero = math::zero<scalar_type>();
895
897 backend::vmul(one, *s, x, zero, x);
898 }
899 else {
900 backend::vmul(one, *Backend::copy_vector(*s, bprm), x, zero, x);
901 }
902 }
903};
904
905/*---------------------------------------------------------------------------*/
906/*---------------------------------------------------------------------------*/
907
908template <class Backend, class Matrix>
909scaled_problem<Backend,
910 std::vector<
911 typename math::scalar_of<
912 typename backend::value_type<Matrix>::type>::type>>
913scale_diagonal(const Matrix& A,
914 const typename Backend::params& bprm = typename Backend::params())
915{
916 typedef typename backend::value_type<Matrix>::type value_type;
917 typedef typename math::scalar_of<value_type>::type scalar_type;
918 ptrdiff_t n = backend::nbRow(A);
919 auto s = std::make_shared<std::vector<scalar_type>>(n);
920
921 arccoreParallelFor(0, n, ForLoopRunInfo{}, [&](Int32 begin, Int32 size) {
922 for (ptrdiff_t i = begin; i < (begin + size); ++i) {
923 for (auto a = backend::row_begin(A, i); a; ++a) {
924 if (a.col() == i) {
925 (*s)[i] = math::inverse(sqrt(math::norm(a.value())));
926 break;
927 }
928 }
929 }
930 });
931
933}
934
935/*---------------------------------------------------------------------------*/
936/*---------------------------------------------------------------------------*/
937
938template <typename Ptr, typename Col, typename Val>
939std::shared_ptr<CSRMatrix<Val>>
940zero_copy(size_t nrows, size_t ncols, const Ptr* ptr, const Col* col, const Val* val)
941{
942 // Check that Ptr and Col types are binary-compatible with ptrdiff_t:
943 static_assert(std::is_integral_v<Ptr>, "Unsupported Ptr type");
944 static_assert(std::is_integral_v<Col>, "Unsupported Col type");
945 static_assert(sizeof(Ptr) == sizeof(ptrdiff_t), "Unsupported Ptr type");
946 static_assert(sizeof(Col) == sizeof(ptrdiff_t), "Unsupported Col type");
947
948 auto A = std::make_shared<CSRMatrix<Val>>();
949 A->setNbRow(nrows);
950 A->ncols = ncols;
951 A->setNbNonZero(nrows ? ptr[nrows] : 0);
952
953 A->ptr.setPointerZeroCopy((ptrdiff_t*)ptr);
954 A->col.setPointerZeroCopy((ptrdiff_t*)col);
955 A->val.setPointerZeroCopy((Val*)val);
956
957 A->own_data = false;
958
959 return A;
960}
961
962/*---------------------------------------------------------------------------*/
963/*---------------------------------------------------------------------------*/
964
965template <typename Ptr, typename Col, typename Val>
966std::shared_ptr<CSRMatrix<Val>>
967zero_copy(size_t n, const Ptr* ptr, const Col* col, const Val* val)
968{
969 return zero_copy(n, n, ptr, col, val);
970}
971
972/*---------------------------------------------------------------------------*/
973/*---------------------------------------------------------------------------*/
974
975template <typename Ptr, typename Col, typename Val>
976std::shared_ptr<CSRMatrix<Val, Col, Ptr>>
977zero_copy_direct(size_t nrows, size_t ncols, const Ptr* ptr, const Col* col, const Val* val)
978{
979 auto A = std::make_shared<CSRMatrix<Val, Col, Ptr>>();
980 A->setNbRow(nrows);
981 A->ncols = ncols;
982 A->setNbNonZero(nrows ? ptr[nrows] : 0);
983
984 A->ptr.setPointerZeroCopy(const_cast<Ptr*>(ptr));
985 A->col.setPointerZeroCopy(const_cast<Col*>(col));
986 A->val.setPointerZeroCopy(const_cast<Val*>(val));
987
988 A->own_data = false;
989
990 return A;
991}
992
993/*---------------------------------------------------------------------------*/
994/*---------------------------------------------------------------------------*/
995
996template <typename Ptr, typename Col, typename Val>
997std::shared_ptr<CSRMatrix<Val, Col, Ptr>>
998zero_copy_direct(size_t n, const Ptr* ptr, const Col* col, const Val* val)
999{
1000 return zero_copy_direct(n, n, ptr, col, val);
1001}
1002
1003/*---------------------------------------------------------------------------*/
1004/*---------------------------------------------------------------------------*/
1005
1006} // namespace Arcane::Alina::adapter
1007
1008/*---------------------------------------------------------------------------*/
1009/*---------------------------------------------------------------------------*/
1010
1011namespace Arcane::Alina::backend
1012{
1013template <class Vector>
1014struct is_builtin_vector<adapter::reordered_vector<Vector>>
1015: is_builtin_vector<typename std::decay<Vector>::type>
1016{};
1017} // namespace Arcane::Alina::backend
1018
1019namespace Arcane::Alina::backend::detail
1020{
1021
1022template <class Matrix, class BlockType>
1023struct use_builtin_matrix_ops<adapter::block_matrix_adapter<Matrix, BlockType>>
1024: std::true_type
1025{};
1026
1027template <class Matrix>
1029: std::true_type
1030{};
1031
1032template <class RowBuilder>
1034: std::true_type
1035{};
1036
1037template <typename N, typename PRng, typename CRng, typename VRng>
1038struct use_builtin_matrix_ops<std::tuple<N, PRng, CRng, VRng>>
1039: std::true_type
1040{};
1041
1042template <class Matrix>
1043struct use_builtin_matrix_ops<adapter::reordered_matrix<Matrix>>
1044: std::true_type
1045{};
1046
1047/*---------------------------------------------------------------------------*/
1048/*---------------------------------------------------------------------------*/
1049
1050} // namespace Arcane::Alina::backend::detail
1051
1052/*---------------------------------------------------------------------------*/
1053/*---------------------------------------------------------------------------*/
1054
1055#endif
NUMA-aware vector container.
Definition NumaVector.h:42
Informations d'exécution d'une boucle.
Matrix class, to be used by user.
Vue d'un tableau d'éléments de type T.
Definition Span.h:633
Classe gérant un vecteur de dimension 2 de type T.
Definition Vector2.h:36
Vector class, to be used by user.
void arccoreParallelFor(const ComplexForLoopRanges< RankValue > &loop_ranges, const ForLoopRunInfo &run_info, const LambdaType &lambda_function, const ReducerArgs &... reducer_args)
Applique en concurrence la fonction lambda lambda_function sur l'intervalle d'itération donné par loo...
Definition ParallelFor.h:85
std::int32_t Int32
Type entier signé sur 32 bits.
Generates matrix rows as needed with help of user-provided functor.
Definition Adapters.h:546
Implementation for function returning the number of columns in a matrix.
Implementation for function returning the number of nonzeros in a matrix.
Metafunction that returns pointer type of a matrix.
Implementation for function returning row iterator for a matrix.
Metafunction returning the row iterator type for a matrix type.
Implementation for function returning the number of nonzeros in a matrix row.
Implementation for function returning the number of rows in a matrix.
Metafunction that returns value type of a matrix or a vector type.
Scalar type of a non-scalar type.
Number of rows for statically sized matrix types.