xref: /petsc/src/mat/tutorials/ex10.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254) !
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   PetscReal      *norms;
24c4762a1bSJed Brown   PetscInt       n,cstart,cend;
25c4762a1bSJed Brown   PetscBool      flg;
26c4762a1bSJed Brown   PetscViewerFormat format;
27c4762a1bSJed Brown 
28*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
29c4762a1bSJed Brown   /*
30c4762a1bSJed Brown      Determine files from which we read the matrix
31c4762a1bSJed Brown   */
325f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&flg));
3328b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option");
34c4762a1bSJed Brown 
35c4762a1bSJed Brown   /*
36c4762a1bSJed Brown      Open binary file.  Note that we use FILE_MODE_READ to indicate
37c4762a1bSJed Brown      reading from this file.
38c4762a1bSJed Brown   */
395f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerCreate(PETSC_COMM_WORLD,&fd));
405f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerSetType(fd,PETSCVIEWERBINARY));
415f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerSetFromOptions(fd));
425f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetEnum(NULL,NULL,"-viewer_format",PetscViewerFormats,(PetscEnum*)&format,&flg));
435f80ce2aSJacob Faibussowitsch   if (flg) CHKERRQ(PetscViewerPushFormat(fd,format));
445f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerFileSetMode(fd,FILE_MODE_READ));
455f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerFileSetName(fd,file));
46c4762a1bSJed Brown 
47c4762a1bSJed Brown   /*
48c4762a1bSJed Brown     Load the matrix; then destroy the viewer.
49c4762a1bSJed Brown     Matrix type is set automatically but you can override it by MatSetType() prior to MatLoad().
50c4762a1bSJed Brown     Do that only if you really insist on the given type.
51c4762a1bSJed Brown   */
525f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
535f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOptionsPrefix(A,"a_"));
545f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) A,"A"));
555f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
565f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLoad(A,fd));
575f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&fd));
58c4762a1bSJed Brown 
595f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetSize(A,NULL,&n));
605f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOwnershipRangeColumn(A,&cstart,&cend));
615f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(n,&norms));
625f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetColumnNorms(A,NORM_2,norms));
635f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRealView(cend-cstart,norms+cstart,PETSC_VIEWER_STDOUT_WORLD));
645f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(norms));
65c4762a1bSJed Brown 
665f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectPrintClassNamePrefixType((PetscObject)A,PETSC_VIEWER_STDOUT_WORLD));
675f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOption(A,MAT_SYMMETRIC,&flg));
685f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"MAT_SYMMETRIC: %" PetscInt_FMT "\n",(PetscInt)flg));
695f80ce2aSJacob Faibussowitsch   CHKERRQ(MatViewFromOptions(A,NULL,"-mat_view"));
70c4762a1bSJed Brown 
715f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
72*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
73*b122ec5aSJacob Faibussowitsch   return 0;
74c4762a1bSJed Brown }
75c4762a1bSJed Brown 
76c4762a1bSJed Brown /*TEST
77c4762a1bSJed Brown 
78c4762a1bSJed Brown    test:
79c4762a1bSJed Brown       suffix: mpiaij
80c4762a1bSJed Brown       nsize: 2
81dfd57a17SPierre Jolivet       requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES)
82c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/small -a_mat_type mpiaij
83c4762a1bSJed Brown       args: -a_matload_symmetric
84c4762a1bSJed Brown 
85c4762a1bSJed Brown    test:
86c4762a1bSJed Brown       suffix: mpiaij_hdf5
87c4762a1bSJed Brown       nsize: 2
88dfd57a17SPierre Jolivet       requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) hdf5 defined(PETSC_HDF5_HAVE_ZLIB)
89c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/matlab/small.mat -a_mat_type mpiaij -viewer_type hdf5 -viewer_format hdf5_mat
90c4762a1bSJed Brown       args: -a_matload_symmetric
91c4762a1bSJed Brown 
92c4762a1bSJed Brown    test:
93c4762a1bSJed Brown       suffix: mpiaij_rect_hdf5
94c4762a1bSJed Brown       nsize: 2
95dfd57a17SPierre Jolivet       requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) hdf5 defined(PETSC_HDF5_HAVE_ZLIB)
96c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/matlab/small_rect.mat -a_mat_type mpiaij -viewer_type hdf5 -viewer_format hdf5_mat
97c4762a1bSJed Brown 
98c4762a1bSJed Brown    test:
99c4762a1bSJed Brown       # test for more processes than rows
100c4762a1bSJed Brown       suffix: mpiaij_hdf5_tiny
101c4762a1bSJed Brown       nsize: 8
102dfd57a17SPierre Jolivet       requires: double !complex !defined(PETSC_USE_64BIT_INDICES) hdf5 defined(PETSC_HDF5_HAVE_ZLIB)
103c4762a1bSJed 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
104c4762a1bSJed Brown       args: -a_matload_symmetric
105c4762a1bSJed Brown 
106c4762a1bSJed Brown    test:
107c4762a1bSJed Brown       # test for more processes than rows, complex
108c4762a1bSJed Brown       TODO: not yet implemented for MATLAB complex format
109c4762a1bSJed Brown       suffix: mpiaij_hdf5_tiny_complex
110c4762a1bSJed Brown       nsize: 8
111dfd57a17SPierre Jolivet       requires: double complex !defined(PETSC_USE_64BIT_INDICES) hdf5 defined(PETSC_HDF5_HAVE_ZLIB)
112c4762a1bSJed 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
113c4762a1bSJed Brown       args: -a_matload_symmetric
114c4762a1bSJed Brown 
115c4762a1bSJed Brown    test:
116c4762a1bSJed Brown       TODO: mpibaij not supported yet
117c4762a1bSJed Brown       suffix: mpibaij_hdf5
118c4762a1bSJed Brown       nsize: 2
119dfd57a17SPierre Jolivet       requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) hdf5 defined(PETSC_HDF5_HAVE_ZLIB)
120c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/matlab/small.mat -a_mat_type mpibaij -a_mat_block_size 2 -viewer_type hdf5 -viewer_format hdf5_mat
121c4762a1bSJed Brown       args: -a_matload_symmetric
122c4762a1bSJed Brown 
123c4762a1bSJed Brown    test:
124c4762a1bSJed Brown       suffix: mpidense
125c4762a1bSJed Brown       nsize: 2
126dfd57a17SPierre Jolivet       requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES)
127c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/small -a_mat_type mpidense
128c4762a1bSJed Brown       args: -a_matload_symmetric
129c4762a1bSJed Brown 
130c4762a1bSJed Brown    test:
131c4762a1bSJed Brown       suffix: seqaij
132dfd57a17SPierre Jolivet       requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES)
133c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/small -a_mat_type seqaij
134c4762a1bSJed Brown       args: -a_matload_symmetric
135c4762a1bSJed Brown 
136c4762a1bSJed Brown    test:
137c4762a1bSJed Brown       suffix: seqaij_hdf5
138dfd57a17SPierre Jolivet       requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) hdf5 defined(PETSC_HDF5_HAVE_ZLIB)
139c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/matlab/small.mat -a_mat_type seqaij -viewer_type hdf5 -viewer_format hdf5_mat
140c4762a1bSJed Brown       args: -a_matload_symmetric
141c4762a1bSJed Brown 
142c4762a1bSJed Brown    test:
143c4762a1bSJed Brown       suffix: seqaij_rect_hdf5
144dfd57a17SPierre Jolivet       requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) hdf5 defined(PETSC_HDF5_HAVE_ZLIB)
145c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/matlab/small_rect.mat -a_mat_type seqaij -viewer_type hdf5 -viewer_format hdf5_mat
146c4762a1bSJed Brown 
147c4762a1bSJed Brown    test:
148c4762a1bSJed Brown       suffix: seqdense
149dfd57a17SPierre Jolivet       requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES)
150c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/small -a_mat_type seqdense
151c4762a1bSJed Brown       args: -a_matload_symmetric
152c4762a1bSJed Brown 
153c4762a1bSJed Brown    test:
154c4762a1bSJed Brown       suffix: seqdense_hdf5
155dfd57a17SPierre Jolivet       requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) hdf5 defined(PETSC_HDF5_HAVE_ZLIB)
156c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/matlab/small_dense.mat -a_mat_type seqdense -viewer_type hdf5 -viewer_format hdf5_mat
157c4762a1bSJed Brown       args: -a_matload_symmetric
158c4762a1bSJed Brown 
159c4762a1bSJed Brown    test:
160c4762a1bSJed Brown       suffix: seqdense_rect_hdf5
161dfd57a17SPierre Jolivet       requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) hdf5 defined(PETSC_HDF5_HAVE_ZLIB)
162c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/matlab/small_rect_dense.mat -a_mat_type seqdense -viewer_type hdf5 -viewer_format hdf5_mat
163c4762a1bSJed Brown 
164c4762a1bSJed Brown    test:
165c4762a1bSJed Brown       suffix: mpidense_hdf5
166c4762a1bSJed Brown       nsize: 2
167dfd57a17SPierre Jolivet       requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) hdf5 defined(PETSC_HDF5_HAVE_ZLIB)
168c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/matlab/small_dense.mat -a_mat_type mpidense -viewer_type hdf5 -viewer_format hdf5_mat
169c4762a1bSJed Brown       args: -a_matload_symmetric
170c4762a1bSJed Brown 
171c4762a1bSJed Brown    test:
172c4762a1bSJed Brown       suffix: mpidense_rect_hdf5
173c4762a1bSJed Brown       nsize: 2
174dfd57a17SPierre Jolivet       requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) hdf5 defined(PETSC_HDF5_HAVE_ZLIB)
175c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/matlab/small_rect_dense.mat -a_mat_type mpidense -viewer_type hdf5 -viewer_format hdf5_mat
176c4762a1bSJed Brown TEST*/
177