11#ifndef ARCCORE_ALINA_ADAPTERS_H
12#define ARCCORE_ALINA_ADAPTERS_H
31#include <boost/iterator/permutation_iterator.hpp>
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"
44namespace Arcane::Alina::backend
50template <
typename N,
typename PRng,
typename CRng,
typename VRng>
53 typedef std::decay_t<decltype(std::declval<VRng>()[0])> type;
56template <
typename N,
typename PRng,
typename CRng,
typename VRng>
59 static size_t get(
const std::tuple<N, PRng, CRng, VRng>& A)
61 return std::get<0>(A);
65template <
typename N,
typename PRng,
typename CRng,
typename VRng>
68 static size_t get(
const std::tuple<N, PRng, CRng, VRng>& A)
70 return std::get<0>(A);
74template <
typename N,
typename PRng,
typename CRng,
typename VRng>
77 static size_t get(
const std::tuple<N, PRng, CRng, VRng>& A)
79 return std::get<1>(A)[std::get<0>(A)];
83template <
typename N,
typename PRng,
typename CRng,
typename VRng>
90 typedef std::decay_t<decltype(std::declval<CRng>()[0])> col_type;
91 typedef std::decay_t<decltype(std::declval<VRng>()[0])> val_type;
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)))
98 typedef std::decay_t<decltype(std::declval<PRng>()[0])>
ptr_type;
100 ptr_type row_begin = std::get<1>(A)[row];
101 ptr_type row_end = std::get<1>(A)[row + 1];
108 operator bool()
const
110 return m_col != m_end;
125 val_type value()
const
132 typedef decltype(std::begin(std::declval<VRng>())) val_iterator;
133 typedef decltype(std::begin(std::declval<CRng>())) col_iterator;
141template <
typename N,
typename PRng,
typename CRng,
typename VRng>
144 typedef std::tuple<N, PRng, CRng, VRng> Matrix;
145 static typename row_iterator<Matrix>::type
146 get(
const Matrix&
matrix,
size_t row)
148 return typename row_iterator<Matrix>::type(
matrix, row);
152template <
typename N,
typename PRng,
typename CRng,
typename VRng>
155 typedef std::tuple<N, PRng, CRng, VRng> Matrix;
157 static size_t get(
const Matrix& A,
size_t row)
159 return std::get<1>(A)[row + 1] - std::get<1>(A)[row];
163template <
typename N,
typename PRng,
typename CRng,
typename VRng>
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)
171 return &std::get<1>(A)[0];
175template <
typename N,
typename PRng,
typename CRng,
typename VRng>
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)
183 return &std::get<2>(A)[0];
187template <
typename N,
typename PRng,
typename CRng,
typename VRng>
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)
195 return &std::get<3>(A)[0];
207namespace Arcane::Alina::adapter
213template <
class Matrix,
class BlockType>
214struct block_matrix_adapter
216 typedef BlockType value_type;
221 block_matrix_adapter(
const Matrix& A)
225 backend::nbRow(A) % BlockSize == 0 &&
226 backend::nbColumn(A) % BlockSize == 0,
227 "Matrix size is not divisible by block size!");
232 return backend::nbRow(A) / BlockSize;
237 return backend::nbColumn(A) / BlockSize;
240 size_t nonzeros()
const
243 return backend::nonzeros(A) / (BlockSize * BlockSize);
248 typedef typename backend::row_iterator<Matrix>::type Base;
249 typedef ptrdiff_t col_type;
250 typedef BlockType val_type;
252 std::array<char,
sizeof(Base) * BlockSize> buf;
259 row_iterator(
const Matrix& A, col_type row)
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));
267 col_type col = base[i].col() / BlockSize;
273 cur_col = std::min<col_type>(cur_col, col);
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();
294 for (
int i = 0; i < BlockSize; ++i)
298 operator bool()
const
303 row_iterator& operator++()
309 col_type end = (cur_col + 1) * BlockSize;
310 for (
int i = 0; i < BlockSize; ++i) {
312 col_type col = base[i].col() / BlockSize;
318 cur_col = std::min<col_type>(cur_col, col);
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();
342 val_type value()
const
358template <
class BlockType,
class Matrix>
359block_matrix_adapter<Matrix, BlockType> block_matrix(
const Matrix& A)
361 return block_matrix_adapter<Matrix, BlockType>(A);
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)
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;
379 const int brows = math::static_rows<Block>::value;
380 const int bcols = math::static_cols<Block>::value;
382 static_assert(brows > 1 || bcols > 1,
"Can not unblock scalar matrix!");
384 auto A = std::make_shared<CSRMatrix<Scalar, Col, Ptr>>();
386 A->set_size(backend::nbRow(B) * brows, backend::nbColumn(B) * bcols);
389 const ptrdiff_t nb = backend::nbRow(B);
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;
404 for (ptrdiff_t ib = begin; ib < (begin + size); ++ib) {
405 for (
auto b = backend::row_begin(B, ib); b; ++b) {
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);
416 A->ptr[ia] = row_head;
422 std::rotate(A->ptr.data(), A->ptr.data() + A->nbRow(), A->ptr.data() + A->nbRow() + 1);
431template <
class Matrix>
432struct complex_adapter
435 "value type should be complex");
441 complex_adapter(
const Matrix& A)
447 return 2 * backend::nbRow(A);
452 return 2 * backend::nbColumn(A);
455 size_t nonzeros()
const
457 return 4 * backend::nonzeros(A);
462 typedef typename backend::row_iterator<Matrix>::type Base;
463 typedef typename Base::col_type col_type;
465 row_iterator(
const Base& base,
bool row_real)
471 operator bool()
const
473 return static_cast<bool>(base);
476 row_iterator& operator++()
478 col_real = !col_real;
488 return base.col() * 2;
490 return base.col() * 2 + 1;
493 value_type value()
const
497 return std::real(base.value());
499 return -std::imag(base.value());
503 return std::imag(base.value());
505 return std::real(base.value());
518 return row_iterator(backend::row_begin(A, i / 2), i % 2 == 0);
522template <
class Matrix>
523complex_adapter<Matrix> complex_matrix(
const Matrix& A)
525 return complex_adapter<Matrix>(A);
528template <
class DataType,
class Range>
531 DataType* b =
reinterpret_cast<DataType*
>(&rng[0]);
532 size_t s = 2 * std::size(rng);
544template <
class RowBuilder>
547 typedef typename RowBuilder::val_type value_type;
548 typedef typename RowBuilder::col_type col_type;
550 RowBuilder build_row;
552 matrix_builder(
const RowBuilder& row_builder)
553 : build_row(row_builder)
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(); }
562 typedef RowBuilder::col_type col_type;
563 typedef RowBuilder::val_type val_type;
565 typedef std::vector<col_type>::const_iterator col_iterator;
566 typedef std::vector<val_type>::const_iterator val_iterator;
568 row_iterator(
const RowBuilder& build_row,
size_t i)
571 build_row(i, m_col, m_val);
574 operator bool()
const
576 return m_col.size() - ptr;
579 row_iterator& operator++()
590 val_type value()
const
598 std::vector<col_type> m_col;
599 std::vector<value_type> m_val;
612template <
class RowBuilder>
613matrix_builder<RowBuilder> make_matrix(
const RowBuilder& row_builder)
615 return matrix_builder<RowBuilder>(row_builder);
621template <
class Matrix>
622struct reordered_matrix
624 typedef backend::value_type<Matrix>::type value_type;
625 typedef backend::row_iterator<Matrix>::type base_iterator;
628 const ptrdiff_t* perm;
629 const ptrdiff_t* iperm;
631 reordered_matrix(
const Matrix& A,
const ptrdiff_t* perm,
const ptrdiff_t* iperm)
639 return backend::nbRow(A);
644 return backend::nbColumn(A);
647 size_t nonzeros()
const
649 return backend::nonzeros(A);
655 const ptrdiff_t* iperm;
657 row_iterator(
const base_iterator& base,
const ptrdiff_t* iperm)
662 operator bool()
const
667 row_iterator& operator++()
673 ptrdiff_t col()
const
675 return iperm[base.col()];
678 value_type value()
const
686 return row_iterator(backend::row_begin(A, perm[i]), iperm);
693template <
class Vector>
694struct reordered_vector
697 typedef std::conditional_t<std::is_const_v<Vector>,
const raw_value_type, raw_value_type> value_type;
700 const ptrdiff_t* perm;
702 reordered_vector(
Vector& x,
const ptrdiff_t* perm)
712 value_type& operator[](
size_t i)
const
717 boost::permutation_iterator<typename std::decay_t<Vector>::iterator,
const ptrdiff_t*>
720 return boost::make_permutation_iterator(std::begin(x), perm);
723 boost::permutation_iterator<typename std::decay_t<Vector>::const_iterator,
const ptrdiff_t*>
726 return boost::make_permutation_iterator(std::begin(x), perm);
729 boost::permutation_iterator<typename std::decay_t<Vector>::iterator,
const ptrdiff_t*>
732 return boost::make_permutation_iterator(std::end(x), perm + size());
735 boost::permutation_iterator<typename std::decay_t<Vector>::const_iterator,
const ptrdiff_t*>
738 return boost::make_permutation_iterator(std::end(x), perm + size());
745template <
class ordering = CuthillMcKeeReorderer<false>>
750 template <
class Matrix>
751 explicit reorder(
const Matrix& A)
752 : n(backend::nbRow(A))
756 ordering::get(A, perm);
758 for (ptrdiff_t i = begin; i < (begin + size); ++i) {
764 template <
class Matrix>
766 operator()(
const Matrix& A)
const
771 template <
class Vector>
773 operator()(
Vector& x)
const
778 template <
class Vector>
780 operator()(
const Vector& x)
const
785 template <
class Vector1,
class Vector2>
786 void forward(
const Vector1& x,
Vector2& y)
const
789 for (ptrdiff_t i = begin; i < (begin + size); ++i) {
795 template <
class Vector1,
class Vector2>
796 void inverse(
const Vector1& x,
Vector2& y)
const
799 for (ptrdiff_t i = begin; i < (begin + size); ++i) {
814template <
class Matrix,
class Scale>
817 typedef typename backend::value_type<Matrix>::type value_type;
818 typedef typename backend::value_type<Scale>::type scale_type;
823 scaled_matrix(
const Matrix& A,
const Scale& s)
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); }
834 typedef typename backend::row_iterator<Matrix>::type Base;
839 row_iterator(
const Matrix& A,
const Scale& s,
size_t i)
845 value_type value()
const
847 return si *
static_cast<const Base*
>(
this)->value() * s[this->col()];
860template <
class Backend,
class Scale>
863 typedef typename Backend::params backend_params;
865 const std::shared_ptr<Scale> s;
866 const backend_params& bprm;
868 scaled_problem(std::shared_ptr<Scale> s,
const backend_params& bprm = backend_params())
873 template <
class Matrix>
879 template <
class Vector>
880 std::shared_ptr<typename Backend::vector> rhs(
const Vector& v)
const
882 auto t = Backend::copy_vector(v, bprm);
887 template <
class Vector>
888 void operator()(
Vector& x)
const
890 typedef typename backend::value_type<Vector>::type value_type;
891 typedef typename math::scalar_of<value_type>::type scalar_type;
893 const auto one = math::identity<scalar_type>();
894 const auto zero = math::zero<scalar_type>();
897 backend::vmul(one, *s, x, zero, x);
900 backend::vmul(one, *Backend::copy_vector(*s, bprm), x, zero, x);
908template <
class Backend,
class Matrix>
912 typename backend::value_type<Matrix>::type>::type>>
913scale_diagonal(
const Matrix& A,
914 const typename Backend::params& bprm =
typename Backend::params())
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);
922 for (ptrdiff_t i = begin; i < (begin + size); ++i) {
923 for (
auto a = backend::row_begin(A, i); a; ++a) {
925 (*s)[i] = math::inverse(sqrt(math::norm(a.value())));
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)
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");
948 auto A = std::make_shared<CSRMatrix<Val>>();
951 A->setNbNonZero(nrows ? ptr[nrows] : 0);
953 A->ptr.setPointerZeroCopy((ptrdiff_t*)ptr);
954 A->col.setPointerZeroCopy((ptrdiff_t*)col);
955 A->val.setPointerZeroCopy((Val*)val);
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)
969 return zero_copy(n, n, ptr, col, val);
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)
979 auto A = std::make_shared<CSRMatrix<Val, Col, Ptr>>();
982 A->setNbNonZero(nrows ? ptr[nrows] : 0);
984 A->ptr.setPointerZeroCopy(
const_cast<Ptr*
>(ptr));
985 A->col.setPointerZeroCopy(
const_cast<Col*
>(col));
986 A->val.setPointerZeroCopy(
const_cast<Val*
>(val));
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)
1000 return zero_copy_direct(n, n, ptr, col, val);
1011namespace Arcane::Alina::backend
1013template <
class Vector>
1019namespace Arcane::Alina::backend::detail
1022template <
class Matrix,
class BlockType>
1027template <
class Matrix>
1032template <
class RowBuilder>
1037template <
typename N,
typename PRng,
typename CRng,
typename VRng>
1042template <
class Matrix>
NUMA-aware vector container.
Informations d'exécution d'une boucle.
Matrix class, to be used by user.
Vue d'un tableau d'éléments de type T.
Classe gérant un vecteur de dimension 2 de type T.
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...
std::int32_t Int32
Type entier signé sur 32 bits.
Generates matrix rows as needed with help of user-provided functor.
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.