PETSc version 3.17.5
MatMPIAIJGetLocalMatMerge
Creates a SeqAIJ from a MATMPIAIJ matrix by taking all its local rows and putting them into a sequential matrix with mlocal rows and n columns. Where n is the sum of the number of columns of the diagonal and offdiagonal part
Synopsis
#include "petscmat.h"
PetscErrorCode MatMPIAIJGetLocalMatMerge(Mat A,MatReuse scall,IS *glob,Mat *A_loc)
Not Collective
Input Parameters
Output Parameters
| glob | - sequential IS with global indices associated with the columns of the local sequential matrix generated (can be NULL)
|
| A_loc | - the local sequential matrix generated
|
Notes
This is different from MatMPIAIJGetLocalMat() since the first columns in the returning matrix are those associated with the diagonal part, then those associated with the offdiagonal part (in its local ordering)
See Also
MatGetOwnershipRange(), MatMPIAIJGetLocalMat(), MatMPIAIJGetLocalMatCondensed()
Level
developer
Location
src/mat/impls/aij/mpi/mpiaij.c
Implementations
MatMPIAIJGetLocalMatMerge_MPIAIJKokkos in src/mat/impls/aij/mpi/kokkos/mpiaijkok.kokkos.cxx
MatMPIAIJGetLocalMatMerge_MPIAIJCUSPARSE in src/mat/impls/aij/mpi/mpicusparse/mpiaijcusparse.cu
Index of all Mat routines
Table of Contents for all manual pages
Index of all manual pages