MatCreateMPIAIJWithSeqAIJ#

creates a MATMPIAIJ matrix using MATSEQAIJ matrices that contain the “diagonal” and “off-diagonal” part of the matrix in CSR format.

Synopsis#

#include "petscmat.h" 
PetscErrorCode MatCreateMPIAIJWithSeqAIJ(MPI_Comm comm, PetscInt M, PetscInt N, Mat A, Mat B, PetscInt *garray, Mat *mat)

Collective

Input Parameters#

  • comm - MPI communicator

  • M - the global row size

  • N - the global column size

  • A - “diagonal” portion of matrix

  • B - if garray is NULL, B should be the offdiag matrix using global col ids and of size N - if garray is not NULL, B should be the offdiag matrix using local col ids and of size garray

  • garray - either NULL or the global index of B columns. If not NULL, it should be allocated by PetscMalloc1() and will be owned by mat thereafter.

Output Parameter#

  • mat - the matrix, with input A as its local diagonal matrix

Notes#

See MatCreateAIJ() for the definition of “diagonal” and “off-diagonal” portion of the matrix.

A and B becomes part of output mat. The user cannot use A and B anymore.

If garray is NULL, B will be compacted to use local indices. In this sense, B’s sparsity pattern (nonzerostate) will be changed. If B is a device matrix, we need to somehow also update B’s copy on device. We do so by increasing B’s nonzerostate. In use of B on device, device matrix types should detect this change (ref. internal routines MatSeqAIJCUSPARSECopyToGPU() or MatAssemblyEnd_SeqAIJKokkos()) and will just destroy and then recreate the device copy of B. It is not optimal, but is easy to implement and less hacky. To avoid this overhead, try to compute garray yourself, see algorithms in the private function MatSetUpMultiply_MPIAIJ().

The NULL-ness of garray doesn’t need to be collective, in other words, garray can be NULL on some processes while not on others.

See Also#

Matrices, Mat, MATMPIAIJ, MATSEQAIJ, MatCreateMPIAIJWithSplitArrays()

Level#

advanced

Location#

src/mat/impls/aij/mpi/mpiaij.c


Index of all Mat routines
Table of Contents for all manual pages
Index of all manual pages