MatCreateMPIAIJWithSplitArrays

creates a MPI AIJ matrix using arrays that contain the "diagonal" and "off-diagonal" part of the matrix in CSR format.

Synopsis

#include "petscmat.h" 
PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithSplitArrays(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt i[],PetscInt j[],PetscScalar a[],
                                                         PetscInt oi[], PetscInt oj[],PetscScalar oa[],Mat *mat)
Collective on MPI_Comm

Input Parameters

comm - MPI communicator
m - number of local rows (Cannot be PETSC_DECIDE)
n - This value should be the same as the local size used in creating the x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have calculated if N is given) For square matrices n is almost always m.
M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
i - row indices for "diagonal" portion of matrix
j - column indices
a - matrix values
oi - row indices for "off-diagonal" portion of matrix
oj - column indices
oa - matrix values

Output Parameter

mat -the matrix

Notes

The i, j, and a arrays ARE NOT copied by this routine into the internal format used by PETSc.

The i and j indices are 0 based

See MatCreateMPIAIJ() for the definition of "diagonal" and "off-diagonal" portion of the matrix

Keywords

matrix, aij, compressed row, sparse, parallel

See Also

MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatMPIAIJSetPreallocation(), MatMPIAIJSetPreallocationCSR(),
MPIAIJ, MatCreateMPIAIJ(), MatCreateMPIAIJWithArrays()

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