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

MatView

Visualizes a matrix object.

Synopsis

#include "petscmat.h" 
PetscErrorCode MatView(Mat mat,PetscViewer viewer)
Collective on Mat

Input Parameters

mat - the matrix
viewer - visualization context

Notes

The available visualization contexts include
PETSC_VIEWER_STDOUT_SELF - for sequential matrices
PETSC_VIEWER_STDOUT_WORLD - for parallel matrices created on PETSC_COMM_WORLD
PETSC_VIEWER_STDOUT_(comm) - for matrices created on MPI communicator comm
PETSC_VIEWER_DRAW_WORLD - graphical display of nonzero structure

The user can open alternative visualization contexts with

PetscViewerASCIIOpen() - Outputs matrix to a specified file
PetscViewerBinaryOpen() - Outputs matrix in binary to a specified file; corresponding input uses MatLoad()
PetscViewerDrawOpen() - Outputs nonzero matrix structure to an X window display
PetscViewerSocketOpen() - Outputs matrix to Socket viewer. Currently only the sequential dense and AIJ matrix types support the Socket viewer.

The user can call PetscViewerPushFormat() to specify the output format of ASCII printed objects (when using PETSC_VIEWER_STDOUT_SELF, PETSC_VIEWER_STDOUT_WORLD and PetscViewerASCIIOpen). Available formats include

PETSC_VIEWER_DEFAULT - default, prints matrix contents
PETSC_VIEWER_ASCII_MATLAB - prints matrix contents in Matlab format
PETSC_VIEWER_ASCII_DENSE - prints entire matrix including zeros
PETSC_VIEWER_ASCII_COMMON - prints matrix contents, using a sparse format common among all matrix types
PETSC_VIEWER_ASCII_IMPL - prints matrix contents, using an implementation-specific format (which is in many cases the same as the default)
PETSC_VIEWER_ASCII_INFO - prints basic information about the matrix size and structure (not the matrix entries)
PETSC_VIEWER_ASCII_INFO_DETAIL - prints more detailed information about the matrix structure

Options Database Keys

-mat_view ::ascii_info - Prints info on matrix at conclusion of MatAssemblyEnd()
-mat_view ::ascii_info_detail - Prints more detailed info
-mat_view - Prints matrix in ASCII format
-mat_view ::ascii_matlab - Prints matrix in Matlab format
-mat_view draw - PetscDraws nonzero structure of matrix, using MatView() and PetscDrawOpenX().
-display <name> - Sets display name (default is host)
-draw_pause <sec> - Sets number of seconds to pause after display
-mat_view socket - Sends matrix to socket, can be accessed from Matlab (see Users-Manual: Chapter 12 Using MATLAB with PETSc for details)
-viewer_socket_machine <machine> -
-viewer_socket_port <port> -
-mat_view binary - save matrix to file in binary format
-viewer_binary_filename <name> -

Notes

The ASCII viewers are only recommended for small matrices on at most a moderate number of processes, the program will seemingly hang and take hours for larger matrices, for larger matrices one should use the binary format.

See the manual page for MatLoad() for the exact format of the binary file when the binary viewer is used.

See share/petsc/matlab/PetscBinaryRead.m for a Matlab code that can read in the binary file when the binary viewer is used.

One can use '-mat_view draw -draw_pause -1' to pause the graphical display of matrix nonzero structure, and then use the following mouse functions.

left mouse: zoom in- . middle mouse: zoom out
right mouse: continue with the simulation-

See Also

PetscViewerPushFormat(), PetscViewerASCIIOpen(), PetscViewerDrawOpen(),
PetscViewerSocketOpen(), PetscViewerBinaryOpen(), MatLoad()

Level

beginner

Location

src/mat/interface/matrix.c

Examples

src/vec/vec/examples/tutorials/ex6.c.html
src/mat/examples/tutorials/ex4.c.html
src/mat/examples/tutorials/ex8.c.html
src/mat/examples/tutorials/ex11.c.html
src/mat/examples/tutorials/ex12.c.html
src/mat/examples/tutorials/ex16.c.html
src/mat/examples/tutorials/ex17.c.html
src/mat/examples/tutorials/ex4f.F90.html
src/mat/examples/tutorials/ex17f.F90.html
src/ksp/pc/examples/tutorials/ex2.c.html
src/ksp/ksp/examples/tutorials/ex28.c.html

Implementations

