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, Mat A, Mat B, const PetscInt garray[], Mat *mat)

Collective

Input Parameters#

  • comm - MPI communicator

  • A - “diagonal” portion of matrix

  • B - “off-diagonal” portion of matrix, may have empty columns, will be destroyed by this routine

  • garray - global index of B columns

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 becomes part of output mat, B is destroyed by this routine. The user cannot use A and B anymore.

See Also#

Matrices, Mat, MATMPIAIJ, MATSEQAIJ, MatCreateMPIAIJWithSplitArrays()

Level#

advanced

Location#

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


Edit on GitLab

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