DOLFINx 0.10.0.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
MatrixCSR.h
1// Copyright (C) 2021-2022 Garth N. Wells and Chris N. Richardson
2//
3// This file is part of DOLFINx (https://www.fenicsproject.org)
4//
5// SPDX-License-Identifier: LGPL-3.0-or-later
6
7#pragma once
8
9#include "SparsityPattern.h"
10#include "Vector.h"
11#include "matrix_csr_impl.h"
12#include <algorithm>
13#include <dolfinx/common/IndexMap.h>
14#include <dolfinx/common/MPI.h>
15#include <dolfinx/graph/AdjacencyList.h>
16#include <mpi.h>
17#include <numeric>
18#include <span>
19#include <utility>
20#include <vector>
21
22namespace dolfinx::la
23{
25enum class BlockMode : int
26{
27 compact = 0,
32};
33
48template <typename Scalar, typename Container = std::vector<Scalar>,
49 typename ColContainer = std::vector<std::int32_t>,
50 typename RowPtrContainer = std::vector<std::int64_t>>
52{
53public:
55 using value_type = Scalar;
56
58 using container_type = Container;
59
61 using column_container_type = ColContainer;
62
64 using rowptr_container_type = RowPtrContainer;
65
66 static_assert(std::is_same_v<value_type, typename container_type::value_type>,
67 "Scalar type and container value type must be the same.");
68 static_assert(std::is_integral_v<typename column_container_type::value_type>);
69 static_assert(std::is_integral_v<typename rowptr_container_type::value_type>);
70
94 template <int BS0 = 1, int BS1 = 1>
96 {
97 if ((BS0 != _bs[0] and BS0 > 1 and _bs[0] > 1)
98 or (BS1 != _bs[1] and BS1 > 1 and _bs[1] > 1))
99 {
100 throw std::runtime_error(
101 "Cannot insert blocks of different size than matrix block size");
102 }
103
104 return [&](std::span<const std::int32_t> rows,
105 std::span<const std::int32_t> cols,
106 std::span<const value_type> data) -> int
107 {
108 this->set<BS0, BS1>(data, rows, cols);
109 return 0;
110 };
111 }
112
136 template <int BS0 = 1, int BS1 = 1>
138 {
139 if ((BS0 != _bs[0] and BS0 > 1 and _bs[0] > 1)
140 or (BS1 != _bs[1] and BS1 > 1 and _bs[1] > 1))
141 {
142 throw std::runtime_error(
143 "Cannot insert blocks of different size than matrix block size");
144 }
145
146 return [&](std::span<const std::int32_t> rows,
147 std::span<const std::int32_t> cols,
148 std::span<const value_type> data) -> int
149 {
150 this->add<BS0, BS1>(data, rows, cols);
151 return 0;
152 };
153 }
154
178 MatrixCSR(const SparsityPattern& p, BlockMode mode = BlockMode::compact);
179
182 MatrixCSR(MatrixCSR&& A) = default;
183
186 MatrixCSR(const MatrixCSR& A) = default;
187
190 template <typename Mat>
191 // requires requires(Mat A) {
192 // //TODO, see examples in la::Vector
193 // }
194 explicit MatrixCSR(const Mat& A)
195 : _index_maps(A.index_maps()), _block_mode(A.block_mode()),
196 _bs(A.block_size()), _data(A.values()), _cols(A.cols()),
197 _row_ptr(A.row_ptr()), _off_diagonal_offset(A.off_diag_offset()),
198 _comm(A.comm())
199 {
200 }
201
205 void set(value_type x) { std::ranges::fill(_data, x); }
206
223 template <int BS0, int BS1>
224 void set(std::span<const value_type> x, std::span<const std::int32_t> rows,
225 std::span<const std::int32_t> cols)
226 {
227 auto set_fn = [](value_type& y, const value_type& x) { y = x; };
228
229 std::int32_t num_rows
230 = _index_maps[0]->size_local() + _index_maps[0]->num_ghosts();
231 assert(x.size() == rows.size() * cols.size() * BS0 * BS1);
232 if (_bs[0] == BS0 and _bs[1] == BS1)
233 {
234 impl::insert_csr<BS0, BS1>(_data, _cols, _row_ptr, x, rows, cols, set_fn,
235 num_rows);
236 }
237 else if (_bs[0] == 1 and _bs[1] == 1)
238 {
239 // Set blocked data in a regular CSR matrix (_bs[0]=1, _bs[1]=1)
240 // with correct sparsity
241 impl::insert_blocked_csr<BS0, BS1>(_data, _cols, _row_ptr, x, rows, cols,
242 set_fn, num_rows);
243 }
244 else
245 {
246 assert(BS0 == 1 and BS1 == 1);
247 // Set non-blocked data in a blocked CSR matrix (BS0=1, BS1=1)
248 impl::insert_nonblocked_csr(_data, _cols, _row_ptr, x, rows, cols, set_fn,
249 num_rows, _bs[0], _bs[1]);
250 }
251 }
252
268 template <int BS0 = 1, int BS1 = 1>
269 void add(std::span<const value_type> x, std::span<const std::int32_t> rows,
270 std::span<const std::int32_t> cols)
271 {
272 auto add_fn = [](value_type& y, const value_type& x) { y += x; };
273
274 assert(x.size() == rows.size() * cols.size() * BS0 * BS1);
275 if (_bs[0] == BS0 and _bs[1] == BS1)
276 {
277 impl::insert_csr<BS0, BS1>(_data, _cols, _row_ptr, x, rows, cols, add_fn,
278 _row_ptr.size());
279 }
280 else if (_bs[0] == 1 and _bs[1] == 1)
281 {
282 // Add blocked data to a regular CSR matrix (_bs[0]=1, _bs[1]=1)
283 impl::insert_blocked_csr<BS0, BS1>(_data, _cols, _row_ptr, x, rows, cols,
284 add_fn, _row_ptr.size());
285 }
286 else
287 {
288 assert(BS0 == 1 and BS1 == 1);
289 // Add non-blocked data to a blocked CSR matrix (BS0=1, BS1=1)
290 impl::insert_nonblocked_csr(_data, _cols, _row_ptr, x, rows, cols, add_fn,
291 _row_ptr.size(), _bs[0], _bs[1]);
292 }
293 }
294
296 std::int32_t num_owned_rows() const { return _index_maps[0]->size_local(); }
297
299 std::int32_t num_all_rows() const { return _row_ptr.size() - 1; }
300
310 std::vector<value_type> to_dense() const
311 {
312 const std::size_t nrows = num_all_rows();
313 const std::size_t ncols = _index_maps[1]->size_global();
314 std::vector<value_type> A(nrows * ncols * _bs[0] * _bs[1], 0.0);
315 for (std::size_t r = 0; r < nrows; ++r)
316 {
317 for (std::int32_t j = _row_ptr[r]; j < _row_ptr[r + 1]; ++j)
318 {
319 for (int i0 = 0; i0 < _bs[0]; ++i0)
320 {
321 for (int i1 = 0; i1 < _bs[1]; ++i1)
322 {
323 std::array<std::int32_t, 1> local_col{_cols[j]};
324 std::array<std::int64_t, 1> global_col{0};
325 _index_maps[1]->local_to_global(local_col, global_col);
326 A[(r * _bs[1] + i0) * ncols * _bs[0] + global_col[0] * _bs[1] + i1]
327 = _data[j * _bs[0] * _bs[1] + i0 * _bs[1] + i1];
328 }
329 }
330 }
331 }
332
333 return A;
334 }
335
343 {
346 }
347
356 {
357 const std::int32_t local_size0 = _index_maps[0]->size_local();
358 const std::int32_t num_ghosts0 = _index_maps[0]->num_ghosts();
359 const int bs2 = _bs[0] * _bs[1];
360
361 // For each ghost row, pack and send values to send to neighborhood
362 std::vector<int> insert_pos = _val_send_disp;
363 _ghost_value_data.resize(_val_send_disp.back());
364 for (int i = 0; i < num_ghosts0; ++i)
365 {
366 const int rank = _ghost_row_to_rank[i];
367
368 // Get position in send buffer to place data to send to this
369 // neighbour
370 const std::int32_t val_pos = insert_pos[rank];
371 std::copy(std::next(_data.data(), _row_ptr[local_size0 + i] * bs2),
372 std::next(_data.data(), _row_ptr[local_size0 + i + 1] * bs2),
373 std::next(_ghost_value_data.begin(), val_pos));
374 insert_pos[rank]
375 += bs2 * (_row_ptr[local_size0 + i + 1] - _row_ptr[local_size0 + i]);
376 }
377
378 _ghost_value_data_in.resize(_val_recv_disp.back());
379
380 // Compute data sizes for send and receive from displacements
381 std::vector<int> val_send_count(_val_send_disp.size() - 1);
382 std::adjacent_difference(std::next(_val_send_disp.begin()),
383 _val_send_disp.end(), val_send_count.begin());
384
385 std::vector<int> val_recv_count(_val_recv_disp.size() - 1);
386 std::adjacent_difference(std::next(_val_recv_disp.begin()),
387 _val_recv_disp.end(), val_recv_count.begin());
388
389 int status = MPI_Ineighbor_alltoallv(
390 _ghost_value_data.data(), val_send_count.data(), _val_send_disp.data(),
391 dolfinx::MPI::mpi_t<value_type>, _ghost_value_data_in.data(),
392 val_recv_count.data(), _val_recv_disp.data(),
393 dolfinx::MPI::mpi_t<value_type>, _comm.comm(), &_request);
394 dolfinx::MPI::check_error(_comm.comm(), status);
395 }
396
403 {
404 int status = MPI_Wait(&_request, MPI_STATUS_IGNORE);
405 dolfinx::MPI::check_error(_comm.comm(), status);
406
407 _ghost_value_data.clear();
408 _ghost_value_data.shrink_to_fit();
409
410 // Add to local rows
411 const int bs2 = _bs[0] * _bs[1];
412 assert(_ghost_value_data_in.size() == _unpack_pos.size() * bs2);
413 for (std::size_t i = 0; i < _unpack_pos.size(); ++i)
414 for (int j = 0; j < bs2; ++j)
415 _data[_unpack_pos[i] * bs2 + j] += _ghost_value_data_in[i * bs2 + j];
416
417 _ghost_value_data_in.clear();
418 _ghost_value_data_in.shrink_to_fit();
419
420 // Set ghost row data to zero
421 const std::int32_t local_size0 = _index_maps[0]->size_local();
422 std::fill(std::next(_data.begin(), _row_ptr[local_size0] * bs2),
423 _data.end(), 0);
424 }
425
428 double squared_norm() const
429 {
430 const std::size_t num_owned_rows = _index_maps[0]->size_local();
431 const int bs2 = _bs[0] * _bs[1];
432 assert(num_owned_rows < _row_ptr.size());
433 double norm_sq_local = std::accumulate(
434 _data.cbegin(),
435 std::next(_data.cbegin(), _row_ptr[num_owned_rows] * bs2), double(0),
436 [](auto norm, value_type y) { return norm + std::norm(y); });
437 double norm_sq;
438 MPI_Allreduce(&norm_sq_local, &norm_sq, 1, MPI_DOUBLE, MPI_SUM,
439 _comm.comm());
440 return norm_sq;
441 }
442
454
456 MPI_Comm comm() const { return _comm.comm(); }
457
465 const std::array<std::shared_ptr<const common::IndexMap>, 2>&
467 {
468 return _index_maps;
469 }
470
478 std::shared_ptr<const common::IndexMap> index_map(int dim) const
479 {
480 return _index_maps.at(dim);
481 }
482
485 container_type& values() { return _data; }
486
489 const container_type& values() const { return _data; }
490
493 const rowptr_container_type& row_ptr() const { return _row_ptr; }
494
497 const column_container_type& cols() const { return _cols; }
498
507 {
508 return _off_diagonal_offset;
509 }
510
513 std::array<int, 2> block_size() const { return _bs; }
514
517 BlockMode block_mode() const { return _block_mode; }
518
519private:
520 // Maps for the distribution of the rows and columns
521 std::array<std::shared_ptr<const common::IndexMap>, 2> _index_maps;
522
523 // Block mode (compact or expanded)
524 BlockMode _block_mode;
525
526 // Block sizes
527 std::array<int, 2> _bs;
528
529 // Matrix data
530 container_type _data;
532 rowptr_container_type _row_ptr;
533
534 // Start of off-diagonal (unowned columns) on each row
535 rowptr_container_type _off_diagonal_offset;
536
537 // Neighborhood communicator (ghost->owner communicator for rows)
538 dolfinx::MPI::Comm _comm;
539
540 // -- Precomputed data for scatter_rev/update
541
542 // Request in non-blocking communication
543 MPI_Request _request;
544
545 // Position in _data to add received data
546 std::vector<int> _unpack_pos;
547
548 // Displacements for alltoall for each neighbor when sending and
549 // receiving
550 std::vector<int> _val_send_disp, _val_recv_disp;
551
552 // Ownership of each row, by neighbor (for the neighbourhood defined
553 // on _comm)
554 std::vector<int> _ghost_row_to_rank;
555
556 // Temporary stores for data during non-blocking communication
557 container_type _ghost_value_data;
558 container_type _ghost_value_data_in;
559};
560//-----------------------------------------------------------------------------
561template <class U, class V, class W, class X>
563 : _index_maps({p.index_map(0),
564 std::make_shared<common::IndexMap>(p.column_index_map())}),
565 _block_mode(mode), _bs({p.block_size(0), p.block_size(1)}),
566 _data(p.num_nonzeros() * _bs[0] * _bs[1], 0),
567 _cols(p.graph().first.begin(), p.graph().first.end()),
568 _row_ptr(p.graph().second.begin(), p.graph().second.end()),
569 _comm(MPI_COMM_NULL)
570{
571 if (_block_mode == BlockMode::expanded)
572 {
573 // Rebuild IndexMaps
574 for (int i = 0; i < 2; ++i)
575 {
576 const auto im = _index_maps[i];
577 const std::int32_t size_local = im->size_local() * _bs[i];
578 std::span ghost_i = im->ghosts();
579 std::vector<std::int64_t> ghosts;
580 const std::vector<int> ghost_owner_i(im->owners().begin(),
581 im->owners().end());
582 std::vector<int> src_rank;
583 for (std::size_t j = 0; j < ghost_i.size(); ++j)
584 {
585 for (int k = 0; k < _bs[i]; ++k)
586 {
587 ghosts.push_back(ghost_i[j] * _bs[i] + k);
588 src_rank.push_back(ghost_owner_i[j]);
589 }
590 }
591
592 std::array<std::vector<int>, 2> src_dest0
593 = {std::vector(_index_maps[i]->src().begin(),
594 _index_maps[i]->src().end()),
595 std::vector(_index_maps[i]->dest().begin(),
596 _index_maps[i]->dest().end())};
597 _index_maps[i] = std::make_shared<common::IndexMap>(
598 _index_maps[i]->comm(), size_local, src_dest0, ghosts, src_rank);
599 }
600
601 // Convert sparsity pattern and set _bs to 1
602
603 column_container_type new_cols;
604 new_cols.reserve(_data.size());
605 rowptr_container_type new_row_ptr{0};
606 new_row_ptr.reserve(_row_ptr.size() * _bs[0]);
607 std::span<const std::int32_t> num_diag_nnz = p.off_diagonal_offsets();
608 for (std::size_t i = 0; i < _row_ptr.size() - 1; ++i)
609 {
610 // Repeat row _bs[0] times
611 for (int q0 = 0; q0 < _bs[0]; ++q0)
612 {
613 _off_diagonal_offset.push_back(new_row_ptr.back()
614 + num_diag_nnz[i] * _bs[1]);
615 for (auto j = _row_ptr[i]; j < _row_ptr[i + 1]; ++j)
616 {
617 for (int q1 = 0; q1 < _bs[1]; ++q1)
618 new_cols.push_back(_cols[j] * _bs[1] + q1);
619 }
620 new_row_ptr.push_back(new_cols.size());
621 }
622 }
623 _cols = new_cols;
624 _row_ptr = new_row_ptr;
625 _bs[0] = 1;
626 _bs[1] = 1;
627 }
628 else
629 {
630 // Compute off-diagonal offset for each row (compact)
631 std::span<const std::int32_t> num_diag_nnz = p.off_diagonal_offsets();
632 _off_diagonal_offset.reserve(num_diag_nnz.size());
633 std::ranges::transform(num_diag_nnz, _row_ptr,
634 std::back_inserter(_off_diagonal_offset),
635 std::plus{});
636 }
637
638 // Some short-hand
639 const std::array local_size
640 = {_index_maps[0]->size_local(), _index_maps[1]->size_local()};
641 const std::array local_range
642 = {_index_maps[0]->local_range(), _index_maps[1]->local_range()};
643 std::span ghosts1 = _index_maps[1]->ghosts();
644
645 std::span ghosts0 = _index_maps[0]->ghosts();
646 std::span src_ranks = _index_maps[0]->src();
647 std::span dest_ranks = _index_maps[0]->dest();
648
649 // Create neighbourhood communicator (owner <- ghost)
650 MPI_Comm comm;
651 MPI_Dist_graph_create_adjacent(_index_maps[0]->comm(), dest_ranks.size(),
652 dest_ranks.data(), MPI_UNWEIGHTED,
653 src_ranks.size(), src_ranks.data(),
654 MPI_UNWEIGHTED, MPI_INFO_NULL, false, &comm);
655 _comm = dolfinx::MPI::Comm(comm, false);
656
657 // Build map from ghost row index position to owning (neighborhood)
658 // rank
659 _ghost_row_to_rank.reserve(_index_maps[0]->owners().size());
660 for (int r : _index_maps[0]->owners())
661 {
662 auto it = std::ranges::lower_bound(src_ranks, r);
663 assert(it != src_ranks.end() and *it == r);
664 std::size_t pos = std::distance(src_ranks.begin(), it);
665 _ghost_row_to_rank.push_back(pos);
666 }
667
668 // Compute size of data to send to each neighbor
669 std::vector<std::int32_t> data_per_proc(src_ranks.size(), 0);
670 for (std::size_t i = 0; i < _ghost_row_to_rank.size(); ++i)
671 {
672 assert(_ghost_row_to_rank[i] < (int)data_per_proc.size());
673 std::size_t pos = local_size[0] + i;
674 data_per_proc[_ghost_row_to_rank[i]] += _row_ptr[pos + 1] - _row_ptr[pos];
675 }
676
677 // Compute send displacements
678 _val_send_disp.resize(src_ranks.size() + 1, 0);
679 std::partial_sum(data_per_proc.begin(), data_per_proc.end(),
680 std::next(_val_send_disp.begin()));
681
682 // For each ghost row, pack and send indices to neighborhood
683 std::vector<std::int64_t> ghost_index_data(2 * _val_send_disp.back());
684 {
685 std::vector<int> insert_pos = _val_send_disp;
686 for (std::size_t i = 0; i < _ghost_row_to_rank.size(); ++i)
687 {
688 const int rank = _ghost_row_to_rank[i];
689 std::int32_t row_id = local_size[0] + i;
690 for (int j = _row_ptr[row_id]; j < _row_ptr[row_id + 1]; ++j)
691 {
692 // Get position in send buffer
693 const std::int32_t idx_pos = 2 * insert_pos[rank];
694
695 // Pack send data (row, col) as global indices
696 ghost_index_data[idx_pos] = ghosts0[i];
697 if (std::int32_t col_local = _cols[j]; col_local < local_size[1])
698 ghost_index_data[idx_pos + 1] = col_local + local_range[1][0];
699 else
700 ghost_index_data[idx_pos + 1] = ghosts1[col_local - local_size[1]];
701
702 insert_pos[rank] += 1;
703 }
704 }
705 }
706
707 // Communicate data with neighborhood
708 std::vector<std::int64_t> ghost_index_array;
709 std::vector<int> recv_disp;
710 {
711 std::vector<int> send_sizes;
712 std::ranges::transform(data_per_proc, std::back_inserter(send_sizes),
713 [](auto x) { return 2 * x; });
714
715 std::vector<int> recv_sizes(dest_ranks.size());
716 send_sizes.reserve(1);
717 recv_sizes.reserve(1);
718 MPI_Neighbor_alltoall(send_sizes.data(), 1, MPI_INT, recv_sizes.data(), 1,
719 MPI_INT, _comm.comm());
720
721 // Build send/recv displacement
722 std::vector<int> send_disp{0};
723 std::partial_sum(send_sizes.begin(), send_sizes.end(),
724 std::back_inserter(send_disp));
725 recv_disp = {0};
726 std::partial_sum(recv_sizes.begin(), recv_sizes.end(),
727 std::back_inserter(recv_disp));
728
729 ghost_index_array.resize(recv_disp.back());
730 MPI_Neighbor_alltoallv(ghost_index_data.data(), send_sizes.data(),
731 send_disp.data(), MPI_INT64_T,
732 ghost_index_array.data(), recv_sizes.data(),
733 recv_disp.data(), MPI_INT64_T, _comm.comm());
734 }
735
736 // Store receive displacements for future use, when transferring
737 // data values
738 _val_recv_disp.resize(recv_disp.size());
739 const int bs2 = _bs[0] * _bs[1];
740 std::ranges::transform(recv_disp, _val_recv_disp.begin(),
741 [&bs2](auto d) { return bs2 * d / 2; });
742 std::ranges::transform(_val_send_disp, _val_send_disp.begin(),
743 [&bs2](auto d) { return d * bs2; });
744
745 // Global-to-local map for ghost columns
746 std::vector<std::pair<std::int64_t, std::int32_t>> global_to_local;
747 global_to_local.reserve(ghosts1.size());
748 for (std::int64_t idx : ghosts1)
749 global_to_local.push_back({idx, global_to_local.size() + local_size[1]});
750 std::ranges::sort(global_to_local);
751
752 // Compute location in which data for each index should be stored
753 // when received
754 for (std::size_t i = 0; i < ghost_index_array.size(); i += 2)
755 {
756 // Row must be on this process
757 const std::int32_t local_row = ghost_index_array[i] - local_range[0][0];
758 assert(local_row >= 0 and local_row < local_size[0]);
759
760 // Column may be owned or unowned
761 std::int32_t local_col = ghost_index_array[i + 1] - local_range[1][0];
762 if (local_col < 0 or local_col >= local_size[1])
763 {
764 auto it = std::ranges::lower_bound(
765 global_to_local, std::pair(ghost_index_array[i + 1], -1),
766 [](auto& a, auto& b) { return a.first < b.first; });
767 assert(it != global_to_local.end()
768 and it->first == ghost_index_array[i + 1]);
769 local_col = it->second;
770 }
771 auto cit0 = std::next(_cols.begin(), _row_ptr[local_row]);
772 auto cit1 = std::next(_cols.begin(), _row_ptr[local_row + 1]);
773
774 // Find position of column index and insert data
775 auto cit = std::lower_bound(cit0, cit1, local_col);
776 assert(cit != cit1);
777 assert(*cit == local_col);
778 std::size_t d = std::distance(_cols.begin(), cit);
779 _unpack_pos.push_back(d);
780 }
781}
782//-----------------------------------------------------------------------------
783
784// The matrix A is distributed across P processes by blocks of rows:
785// A = | A_0 |
786// | A_1 |
787// | ... |
788// | A_P-1 |
789//
790// Each submatrix A_i is owned by a single process "i" and can be further
791// decomposed into diagonal (Ai[0]) and off diagonal (Ai[1]) blocks:
792// Ai = |Ai[0] Ai[1]|
793//
794// If A is square, the diagonal block Ai[0] is also square and contains
795// only owned columns and rows. The block Ai[1] contains ghost columns
796// (unowned dofs).
797
798// Likewise, a local vector x can be decomposed into owned and ghost blocks:
799// xi = | x[0] |
800// | x[1] |
801//
802// So the product y = Ax can be computed into two separate steps:
803// y[0] = |Ai[0] Ai[1]| | x[0] | = Ai[0] x[0] + Ai[1] x[1]
804// | x[1] |
805//
808template <typename Scalar, typename V, typename W, typename X>
811{
812 // start communication (update ghosts)
814
815 const std::int32_t nrowslocal = num_owned_rows();
816 std::span<const std::int64_t> Arow_ptr(row_ptr().data(), nrowslocal + 1);
817 std::span<const std::int32_t> Acols(cols().data(), Arow_ptr[nrowslocal]);
818 std::span<const std::int64_t> Aoff_diag_offset(off_diag_offset().data(),
819 nrowslocal);
820 std::span<const Scalar> Avalues(values().data(), Arow_ptr[nrowslocal]);
821
822 std::span<const Scalar> _x = x.array();
823 std::span<Scalar> _y = y.array();
824
825 std::span<const std::int64_t> Arow_begin(Arow_ptr.data(), nrowslocal);
826 std::span<const std::int64_t> Arow_end(Arow_ptr.data() + 1, nrowslocal);
827
828 // First stage: spmv - diagonal
829 // yi[0] += Ai[0] * xi[0]
830 if (_bs[1] == 1)
831 {
832 impl::spmv<Scalar, 1>(Avalues, Arow_begin, Aoff_diag_offset, Acols, _x, _y,
833 _bs[0], 1);
834 }
835 else
836 {
837 impl::spmv<Scalar, -1>(Avalues, Arow_begin, Aoff_diag_offset, Acols, _x, _y,
838 _bs[0], _bs[1]);
839 }
840
841 // finalize ghost update
842 x.scatter_fwd_end();
843
844 // Second stage: spmv - off-diagonal
845 // yi[0] += Ai[1] * xi[1]
846 if (_bs[1] == 1)
847 {
848 impl::spmv<Scalar, 1>(Avalues, Aoff_diag_offset, Arow_end, Acols, _x, _y,
849 _bs[0], 1);
850 }
851 else
852 {
853 impl::spmv<Scalar, -1>(Avalues, Aoff_diag_offset, Arow_end, Acols, _x, _y,
854 _bs[0], _bs[1]);
855 }
856}
857
858} // namespace dolfinx::la
A duplicate MPI communicator and manage lifetime of the communicator.
Definition MPI.h:42
const container_type & values() const
Definition MatrixCSR.h:489
std::shared_ptr< const common::IndexMap > index_map(int dim) const
Index map for the row or column space.
Definition MatrixCSR.h:478
const rowptr_container_type & off_diag_offset() const
Definition MatrixCSR.h:506
void set(std::span< const value_type > x, std::span< const std::int32_t > rows, std::span< const std::int32_t > cols)
Set values in the matrix.
Definition MatrixCSR.h:224
RowPtrContainer rowptr_container_type
Row pointer container type.
Definition MatrixCSR.h:64
void scatter_rev_end()
End transfer of ghost row data to owning ranks.
Definition MatrixCSR.h:402
container_type & values()
Definition MatrixCSR.h:485
auto mat_add_values()
Insertion functor for adding values to a matrix. It is typically used in finite element assembly func...
Definition MatrixCSR.h:137
BlockMode block_mode() const
Get 'block mode'.
Definition MatrixCSR.h:517
void add(std::span< const value_type > x, std::span< const std::int32_t > rows, std::span< const std::int32_t > cols)
Accumulate values in the matrix.
Definition MatrixCSR.h:269
std::int32_t num_owned_rows() const
Number of local rows excluding ghost rows.
Definition MatrixCSR.h:296
MatrixCSR(const SparsityPattern &p, BlockMode mode=BlockMode::compact)
Create a distributed matrix.
Definition MatrixCSR.h:562
ColContainer column_container_type
Column index container type.
Definition MatrixCSR.h:61
MatrixCSR(MatrixCSR &&A)=default
void mult(Vector< value_type > &x, Vector< value_type > &y)
Compute the product y += Ax.
Definition MatrixCSR.h:809
double squared_norm() const
Compute the Frobenius norm squared across all processes.
Definition MatrixCSR.h:428
void scatter_rev()
Transfer ghost row data to the owning ranks accumulating received values on the owned rows,...
Definition MatrixCSR.h:342
Container container_type
Matrix entries container type.
Definition MatrixCSR.h:58
Scalar value_type
Scalar type.
Definition MatrixCSR.h:55
void scatter_rev_begin()
Begin transfer of ghost row data to owning ranks, where it will be accumulated into existing owned ro...
Definition MatrixCSR.h:355
const column_container_type & cols() const
Definition MatrixCSR.h:497
void set(value_type x)
Set all non-zero local entries to a value including entries in ghost rows.
Definition MatrixCSR.h:205
std::array< int, 2 > block_size() const
Get block sizes.
Definition MatrixCSR.h:513
std::int32_t num_all_rows() const
Number of local rows including ghost rows.
Definition MatrixCSR.h:299
const rowptr_container_type & row_ptr() const
Definition MatrixCSR.h:493
std::vector< value_type > to_dense() const
Copy to a dense matrix.
Definition MatrixCSR.h:310
MPI_Comm comm() const
Get MPI communicator that matrix is defined on.
Definition MatrixCSR.h:456
const std::array< std::shared_ptr< const common::IndexMap >, 2 > & index_maps() const
Index maps for the row and column space.
Definition MatrixCSR.h:466
MatrixCSR(const Mat &A)
Copy-convert to different container types.
Definition MatrixCSR.h:194
MatrixCSR(const MatrixCSR &A)=default
auto mat_set_values()
Insertion functor for setting values in a matrix. It is typically used in finite element assembly fun...
Definition MatrixCSR.h:95
Definition SparsityPattern.h:26
A vector that can be distributed across processes.
Definition Vector.h:50
container_type & array()
Get the process-local part of the vector.
Definition Vector.h:344
void scatter_fwd_end(U unpack)
End scatter (send) of local data values that are ghosted on other processes.
Definition Vector.h:213
void scatter_fwd_begin(U pack, GetPtr get_ptr)
Begin scatter (send) of local data that is ghosted on other processes.
Definition Vector.h:176
MPI_Datatype mpi_t
Retrieves the MPI data type associated to the provided type.
Definition MPI.h:280
void check_error(MPI_Comm comm, int code)
Check MPI error code. If the error code is not equal to MPI_SUCCESS, then std::abort is called.
Definition MPI.cpp:80
int size(MPI_Comm comm)
Definition MPI.cpp:72
int rank(MPI_Comm comm)
Return process rank for the communicator.
Definition MPI.cpp:64
constexpr std::array< std::int64_t, 2 > local_range(int rank, std::int64_t N, int size)
Return local range for the calling process, partitioning the global [0, N - 1] range across all ranks...
Definition MPI.h:89
Linear algebra interface.
Definition sparsitybuild.h:15
BlockMode
Modes for representing block structured matrices.
Definition MatrixCSR.h:26
@ expanded
Definition MatrixCSR.h:29
auto norm(const V &x, Norm type=Norm::l2)
Compute the norm of the vector.
Definition Vector.h:434