MatView_MPI_DA in src/dm/impls/da/fdda.c
MatView_MPIAdj_ASCII in src/mat/impls/adj/mpi/mpiadj.c
MatView_MPIAdj in src/mat/impls/adj/mpi/mpiadj.c
MatView_SparseElemental in src/mat/impls/aij/mpi/clique/clique.cxx
MatView_MKL_CPARDISO in src/mat/impls/aij/mpi/mkl_cpardiso/mkl_cpardiso.c
MatView_MPIAIJ_Binary in src/mat/impls/aij/mpi/mpiaij.c
MatView_MPIAIJ_ASCIIorDraworSocket in src/mat/impls/aij/mpi/mpiaij.c
MatView_MPIAIJ in src/mat/impls/aij/mpi/mpiaij.c
MatView_MPIAIJ_PtAP in src/mat/impls/aij/mpi/mpiptap.c
MatView_MUMPS in src/mat/impls/aij/mpi/mumps/mumps.c
MatView_PaStiX in src/mat/impls/aij/mpi/pastix/pastix.c
MatView_Info_STRUMPACK in src/mat/impls/aij/mpi/strumpack/strumpack.c
MatView_STRUMPACK in src/mat/impls/aij/mpi/strumpack/strumpack.c
MatView_Info_SuperLU_DIST in src/mat/impls/aij/mpi/superlu_dist/superlu_dist.c
MatView_SuperLU_DIST in src/mat/impls/aij/mpi/superlu_dist/superlu_dist.c
MatView_SeqAIJ_Binary in src/mat/impls/aij/seq/aij.c
MatView_SeqAIJ_ASCII_structonly in src/mat/impls/aij/seq/aij.c
MatView_SeqAIJ_ASCII in src/mat/impls/aij/seq/aij.c
MatView_SeqAIJ_Draw_Zoom in src/mat/impls/aij/seq/aij.c
MatView_SeqAIJ_Draw in src/mat/impls/aij/seq/aij.c
MatView_SeqAIJ in src/mat/impls/aij/seq/aij.c
MatView_SeqAIJ_Inode in src/mat/impls/aij/seq/inode2.c
MatView_Info_KLU in src/mat/impls/aij/seq/klu/klu.c
MatView_KLU in src/mat/impls/aij/seq/klu/klu.c
MatView_Info_Matlab in src/mat/impls/aij/seq/matlab/aijmatlab.c
MatView_Matlab in src/mat/impls/aij/seq/matlab/aijmatlab.c
MatView_MKL_PARDISO in src/mat/impls/aij/seq/mkl_pardiso/mkl_pardiso.c
MatView_Info_SuperLU in src/mat/impls/aij/seq/superlu/superlu.c
MatView_SuperLU in src/mat/impls/aij/seq/superlu/superlu.c
MatView_Info_UMFPACK in src/mat/impls/aij/seq/umfpack/umfpack.c
MatView_UMFPACK in src/mat/impls/aij/seq/umfpack/umfpack.c
MatView_MPIBAIJ_ASCIIorDraworSocket in src/mat/impls/baij/mpi/mpibaij.c
MatView_MPIBAIJ_Binary in src/mat/impls/baij/mpi/mpibaij.c
MatView_MPIBAIJ in src/mat/impls/baij/mpi/mpibaij.c
MatView_SeqBAIJ_Binary in src/mat/impls/baij/seq/baij.c
MatView_SeqBAIJ_ASCII_structonly in src/mat/impls/baij/seq/baij.c
MatView_SeqBAIJ_ASCII in src/mat/impls/baij/seq/baij.c
MatView_SeqBAIJ_Draw_Zoom in src/mat/impls/baij/seq/baij.c
MatView_SeqBAIJ_Draw in src/mat/impls/baij/seq/baij.c
MatView_SeqBAIJ in src/mat/impls/baij/seq/baij.c
MatView_BlockMat in src/mat/impls/blockmat/seq/blockmat.c
MatView_MPIDense_Binary in src/mat/impls/dense/mpi/mpidense.c
MatView_MPIDense_ASCIIorDraworSocket in src/mat/impls/dense/mpi/mpidense.c
MatView_MPIDense in src/mat/impls/dense/mpi/mpidense.c
MatView_SeqDense_ASCII in src/mat/impls/dense/seq/dense.c
MatView_SeqDense_Binary in src/mat/impls/dense/seq/dense.c
MatView_SeqDense_Draw_Zoom in src/mat/impls/dense/seq/dense.c
MatView_SeqDense_Draw in src/mat/impls/dense/seq/dense.c
MatView_SeqDense in src/mat/impls/dense/seq/dense.c
MatView_Elemental in src/mat/impls/elemental/matelem.cxx
MatView_HYPRE in src/mat/impls/hypre/mhypre.c
MatView_IS in src/mat/impls/is/matis.c
MatView_SeqMAIJ in src/mat/impls/maij/maij.c
MatView_MPIMAIJ in src/mat/impls/maij/maij.c
MatView_MFFD in src/mat/impls/mffd/mffd.c
MatView_Nest in src/mat/impls/nest/matnest.c
MatView_Preallocator in src/mat/impls/preallocator/matpreallocator.c
MatView_MPISBAIJ_ASCIIorDraworSocket in src/mat/impls/sbaij/mpi/mpisbaij.c
MatView_MPISBAIJ_Binary in src/mat/impls/sbaij/mpi/mpisbaij.c
MatView_MPISBAIJ in src/mat/impls/sbaij/mpi/mpisbaij.c
MatView_Info_CHOLMOD in src/mat/impls/sbaij/seq/cholmod/sbaijcholmod.c
MatView_CHOLMOD in src/mat/impls/sbaij/seq/cholmod/sbaijcholmod.c
MatView_SeqSBAIJ_ASCII in src/mat/impls/sbaij/seq/sbaij.c
MatView_SeqSBAIJ_Draw_Zoom in src/mat/impls/sbaij/seq/sbaij.c
MatView_SeqSBAIJ_Draw in src/mat/impls/sbaij/seq/sbaij.c
MatView_SeqSBAIJ in src/mat/impls/sbaij/seq/sbaij.c
MatView_MPISELL_ASCIIorDraworSocket in src/mat/impls/sell/mpi/mpisell.c
MatView_MPISELL in src/mat/impls/sell/mpi/mpisell.c
MatView_SeqSELL_ASCII in src/mat/impls/sell/seq/sell.c
MatView_SeqSELL_Draw_Zoom in src/mat/impls/sell/seq/sell.c
MatView_SeqSELL_Draw in src/mat/impls/sell/seq/sell.c
MatView_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