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