petsc-3.11.4 2019-09-28
Report Typos and Errors

MatCreateMPIBAIJWithArrays

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

Synopsis

#include "petscmat.h"  
PetscErrorCode  MatCreateMPIBAIJWithArrays(MPI_Comm comm,PetscInt bs,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
bs - the block size, only a block size of 1 is supported
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; that is i[0] = 0, i[row] = i[row-1] + number of block elements in that rowth block row of the matrix
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 order of the entries in values is the same as the block compressed sparse row storage format; that is, it is the same as a three dimensional array in Fortran values(bs,bs,nnz) that contains the first column of the first block, followed by the second column of the first block etc etc. That is, the blocks are contiguous in memory with column-major ordering within blocks.

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

Keywords

matrix, aij, compressed row, sparse, parallel

See Also

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

Level

intermediate

Location

src/mat/impls/baij/mpi/mpibaij.c
Index of all Mat routines
Table of Contents for all manual pages
Index of all manual pages