Actual source code: ex206.c

  1: static char help[] = "Reads binary matrix - twice\n";

  3: #include <petscmat.h>
  4: int main(int argc, char **args)
  5: {
  6:   Mat         A;
  7:   PetscViewer fd;                       /* viewer */
  8:   char        file[PETSC_MAX_PATH_LEN]; /* input file name */
  9:   PetscBool   flg;

 11:   PetscInitialize(&argc, &args, (char *)0, help);

 13:   PetscOptionsGetString(NULL, NULL, "-f", file, sizeof(file), &flg);

 16:   MatCreate(PETSC_COMM_WORLD, &A);
 17:   MatSetFromOptions(A);
 18:   PetscPrintf(PETSC_COMM_WORLD, "First MatLoad! \n");
 19:   PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &fd);
 20:   MatLoad(A, fd);
 21:   PetscViewerDestroy(&fd);
 22:   MatView(A, PETSC_VIEWER_STDOUT_WORLD);

 24:   PetscOptionsGetString(NULL, NULL, "-f2", file, sizeof(file), &flg);
 25:   PetscPrintf(PETSC_COMM_WORLD, "Second MatLoad! \n");
 26:   PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &fd);
 27:   MatLoad(A, fd);
 28:   PetscViewerDestroy(&fd);
 29:   MatView(A, PETSC_VIEWER_STDOUT_WORLD);

 31:   MatDestroy(&A);
 32:   PetscFinalize();
 33:   return 0;
 34: }

 36: /*TEST

 38:    test:
 39:       suffix: aij_1
 40:       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
 41:       args: -f ${DATAFILESPATH}/matrices/small -mat_type aij -mat_block_size 1
 42:       filter: grep -v Mat_

 44:    test:
 45:       suffix: aij_2
 46:       nsize: 2
 47:       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
 48:       args: -f ${DATAFILESPATH}/matrices/small -mat_type aij -mat_block_size 1
 49:       filter: grep -v Mat_

 51:    test:
 52:       suffix: aij_2_d
 53:       nsize: 2
 54:       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
 55:       args: -f ${DATAFILESPATH}/matrices/small -f2 ${DATAFILESPATH}/matrices/smallbs2 -mat_type aij -mat_block_size 1
 56:       filter: grep -v Mat_

 58:    test:
 59:       suffix: baij_1_2
 60:       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
 61:       args: -f ${DATAFILESPATH}/matrices/small -mat_type baij -mat_block_size 2
 62:       filter: grep -v Mat_

 64:    test:
 65:       suffix: baij_2_1_d
 66:       nsize: 2
 67:       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
 68:       args: -f ${DATAFILESPATH}/matrices/small -f2 ${DATAFILESPATH}/matrices/smallbs2 -mat_type baij -mat_block_size 1
 69:       filter: grep -v Mat_

 71:    test:
 72:       suffix: baij_2_2
 73:       nsize: 2
 74:       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
 75:       args: -f ${DATAFILESPATH}/matrices/small -mat_type baij -mat_block_size 2
 76:       filter: grep -v Mat_

 78:    test:
 79:       suffix: baij_2_2_d
 80:       nsize: 2
 81:       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
 82:       args: -f ${DATAFILESPATH}/matrices/small -f2 ${DATAFILESPATH}/matrices/smallbs2 -mat_type baij -mat_block_size 2
 83:       filter: grep -v Mat_

 85:    test:
 86:       suffix: sbaij_1_1
 87:       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
 88:       args: -f ${DATAFILESPATH}/matrices/small -mat_type sbaij -mat_block_size 1
 89:       filter: grep -v Mat_

 91:    test:
 92:       suffix: sbaij_1_2
 93:       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
 94:       args: -f ${DATAFILESPATH}/matrices/small -mat_type sbaij -mat_block_size 2
 95:       filter: grep -v Mat_

 97:    test:
 98:       suffix: sbaij_2_1_d
 99:       nsize: 2
100:       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
101:       args: -f ${DATAFILESPATH}/matrices/small -f2 ${DATAFILESPATH}/matrices/smallbs2 -mat_type sbaij -mat_block_size 1
102:       filter: grep -v Mat_

104:    test:
105:       suffix: sbaij_2_2
106:       nsize: 2
107:       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
108:       args: -f ${DATAFILESPATH}/matrices/small -mat_type sbaij -mat_block_size 2
109:       filter: grep -v Mat_

111: TEST*/