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