Actual source code: ex10.c


  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>

 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 */
 24:   PetscReal      *norms;
 25:   PetscInt       n,cstart,cend;
 26:   PetscBool      flg;
 27:   PetscViewerFormat format;

 29:   PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
 30:   /*
 31:      Determine files from which we read the matrix
 32:   */
 33:   PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&flg);
 34:   if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option");

 36:   /*
 37:      Open binary file.  Note that we use FILE_MODE_READ to indicate
 38:      reading from this file.
 39:   */
 40:   PetscViewerCreate(PETSC_COMM_WORLD,&fd);
 41:   PetscViewerSetType(fd,PETSCVIEWERBINARY);
 42:   PetscViewerSetFromOptions(fd);
 43:   PetscOptionsGetEnum(NULL,NULL,"-viewer_format",PetscViewerFormats,(PetscEnum*)&format,&flg);
 44:   if (flg) {PetscViewerPushFormat(fd,format);}
 45:   PetscViewerFileSetMode(fd,FILE_MODE_READ);
 46:   PetscViewerFileSetName(fd,file);

 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:   MatCreate(PETSC_COMM_WORLD,&A);
 54:   MatSetOptionsPrefix(A,"a_");
 55:   PetscObjectSetName((PetscObject) A,"A");
 56:   MatSetFromOptions(A);
 57:   MatLoad(A,fd);
 58:   PetscViewerDestroy(&fd);

 60:   MatGetSize(A,NULL,&n);
 61:   MatGetOwnershipRangeColumn(A,&cstart,&cend);
 62:   PetscMalloc1(n,&norms);
 63:   MatGetColumnNorms(A,NORM_2,norms);
 64:   PetscRealView(cend-cstart,norms+cstart,PETSC_VIEWER_STDOUT_WORLD);
 65:   PetscFree(norms);

 67:   PetscObjectPrintClassNamePrefixType((PetscObject)A,PETSC_VIEWER_STDOUT_WORLD);
 68:   MatGetOption(A,MAT_SYMMETRIC,&flg);
 69:   PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"MAT_SYMMETRIC: %D\n",flg);
 70:   MatViewFromOptions(A,NULL,"-mat_view");

 72:   MatDestroy(&A);
 73:   PetscFinalize();
 74:   return ierr;
 75: }

 77: /*TEST

 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

 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

 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

 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

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

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

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

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

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

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

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

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

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

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

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*/