xref: /petsc/src/mat/tutorials/ex10.c (revision 28b400f66ebc7ae0049166a2294dfcd3df27e64b)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Reads a PETSc matrix and computes the 2 norm of the columns\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown /*T
5c4762a1bSJed Brown    Concepts: Mat^loading a binary matrix;
6c4762a1bSJed Brown    Processors: n
7c4762a1bSJed Brown T*/
8c4762a1bSJed Brown 
9c4762a1bSJed Brown /*
10c4762a1bSJed Brown   Include "petscmat.h" so that we can use matrices.
11c4762a1bSJed Brown   automatically includes:
12c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h    - vectors
13c4762a1bSJed Brown      petscmat.h    - matrices
14c4762a1bSJed Brown      petscis.h     - index sets            petscviewer.h - viewers
15c4762a1bSJed Brown */
16c4762a1bSJed Brown #include <petscmat.h>
17c4762a1bSJed Brown 
18c4762a1bSJed Brown int main(int argc,char **args)
19c4762a1bSJed Brown {
20c4762a1bSJed Brown   Mat            A;                       /* matrix */
21c4762a1bSJed Brown   PetscViewer    fd;                      /* viewer */
22c4762a1bSJed Brown   char           file[PETSC_MAX_PATH_LEN];            /* input file name */
23c4762a1bSJed Brown   PetscErrorCode ierr;
24c4762a1bSJed Brown   PetscReal      *norms;
25c4762a1bSJed Brown   PetscInt       n,cstart,cend;
26c4762a1bSJed Brown   PetscBool      flg;
27c4762a1bSJed Brown   PetscViewerFormat format;
28c4762a1bSJed Brown 
29c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
30c4762a1bSJed Brown   /*
31c4762a1bSJed Brown      Determine files from which we read the matrix
32c4762a1bSJed Brown   */
335f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&flg));
34*28b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option");
35c4762a1bSJed Brown 
36c4762a1bSJed Brown   /*
37c4762a1bSJed Brown      Open binary file.  Note that we use FILE_MODE_READ to indicate
38c4762a1bSJed Brown      reading from this file.
39c4762a1bSJed Brown   */
405f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerCreate(PETSC_COMM_WORLD,&fd));
415f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerSetType(fd,PETSCVIEWERBINARY));
425f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerSetFromOptions(fd));
435f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetEnum(NULL,NULL,"-viewer_format",PetscViewerFormats,(PetscEnum*)&format,&flg));
445f80ce2aSJacob Faibussowitsch   if (flg) CHKERRQ(PetscViewerPushFormat(fd,format));
455f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerFileSetMode(fd,FILE_MODE_READ));
465f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerFileSetName(fd,file));
47c4762a1bSJed Brown 
48c4762a1bSJed Brown   /*
49c4762a1bSJed Brown     Load the matrix; then destroy the viewer.
50c4762a1bSJed Brown     Matrix type is set automatically but you can override it by MatSetType() prior to MatLoad().
51c4762a1bSJed Brown     Do that only if you really insist on the given type.
52c4762a1bSJed Brown   */
535f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
545f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOptionsPrefix(A,"a_"));
555f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) A,"A"));
565f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
575f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLoad(A,fd));
585f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&fd));
59c4762a1bSJed Brown 
605f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetSize(A,NULL,&n));
615f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOwnershipRangeColumn(A,&cstart,&cend));
625f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(n,&norms));
635f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetColumnNorms(A,NORM_2,norms));
645f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRealView(cend-cstart,norms+cstart,PETSC_VIEWER_STDOUT_WORLD));
655f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(norms));
66c4762a1bSJed Brown 
675f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectPrintClassNamePrefixType((PetscObject)A,PETSC_VIEWER_STDOUT_WORLD));
685f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOption(A,MAT_SYMMETRIC,&flg));
695f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"MAT_SYMMETRIC: %" PetscInt_FMT "\n",(PetscInt)flg));
705f80ce2aSJacob Faibussowitsch   CHKERRQ(MatViewFromOptions(A,NULL,"-mat_view"));
71c4762a1bSJed Brown 
725f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
73c4762a1bSJed Brown   ierr = PetscFinalize();
74c4762a1bSJed Brown   return ierr;
75c4762a1bSJed Brown }
76c4762a1bSJed Brown 
77c4762a1bSJed Brown /*TEST
78c4762a1bSJed Brown 
79c4762a1bSJed Brown    test:
80c4762a1bSJed Brown       suffix: mpiaij
81c4762a1bSJed Brown       nsize: 2
82dfd57a17SPierre Jolivet       requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES)
83c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/small -a_mat_type mpiaij
84c4762a1bSJed Brown       args: -a_matload_symmetric
85c4762a1bSJed Brown 
86c4762a1bSJed Brown    test:
87c4762a1bSJed Brown       suffix: mpiaij_hdf5
88c4762a1bSJed Brown       nsize: 2
89dfd57a17SPierre Jolivet       requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) hdf5 defined(PETSC_HDF5_HAVE_ZLIB)
90c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/matlab/small.mat -a_mat_type mpiaij -viewer_type hdf5 -viewer_format hdf5_mat
91c4762a1bSJed Brown       args: -a_matload_symmetric
92c4762a1bSJed Brown 
93c4762a1bSJed Brown    test:
94c4762a1bSJed Brown       suffix: mpiaij_rect_hdf5
95c4762a1bSJed Brown       nsize: 2
96dfd57a17SPierre Jolivet       requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) hdf5 defined(PETSC_HDF5_HAVE_ZLIB)
97c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/matlab/small_rect.mat -a_mat_type mpiaij -viewer_type hdf5 -viewer_format hdf5_mat
98c4762a1bSJed Brown 
99c4762a1bSJed Brown    test:
100c4762a1bSJed Brown       # test for more processes than rows
101c4762a1bSJed Brown       suffix: mpiaij_hdf5_tiny
102c4762a1bSJed Brown       nsize: 8
103dfd57a17SPierre Jolivet       requires: double !complex !defined(PETSC_USE_64BIT_INDICES) hdf5 defined(PETSC_HDF5_HAVE_ZLIB)
104c4762a1bSJed Brown       args: -f ${wPETSC_DIR}/share/petsc/datafiles/matrices/tiny_system_with_x0.mat -a_mat_type mpiaij -viewer_type hdf5 -viewer_format hdf5_mat
105c4762a1bSJed Brown       args: -a_matload_symmetric
106c4762a1bSJed Brown 
107c4762a1bSJed Brown    test:
108c4762a1bSJed Brown       # test for more processes than rows, complex
109c4762a1bSJed Brown       TODO: not yet implemented for MATLAB complex format
110c4762a1bSJed Brown       suffix: mpiaij_hdf5_tiny_complex
111c4762a1bSJed Brown       nsize: 8
112dfd57a17SPierre Jolivet       requires: double complex !defined(PETSC_USE_64BIT_INDICES) hdf5 defined(PETSC_HDF5_HAVE_ZLIB)
113c4762a1bSJed Brown       args: -f ${wPETSC_DIR}/share/petsc/datafiles/matrices/tiny_system_with_x0_complex.mat -a_mat_type mpiaij -viewer_type hdf5 -viewer_format hdf5_mat
114c4762a1bSJed Brown       args: -a_matload_symmetric
115c4762a1bSJed Brown 
116c4762a1bSJed Brown    test:
117c4762a1bSJed Brown       TODO: mpibaij not supported yet
118c4762a1bSJed Brown       suffix: mpibaij_hdf5
119c4762a1bSJed Brown       nsize: 2
120dfd57a17SPierre Jolivet       requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) hdf5 defined(PETSC_HDF5_HAVE_ZLIB)
121c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/matlab/small.mat -a_mat_type mpibaij -a_mat_block_size 2 -viewer_type hdf5 -viewer_format hdf5_mat
122c4762a1bSJed Brown       args: -a_matload_symmetric
123c4762a1bSJed Brown 
124c4762a1bSJed Brown    test:
125c4762a1bSJed Brown       suffix: mpidense
126c4762a1bSJed Brown       nsize: 2
127dfd57a17SPierre Jolivet       requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES)
128c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/small -a_mat_type mpidense
129c4762a1bSJed Brown       args: -a_matload_symmetric
130c4762a1bSJed Brown 
131c4762a1bSJed Brown    test:
132c4762a1bSJed Brown       suffix: seqaij
133dfd57a17SPierre Jolivet       requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES)
134c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/small -a_mat_type seqaij
135c4762a1bSJed Brown       args: -a_matload_symmetric
136c4762a1bSJed Brown 
137c4762a1bSJed Brown    test:
138c4762a1bSJed Brown       suffix: seqaij_hdf5
139dfd57a17SPierre Jolivet       requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) hdf5 defined(PETSC_HDF5_HAVE_ZLIB)
140c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/matlab/small.mat -a_mat_type seqaij -viewer_type hdf5 -viewer_format hdf5_mat
141c4762a1bSJed Brown       args: -a_matload_symmetric
142c4762a1bSJed Brown 
143c4762a1bSJed Brown    test:
144c4762a1bSJed Brown       suffix: seqaij_rect_hdf5
145dfd57a17SPierre Jolivet       requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) hdf5 defined(PETSC_HDF5_HAVE_ZLIB)
146c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/matlab/small_rect.mat -a_mat_type seqaij -viewer_type hdf5 -viewer_format hdf5_mat
147c4762a1bSJed Brown 
148c4762a1bSJed Brown    test:
149c4762a1bSJed Brown       suffix: seqdense
150dfd57a17SPierre Jolivet       requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES)
151c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/small -a_mat_type seqdense
152c4762a1bSJed Brown       args: -a_matload_symmetric
153c4762a1bSJed Brown 
154c4762a1bSJed Brown    test:
155c4762a1bSJed Brown       suffix: seqdense_hdf5
156dfd57a17SPierre Jolivet       requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) hdf5 defined(PETSC_HDF5_HAVE_ZLIB)
157c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/matlab/small_dense.mat -a_mat_type seqdense -viewer_type hdf5 -viewer_format hdf5_mat
158c4762a1bSJed Brown       args: -a_matload_symmetric
159c4762a1bSJed Brown 
160c4762a1bSJed Brown    test:
161c4762a1bSJed Brown       suffix: seqdense_rect_hdf5
162dfd57a17SPierre Jolivet       requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) hdf5 defined(PETSC_HDF5_HAVE_ZLIB)
163c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/matlab/small_rect_dense.mat -a_mat_type seqdense -viewer_type hdf5 -viewer_format hdf5_mat
164c4762a1bSJed Brown 
165c4762a1bSJed Brown    test:
166c4762a1bSJed Brown       suffix: mpidense_hdf5
167c4762a1bSJed Brown       nsize: 2
168dfd57a17SPierre Jolivet       requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) hdf5 defined(PETSC_HDF5_HAVE_ZLIB)
169c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/matlab/small_dense.mat -a_mat_type mpidense -viewer_type hdf5 -viewer_format hdf5_mat
170c4762a1bSJed Brown       args: -a_matload_symmetric
171c4762a1bSJed Brown 
172c4762a1bSJed Brown    test:
173c4762a1bSJed Brown       suffix: mpidense_rect_hdf5
174c4762a1bSJed Brown       nsize: 2
175dfd57a17SPierre Jolivet       requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) hdf5 defined(PETSC_HDF5_HAVE_ZLIB)
176c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/matlab/small_rect_dense.mat -a_mat_type mpidense -viewer_type hdf5 -viewer_format hdf5_mat
177c4762a1bSJed Brown TEST*/
178