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