MatRestoreRow#
Frees any temporary space allocated by MatGetRow()
.
Synopsis#
#include "petscmat.h"
PetscErrorCode MatRestoreRow(Mat mat, PetscInt row, PetscInt *ncols, const PetscInt *cols[], const PetscScalar *vals[])
Not Collective
Input Parameters#
mat - the matrix
row - the row to get
ncols - the number of nonzeros
cols - the columns of the nonzeros
vals - if nonzero the column values
Notes#
This routine should be called after you have finished examining the entries.
This routine zeros out ncols
, cols
, and vals
. This is to prevent accidental
us of the array after it has been restored. If you pass NULL
, it will
not zero the pointers. Use of cols
or vals
after MatRestoreRow()
is invalid.
Fortran Note#
MatRestoreRow()
MUST be called after MatGetRow()
before another call to MatGetRow()
can be made.
See Also#
Level#
advanced
Location#
Examples#
src/mat/tutorials/ex16.c
src/mat/tutorials/ex12.c
src/ksp/ksp/tutorials/ex59.c
Implementations#
MatRestoreRow_MPIAdj() in src/mat/impls/adj/mpi/mpiadj.c
MatRestoreRow_MPIAIJ() in src/mat/impls/aij/mpi/mpiaij.c
MatRestoreRow_SeqAIJ() in src/mat/impls/aij/seq/aij.c
MatRestoreRow_MPIBAIJ() in src/mat/impls/baij/mpi/mpibaij.c
MatRestoreRow_SeqBAIJ() in src/mat/impls/baij/seq/baij.c
MatRestoreRow_ConstantDiagonal() in src/mat/impls/cdiagonal/cdiagonal.c
MatRestoreRow_MPIDense() in src/mat/impls/dense/mpi/mpidense.c
MatRestoreRow_SeqDense() in src/mat/impls/dense/seq/dense.c
MatRestoreRow_Diagonal() in src/mat/impls/diagonal/diagonal.c
MatRestoreRow_Htool() in src/mat/impls/htool/htool.cxx
MatRestoreRow_HYPRE() in src/mat/impls/hypre/mhypre.c
MatRestoreRow_SeqKAIJ() in src/mat/impls/kaij/kaij.c
MatRestoreRow_MPIKAIJ() in src/mat/impls/kaij/kaij.c
MatRestoreRow_MPISBAIJ() in src/mat/impls/sbaij/mpi/mpisbaij.c
MatRestoreRow_SeqSBAIJ() in src/mat/impls/sbaij/seq/sbaij.c
MatRestoreRow_SeqSELL() in src/mat/impls/sell/seq/sell.c
Index of all Mat routines
Table of Contents for all manual pages
Index of all manual pages