:orphan: # MatInvertVariableBlockDiagonal Inverts the point block diagonal entries. ## Synopsis ``` #include "petscmat.h" PetscErrorCode MatInvertVariableBlockDiagonal(Mat mat, PetscInt nblocks, const PetscInt *bsizes, PetscScalar *values) ``` Collective; No Fortran Support ## Input Parameters - ***mat -*** the matrix - ***nblocks -*** the number of blocks on the process, set with `MatSetVariableBlockSizes()` - ***bsizes -*** the size of each block on the process, set with `MatSetVariableBlockSizes()` ## Output Parameter - ***values -*** the block inverses in column major order (FORTRAN-like) ## Notes Use `MatInvertBlockDiagonal()` if all blocks have the same size The blocks never overlap between two MPI ranks, use `MatInvertVariableBlockEnvelope()` for that case ## See Also [](ch_matrices), `Mat`, `MatInvertBlockDiagonal()`, `MatSetVariableBlockSizes()`, `MatInvertVariableBlockEnvelope()` ## Level advanced ## Location src/mat/interface/matrix.c ## Implementations MatInvertVariableBlockDiagonal_MPIAIJ in src/mat/impls/aij/mpi/mpiaij.c
MatInvertVariableBlockDiagonal_SeqAIJ in src/mat/impls/aij/seq/aij.c
--- [Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/mat/interface/matrix.c) [Index of all Mat routines](index.md) [Table of Contents for all manual pages](/manualpages/index.md) [Index of all manual pages](/manualpages/singleindex.md)