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