petsc-3.7.5 2017-01-01
Report Typos and Errors

MatCreateMPIAIJWithArrays

creates a MPI AIJ matrix using arrays that contain in standard CSR format the local rows.

Synopsis

#include "petscmat.h" 
PetscErrorCode  MatCreateMPIAIJWithArrays(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt i[],const PetscInt j[],const PetscScalar a[],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
j - column indices
a - matrix values

Output Parameter

mat -the matrix

Notes

The i, j, and a arrays ARE copied by this routine into the internal format used by PETSc; thus you CANNOT change the matrix entries by changing the values of a[] after you have called this routine. Use MatCreateMPIAIJWithSplitArrays() to avoid needing to copy the arrays.

The i and j indices are 0 based, and i indices are indices corresponding to the local j array.

The format which is used for the sparse matrix input, is equivalent to a row-major ordering.. i.e for the following matrix, the input data expected is as shown

       1 0 0
       2 0 3     P0
      -------
       4 5 6     P1

    Process0 [P0]: rows_owned=[0,1]
       i =  {0,1,3}  [size = nrow+1  = 2+1]
       j =  {0,0,2}  [size = 3]
       v =  {1,2,3}  [size = 3]

    Process1 [P1]: rows_owned=[2]
       i =  {0,3}    [size = nrow+1  = 1+1]
       j =  {0,1,2}  [size = 3]
       v =  {4,5,6}  [size = 3]

Keywords

matrix, aij, compressed row, sparse, parallel

See Also

MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatMPIAIJSetPreallocation(), MatMPIAIJSetPreallocationCSR(),
MPIAIJ, MatCreateAIJ(), MatCreateMPIAIJWithSplitArrays()

Level:intermediate
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