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 "matrix_csr_impl.h"
11#include <algorithm>
12#include <dolfinx/common/IndexMap.h>
13#include <dolfinx/common/MPI.h>
14#include <dolfinx/graph/AdjacencyList.h>
15#include <mpi.h>
16#include <numeric>
17#include <span>
18#include <utility>
19#include <vector>
20
21namespace dolfinx::la
22{
24enum class BlockMode : int
25{
26 compact = 0,
28 expanded = 1
31};
32
47template <class Scalar, class Container = std::vector<Scalar>,
48 class ColContainer = std::vector<std::int32_t>,
49 class RowPtrContainer = std::vector<std::int64_t>>
51{
52 static_assert(std::is_same_v<typename Container::value_type, Scalar>);
53
54public:
56 using value_type = Scalar;
57
59 using container_type = Container;
60
62 using column_container_type = ColContainer;
63
65 using rowptr_container_type = RowPtrContainer;
66
67 static_assert(std::is_same_v<value_type, typename container_type::value_type>,
68 "Scalar type and container value type must be the same.");
69
93 template <int BS0 = 1, int BS1 = 1>
95 {
96 if ((BS0 != _bs[0] and BS0 > 1 and _bs[0] > 1)
97 or (BS1 != _bs[1] and BS1 > 1 and _bs[1] > 1))
98 {
99 throw std::runtime_error(
100 "Cannot insert blocks of different size than matrix block size");
101 }
102
103 return [&](std::span<const std::int32_t> rows,
104 std::span<const std::int32_t> cols,
105 std::span<const value_type> data) -> int
106 {
107 this->set<BS0, BS1>(data, rows, cols);
108 return 0;
109 };
110 }
111
135 template <int BS0 = 1, int BS1 = 1>
137 {
138 if ((BS0 != _bs[0] and BS0 > 1 and _bs[0] > 1)
139 or (BS1 != _bs[1] and BS1 > 1 and _bs[1] > 1))
140 {
141 throw std::runtime_error(
142 "Cannot insert blocks of different size than matrix block size");
143 }
144
145 return [&](std::span<const std::int32_t> rows,
146 std::span<const std::int32_t> cols,
147 std::span<const value_type> data) -> int
148 {
149 this->add<BS0, BS1>(data, rows, cols);
150 return 0;
151 };
152 }
153
177 MatrixCSR(const SparsityPattern& p, BlockMode mode = BlockMode::compact);
178
181 MatrixCSR(MatrixCSR&& A) = default;
182
186 void set(value_type x) { std::ranges::fill(_data, x); }
187
204 template <int BS0, int BS1>
205 void set(std::span<const value_type> x, std::span<const std::int32_t> rows,
206 std::span<const std::int32_t> cols)
207 {
208 auto set_fn = [](value_type& y, const value_type& x) { y = x; };
209
210 std::int32_t num_rows
211 = _index_maps[0]->size_local() + _index_maps[0]->num_ghosts();
212 assert(x.size() == rows.size() * cols.size() * BS0 * BS1);
213 if (_bs[0] == BS0 and _bs[1] == BS1)
214 {
215 impl::insert_csr<BS0, BS1>(_data, _cols, _row_ptr, x, rows, cols, set_fn,
216 num_rows);
217 }
218 else if (_bs[0] == 1 and _bs[1] == 1)
219 {
220 // Set blocked data in a regular CSR matrix (_bs[0]=1, _bs[1]=1)
221 // with correct sparsity
222 impl::insert_blocked_csr<BS0, BS1>(_data, _cols, _row_ptr, x, rows, cols,
223 set_fn, num_rows);
224 }
225 else
226 {
227 assert(BS0 == 1 and BS1 == 1);
228 // Set non-blocked data in a blocked CSR matrix (BS0=1, BS1=1)
229 impl::insert_nonblocked_csr(_data, _cols, _row_ptr, x, rows, cols, set_fn,
230 num_rows, _bs[0], _bs[1]);
231 }
232 }
233
249 template <int BS0 = 1, int BS1 = 1>
250 void add(std::span<const value_type> x, std::span<const std::int32_t> rows,
251 std::span<const std::int32_t> cols)
252 {
253 auto add_fn = [](value_type& y, const value_type& x) { y += x; };
254
255 assert(x.size() == rows.size() * cols.size() * BS0 * BS1);
256 if (_bs[0] == BS0 and _bs[1] == BS1)
257 {
258 impl::insert_csr<BS0, BS1>(_data, _cols, _row_ptr, x, rows, cols, add_fn,
259 _row_ptr.size());
260 }
261 else if (_bs[0] == 1 and _bs[1] == 1)
262 {
263 // Add blocked data to a regular CSR matrix (_bs[0]=1, _bs[1]=1)
264 impl::insert_blocked_csr<BS0, BS1>(_data, _cols, _row_ptr, x, rows, cols,
265 add_fn, _row_ptr.size());
266 }
267 else
268 {
269 assert(BS0 == 1 and BS1 == 1);
270 // Add non-blocked data to a blocked CSR matrix (BS0=1, BS1=1)
271 impl::insert_nonblocked_csr(_data, _cols, _row_ptr, x, rows, cols, add_fn,
272 _row_ptr.size(), _bs[0], _bs[1]);
273 }
274 }
275
277 std::int32_t num_owned_rows() const { return _index_maps[0]->size_local(); }
278
280 std::int32_t num_all_rows() const { return _row_ptr.size() - 1; }
281
291 std::vector<value_type> to_dense() const;
292
300 {
303 }
304
312 void scatter_rev_begin();
313
319 void scatter_rev_end();
320
323 double squared_norm() const;
324
332 std::shared_ptr<const common::IndexMap> index_map(int dim) const
333 {
334 return _index_maps.at(dim);
335 }
336
339 container_type& values() { return _data; }
340
343 const container_type& values() const { return _data; }
344
347 const rowptr_container_type& row_ptr() const { return _row_ptr; }
348
351 const column_container_type& cols() const { return _cols; }
352
361 {
362 return _off_diagonal_offset;
363 }
364
367 std::array<int, 2> block_size() const { return _bs; }
368
369private:
370 // Maps for the distribution of the ows and columns
371 std::array<std::shared_ptr<const common::IndexMap>, 2> _index_maps;
372
373 // Block mode (compact or expanded)
374 BlockMode _block_mode;
375
376 // Block sizes
377 std::array<int, 2> _bs;
378
379 // Matrix data
380 container_type _data;
382 rowptr_container_type _row_ptr;
383
384 // Start of off-diagonal (unowned columns) on each row
385 rowptr_container_type _off_diagonal_offset;
386
387 // Neighborhood communicator (ghost->owner communicator for rows)
388 dolfinx::MPI::Comm _comm;
389
390 // -- Precomputed data for scatter_rev/update
391
392 // Request in non-blocking communication
393 MPI_Request _request;
394
395 // Position in _data to add received data
396 std::vector<int> _unpack_pos;
397
398 // Displacements for alltoall for each neighbor when sending and
399 // receiving
400 std::vector<int> _val_send_disp, _val_recv_disp;
401
402 // Ownership of each row, by neighbor (for the neighbourhood defined
403 // on _comm)
404 std::vector<int> _ghost_row_to_rank;
405
406 // Temporary stores for data during non-blocking communication
407 container_type _ghost_value_data;
408 container_type _ghost_value_data_in;
409};
410//-----------------------------------------------------------------------------
411template <class U, class V, class W, class X>
413 : _index_maps({p.index_map(0),
414 std::make_shared<common::IndexMap>(p.column_index_map())}),
415 _block_mode(mode), _bs({p.block_size(0), p.block_size(1)}),
416 _data(p.num_nonzeros() * _bs[0] * _bs[1], 0),
417 _cols(p.graph().first.begin(), p.graph().first.end()),
418 _row_ptr(p.graph().second.begin(), p.graph().second.end()),
419 _comm(MPI_COMM_NULL)
420{
421 if (_block_mode == BlockMode::expanded)
422 {
423 // Rebuild IndexMaps
424 for (int i = 0; i < 2; ++i)
425 {
426 const auto im = _index_maps[i];
427 const int size_local = im->size_local() * _bs[i];
428 std::span ghost_i = im->ghosts();
429 std::vector<std::int64_t> ghosts;
430 const std::vector<int> ghost_owner_i(im->owners().begin(),
431 im->owners().end());
432 std::vector<int> src_rank;
433 for (std::size_t j = 0; j < ghost_i.size(); ++j)
434 {
435 for (int k = 0; k < _bs[i]; ++k)
436 {
437 ghosts.push_back(ghost_i[j] * _bs[i] + k);
438 src_rank.push_back(ghost_owner_i[j]);
439 }
440 }
441 const std::array<std::vector<int>, 2> src_dest0
442 = {std::vector(_index_maps[i]->src().begin(),
443 _index_maps[i]->src().end()),
444 std::vector(_index_maps[i]->dest().begin(),
445 _index_maps[i]->dest().end())};
446 _index_maps[i] = std::make_shared<common::IndexMap>(
447 _index_maps[i]->comm(), size_local, src_dest0, ghosts, src_rank);
448 }
449
450 // Convert sparsity pattern and set _bs to 1
451
452 column_container_type new_cols;
453 new_cols.reserve(_data.size());
454 rowptr_container_type new_row_ptr = {0};
455 new_row_ptr.reserve(_row_ptr.size() * _bs[0]);
456 std::span<const std::int32_t> num_diag_nnz = p.off_diagonal_offsets();
457 for (std::size_t i = 0; i < _row_ptr.size() - 1; ++i)
458 {
459 // Repeat row _bs[0] times
460 for (int q0 = 0; q0 < _bs[0]; ++q0)
461 {
462 _off_diagonal_offset.push_back(new_row_ptr.back()
463 + num_diag_nnz[i] * _bs[1]);
464 for (auto j = _row_ptr[i]; j < _row_ptr[i + 1]; ++j)
465 {
466 for (int q1 = 0; q1 < _bs[1]; ++q1)
467 new_cols.push_back(_cols[j] * _bs[1] + q1);
468 }
469 new_row_ptr.push_back(new_cols.size());
470 }
471 }
472 _cols = new_cols;
473 _row_ptr = new_row_ptr;
474 _bs[0] = 1;
475 _bs[1] = 1;
476 }
477 else
478 {
479 // Compute off-diagonal offset for each row (compact)
480 std::span<const std::int32_t> num_diag_nnz = p.off_diagonal_offsets();
481 _off_diagonal_offset.reserve(num_diag_nnz.size());
482 std::ranges::transform(num_diag_nnz, _row_ptr,
483 std::back_inserter(_off_diagonal_offset),
484 std::plus{});
485 }
486
487 // Some short-hand
488 const std::array local_size
489 = {_index_maps[0]->size_local(), _index_maps[1]->size_local()};
490 const std::array local_range
491 = {_index_maps[0]->local_range(), _index_maps[1]->local_range()};
492 std::span ghosts1 = _index_maps[1]->ghosts();
493
494 std::span ghosts0 = _index_maps[0]->ghosts();
495 std::span src_ranks = _index_maps[0]->src();
496 std::span dest_ranks = _index_maps[0]->dest();
497
498 // Create neighbourhood communicator (owner <- ghost)
499 MPI_Comm comm;
500 MPI_Dist_graph_create_adjacent(_index_maps[0]->comm(), dest_ranks.size(),
501 dest_ranks.data(), MPI_UNWEIGHTED,
502 src_ranks.size(), src_ranks.data(),
503 MPI_UNWEIGHTED, MPI_INFO_NULL, false, &comm);
504 _comm = dolfinx::MPI::Comm(comm, false);
505
506 // Build map from ghost row index position to owning (neighborhood)
507 // rank
508 _ghost_row_to_rank.reserve(_index_maps[0]->owners().size());
509 for (int r : _index_maps[0]->owners())
510 {
511 auto it = std::ranges::lower_bound(src_ranks, r);
512 assert(it != src_ranks.end() and *it == r);
513 int pos = std::distance(src_ranks.begin(), it);
514 _ghost_row_to_rank.push_back(pos);
515 }
516
517 // Compute size of data to send to each neighbor
518 std::vector<std::int32_t> data_per_proc(src_ranks.size(), 0);
519 for (std::size_t i = 0; i < _ghost_row_to_rank.size(); ++i)
520 {
521 assert(_ghost_row_to_rank[i] < data_per_proc.size());
522 std::size_t pos = local_size[0] + i;
523 data_per_proc[_ghost_row_to_rank[i]] += _row_ptr[pos + 1] - _row_ptr[pos];
524 }
525
526 // Compute send displacements
527 _val_send_disp.resize(src_ranks.size() + 1, 0);
528 std::partial_sum(data_per_proc.begin(), data_per_proc.end(),
529 std::next(_val_send_disp.begin()));
530
531 // For each ghost row, pack and send indices to neighborhood
532 std::vector<std::int64_t> ghost_index_data(2 * _val_send_disp.back());
533 {
534 std::vector<int> insert_pos = _val_send_disp;
535 for (std::size_t i = 0; i < _ghost_row_to_rank.size(); ++i)
536 {
537 const int rank = _ghost_row_to_rank[i];
538 int row_id = local_size[0] + i;
539 for (int j = _row_ptr[row_id]; j < _row_ptr[row_id + 1]; ++j)
540 {
541 // Get position in send buffer
542 const std::int32_t idx_pos = 2 * insert_pos[rank];
543
544 // Pack send data (row, col) as global indices
545 ghost_index_data[idx_pos] = ghosts0[i];
546 if (std::int32_t col_local = _cols[j]; col_local < local_size[1])
547 ghost_index_data[idx_pos + 1] = col_local + local_range[1][0];
548 else
549 ghost_index_data[idx_pos + 1] = ghosts1[col_local - local_size[1]];
550
551 insert_pos[rank] += 1;
552 }
553 }
554 }
555
556 // Communicate data with neighborhood
557 std::vector<std::int64_t> ghost_index_array;
558 std::vector<int> recv_disp;
559 {
560 std::vector<int> send_sizes;
561 std::ranges::transform(data_per_proc, std::back_inserter(send_sizes),
562 [](auto x) { return 2 * x; });
563
564 std::vector<int> recv_sizes(dest_ranks.size());
565 send_sizes.reserve(1);
566 recv_sizes.reserve(1);
567 MPI_Neighbor_alltoall(send_sizes.data(), 1, MPI_INT, recv_sizes.data(), 1,
568 MPI_INT, _comm.comm());
569
570 // Build send/recv displacement
571 std::vector<int> send_disp = {0};
572 std::partial_sum(send_sizes.begin(), send_sizes.end(),
573 std::back_inserter(send_disp));
574 recv_disp = {0};
575 std::partial_sum(recv_sizes.begin(), recv_sizes.end(),
576 std::back_inserter(recv_disp));
577
578 ghost_index_array.resize(recv_disp.back());
579 MPI_Neighbor_alltoallv(ghost_index_data.data(), send_sizes.data(),
580 send_disp.data(), MPI_INT64_T,
581 ghost_index_array.data(), recv_sizes.data(),
582 recv_disp.data(), MPI_INT64_T, _comm.comm());
583 }
584
585 // Store receive displacements for future use, when transferring
586 // data values
587 _val_recv_disp.resize(recv_disp.size());
588 const int bs2 = _bs[0] * _bs[1];
589 std::ranges::transform(recv_disp, _val_recv_disp.begin(),
590 [&bs2](auto d) { return bs2 * d / 2; });
591 std::ranges::transform(_val_send_disp, _val_send_disp.begin(),
592 [&bs2](auto d) { return d * bs2; });
593
594 // Global-to-local map for ghost columns
595 std::vector<std::pair<std::int64_t, std::int32_t>> global_to_local;
596 global_to_local.reserve(ghosts1.size());
597 for (std::int64_t idx : ghosts1)
598 global_to_local.push_back({idx, global_to_local.size() + local_size[1]});
599 std::ranges::sort(global_to_local);
600
601 // Compute location in which data for each index should be stored
602 // when received
603 for (std::size_t i = 0; i < ghost_index_array.size(); i += 2)
604 {
605 // Row must be on this process
606 const std::int32_t local_row = ghost_index_array[i] - local_range[0][0];
607 assert(local_row >= 0 and local_row < local_size[0]);
608
609 // Column may be owned or unowned
610 std::int32_t local_col = ghost_index_array[i + 1] - local_range[1][0];
611 if (local_col < 0 or local_col >= local_size[1])
612 {
613 auto it = std::ranges::lower_bound(
614 global_to_local, std::pair(ghost_index_array[i + 1], -1),
615 [](auto& a, auto& b) { return a.first < b.first; });
616 assert(it != global_to_local.end()
617 and it->first == ghost_index_array[i + 1]);
618 local_col = it->second;
619 }
620 auto cit0 = std::next(_cols.begin(), _row_ptr[local_row]);
621 auto cit1 = std::next(_cols.begin(), _row_ptr[local_row + 1]);
622
623 // Find position of column index and insert data
624 auto cit = std::lower_bound(cit0, cit1, local_col);
625 assert(cit != cit1);
626 assert(*cit == local_col);
627 std::size_t d = std::distance(_cols.begin(), cit);
628 _unpack_pos.push_back(d);
629 }
630}
631//-----------------------------------------------------------------------------
632template <typename U, typename V, typename W, typename X>
633std::vector<typename MatrixCSR<U, V, W, X>::value_type>
635{
636 const std::size_t nrows = num_all_rows();
637 const std::size_t ncols = _index_maps[1]->size_global();
638 std::vector<value_type> A(nrows * ncols * _bs[0] * _bs[1], 0.0);
639 for (std::size_t r = 0; r < nrows; ++r)
640 for (std::int32_t j = _row_ptr[r]; j < _row_ptr[r + 1]; ++j)
641 for (int i0 = 0; i0 < _bs[0]; ++i0)
642 for (int i1 = 0; i1 < _bs[1]; ++i1)
643 {
644 std::array<std::int32_t, 1> local_col{_cols[j]};
645 std::array<std::int64_t, 1> global_col{0};
646 _index_maps[1]->local_to_global(local_col, global_col);
647 A[(r * _bs[1] + i0) * ncols * _bs[0] + global_col[0] * _bs[1] + i1]
648 = _data[j * _bs[0] * _bs[1] + i0 * _bs[1] + i1];
649 }
650
651 return A;
652}
653//-----------------------------------------------------------------------------
654template <typename U, typename V, typename W, typename X>
656{
657 const std::int32_t local_size0 = _index_maps[0]->size_local();
658 const std::int32_t num_ghosts0 = _index_maps[0]->num_ghosts();
659 const int bs2 = _bs[0] * _bs[1];
660
661 // For each ghost row, pack and send values to send to neighborhood
662 std::vector<int> insert_pos = _val_send_disp;
663 _ghost_value_data.resize(_val_send_disp.back());
664 for (int i = 0; i < num_ghosts0; ++i)
665 {
666 const int rank = _ghost_row_to_rank[i];
667
668 // Get position in send buffer to place data to send to this
669 // neighbour
670 const std::int32_t val_pos = insert_pos[rank];
671 std::copy(std::next(_data.data(), _row_ptr[local_size0 + i] * bs2),
672 std::next(_data.data(), _row_ptr[local_size0 + i + 1] * bs2),
673 std::next(_ghost_value_data.begin(), val_pos));
674 insert_pos[rank]
675 += bs2 * (_row_ptr[local_size0 + i + 1] - _row_ptr[local_size0 + i]);
676 }
677
678 _ghost_value_data_in.resize(_val_recv_disp.back());
679
680 // Compute data sizes for send and receive from displacements
681 std::vector<int> val_send_count(_val_send_disp.size() - 1);
682 std::adjacent_difference(std::next(_val_send_disp.begin()),
683 _val_send_disp.end(), val_send_count.begin());
684
685 std::vector<int> val_recv_count(_val_recv_disp.size() - 1);
686 std::adjacent_difference(std::next(_val_recv_disp.begin()),
687 _val_recv_disp.end(), val_recv_count.begin());
688
689 int status = MPI_Ineighbor_alltoallv(
690 _ghost_value_data.data(), val_send_count.data(), _val_send_disp.data(),
691 dolfinx::MPI::mpi_t<value_type>, _ghost_value_data_in.data(),
692 val_recv_count.data(), _val_recv_disp.data(),
693 dolfinx::MPI::mpi_t<value_type>, _comm.comm(), &_request);
694 assert(status == MPI_SUCCESS);
695}
696//-----------------------------------------------------------------------------
697template <typename U, typename V, typename W, typename X>
699{
700 int status = MPI_Wait(&_request, MPI_STATUS_IGNORE);
701 assert(status == MPI_SUCCESS);
702
703 _ghost_value_data.clear();
704 _ghost_value_data.shrink_to_fit();
705
706 // Add to local rows
707 const int bs2 = _bs[0] * _bs[1];
708 assert(_ghost_value_data_in.size() == _unpack_pos.size() * bs2);
709 for (std::size_t i = 0; i < _unpack_pos.size(); ++i)
710 for (int j = 0; j < bs2; ++j)
711 _data[_unpack_pos[i] * bs2 + j] += _ghost_value_data_in[i * bs2 + j];
712
713 _ghost_value_data_in.clear();
714 _ghost_value_data_in.shrink_to_fit();
715
716 // Set ghost row data to zero
717 const std::int32_t local_size0 = _index_maps[0]->size_local();
718 std::fill(std::next(_data.begin(), _row_ptr[local_size0] * bs2), _data.end(),
719 0);
720}
721//-----------------------------------------------------------------------------
722template <typename U, typename V, typename W, typename X>
724{
725 const std::size_t num_owned_rows = _index_maps[0]->size_local();
726 const int bs2 = _bs[0] * _bs[1];
727 assert(num_owned_rows < _row_ptr.size());
728 double norm_sq_local = std::accumulate(
729 _data.cbegin(), std::next(_data.cbegin(), _row_ptr[num_owned_rows] * bs2),
730 double(0), [](auto norm, value_type y) { return norm + std::norm(y); });
731 double norm_sq;
732 MPI_Allreduce(&norm_sq_local, &norm_sq, 1, MPI_DOUBLE, MPI_SUM, _comm.comm());
733 return norm_sq;
734}
735//-----------------------------------------------------------------------------
736
737} // namespace dolfinx::la
A duplicate MPI communicator and manage lifetime of the communicator.
Definition MPI.h:44
Distributed sparse matrix.
Definition MatrixCSR.h:51
const container_type & values() const
Definition MatrixCSR.h:343
std::shared_ptr< const common::IndexMap > index_map(int dim) const
Index maps for the row and column space.
Definition MatrixCSR.h:332
const rowptr_container_type & off_diag_offset() const
Definition MatrixCSR.h:360
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:205
RowPtrContainer rowptr_container_type
Row pointer container type.
Definition MatrixCSR.h:65
void scatter_rev_end()
End transfer of ghost row data to owning ranks.
Definition MatrixCSR.h:698
container_type & values()
Definition MatrixCSR.h:339
auto mat_add_values()
Insertion functor for adding values to a matrix. It is typically used in finite element assembly func...
Definition MatrixCSR.h:136
std::vector< value_type > to_dense() const
Copy to a dense matrix.
Definition MatrixCSR.h:634
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:250
std::int32_t num_owned_rows() const
Number of local rows excluding ghost rows.
Definition MatrixCSR.h:277
MatrixCSR(const SparsityPattern &p, BlockMode mode=BlockMode::compact)
Create a distributed matrix.
Definition MatrixCSR.h:412
ColContainer column_container_type
Column index container type.
Definition MatrixCSR.h:62
MatrixCSR(MatrixCSR &&A)=default
double squared_norm() const
Compute the Frobenius norm squared across all processes.
Definition MatrixCSR.h:723
void scatter_rev()
Transfer ghost row data to the owning ranks accumulating received values on the owned rows,...
Definition MatrixCSR.h:299
Container container_type
Matrix entries container type.
Definition MatrixCSR.h:59
Scalar value_type
Scalar type.
Definition MatrixCSR.h:56
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:655
const column_container_type & cols() const
Definition MatrixCSR.h:351
void set(value_type x)
Set all non-zero local entries to a value including entries in ghost rows.
Definition MatrixCSR.h:186
std::array< int, 2 > block_size() const
Definition MatrixCSR.h:367
std::int32_t num_all_rows() const
Number of local rows including ghost rows.
Definition MatrixCSR.h:280
const rowptr_container_type & row_ptr() const
Definition MatrixCSR.h:347
auto mat_set_values()
Insertion functor for setting values in a matrix. It is typically used in finite element assembly fun...
Definition MatrixCSR.h:94
Definition SparsityPattern.h:26
MPI_Datatype mpi_t
Retrieves the MPI data type associated to the provided type.
Definition MPI.h:282
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:91
Linear algebra interface.
Definition sparsitybuild.h:15
BlockMode
Modes for representing block structured matrices.
Definition MatrixCSR.h:25
auto norm(const V &x, Norm type=Norm::l2)
Definition Vector.h:268