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