MatMPIAIJGetSeqAIJ#

Returns the local pieces of this distributed matrix

Synopsis#

#include "petscmat.h" 
PetscErrorCode MatMPIAIJGetSeqAIJ(Mat A, Mat *Ad, Mat *Ao, const PetscInt *colmap[])

Not Collective

Input Parameter#

Output Parameters#

  • Ad - The local diagonal block as a MATSEQAIJ matrix

  • Ao - The local off-diagonal block as a MATSEQAIJ matrix

  • colmap - An array mapping local column numbers of Ao to global column numbers of the parallel matrix

Note#

The rows in Ad and Ao are in [0, Nr), where Nr is the number of local rows on this process. The columns in Ad are in [0, Nc) where Nc is the number of local columns. The columns are Ao are in [0, Nco), where Nco is the number of nonzero columns in the local off-diagonal piece of the matrix A. The array colmap maps these local column numbers to global column numbers in the original matrix.

Fortran Notes#

MatMPIAIJGetSeqAIJ() Fortran binding is deprecated (since PETSc 3.19), use MatMPIAIJGetSeqAIJF90()

See Also#

Matrices, Mat, MATMPIAIJ, MatMPIAIJGetSeqAIJF90(), MatMPIAIJRestoreSeqAIJF90(), MatMPIAIJGetLocalMat(), MatMPIAIJGetLocalMatCondensed(), MatCreateAIJ(), MATSEQAIJ

Level#

intermediate

Location#

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

Examples#

src/mat/tutorials/ex8.c
src/ksp/ksp/tutorials/ex19.c


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