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 notNULL
, B should be the offdiag matrix using local col ids and of size garraygarray - either
NULL
or the global index ofB
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
and B
becomes part of output mat. 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
Index of all Mat routines
Table of Contents for all manual pages
Index of all manual pages