Actual source code: ex47.c


  2: static char help[] = "Tests the various routines in MatBAIJ format.\n\
  3: Input arguments are:\n\
  4:   -f <input_file> : file to load.  For example see $PETSC_DIR/share/petsc/datafiles/matrices\n\n";

  6: #include <petscmat.h>

  8: int main(int argc, char **args)
  9: {
 10:   Mat                A, B, C;
 11:   PetscViewer        va, vb, vc;
 12:   Vec                x, y;
 13:   PetscInt           i, j, row, m, n, ncols1, ncols2, ct, m2, n2;
 14:   const PetscInt    *cols1, *cols2;
 15:   char               file[PETSC_MAX_PATH_LEN];
 16:   PetscBool          tflg;
 17:   PetscScalar        rval;
 18:   const PetscScalar *vals1, *vals2;
 19:   PetscReal          norm1, norm2, rnorm;
 20:   PetscRandom        r;

 23:   PetscInitialize(&argc, &args, (char *)0, help);
 24:   PetscOptionsGetString(NULL, NULL, "-f", file, sizeof(file), NULL);

 26:   /* Load the matrix as AIJ format */
 27:   PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &va);
 28:   MatCreate(PETSC_COMM_WORLD, &A);
 29:   MatSetType(A, MATSEQAIJ);
 30:   MatLoad(A, va);
 31:   PetscViewerDestroy(&va);

 33:   /* Load the matrix as BAIJ format */
 34:   PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &vb);
 35:   MatCreate(PETSC_COMM_WORLD, &B);
 36:   MatSetType(B, MATSEQBAIJ);
 37:   MatLoad(B, vb);
 38:   PetscViewerDestroy(&vb);

 40:   /* Load the matrix as BAIJ format */
 41:   PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &vc);
 42:   MatCreate(PETSC_COMM_WORLD, &C);
 43:   MatSetType(C, MATSEQBAIJ);
 44:   MatSetFromOptions(C);
 45:   MatLoad(C, vc);
 46:   PetscViewerDestroy(&vc);

 48:   MatGetSize(A, &m, &n);
 49:   MatGetSize(B, &m2, &n2);

 52:   /* Test MatEqual() */
 53:   MatEqual(B, C, &tflg);

 56:   /* Test MatGetDiagonal() */
 57:   VecCreateSeq(PETSC_COMM_SELF, m, &x);
 58:   VecCreateSeq(PETSC_COMM_SELF, m, &y);

 60:   MatGetDiagonal(A, x);
 61:   MatGetDiagonal(B, y);

 63:   VecEqual(x, y, &tflg);

 66:   /* Test MatDiagonalScale() */
 67:   PetscRandomCreate(PETSC_COMM_SELF, &r);
 68:   PetscRandomSetFromOptions(r);
 69:   VecSetRandom(x, r);
 70:   VecSetRandom(y, r);

 72:   MatDiagonalScale(A, x, y);
 73:   MatDiagonalScale(B, x, y);
 74:   MatMult(A, x, y);
 75:   VecNorm(y, NORM_2, &norm1);
 76:   MatMult(B, x, y);
 77:   VecNorm(y, NORM_2, &norm2);
 78:   rnorm = ((norm1 - norm2) * 100) / norm1;

 81:   /* Test MatGetRow()/ MatRestoreRow() */
 82:   for (ct = 0; ct < 100; ct++) {
 83:     PetscRandomGetValue(r, &rval);
 84:     row = (int)(rval * m);
 85:     MatGetRow(A, row, &ncols1, &cols1, &vals1);
 86:     MatGetRow(B, row, &ncols2, &cols2, &vals2);

 88:     for (i = 0, j = 0; i < ncols1 && j < ncols2; i++) {
 89:       while (cols2[j] != cols1[i]) j++;
 91:     }

 94:     MatRestoreRow(A, row, &ncols1, &cols1, &vals1);
 95:     MatRestoreRow(B, row, &ncols2, &cols2, &vals2);
 96:   }

 98:   MatDestroy(&A);
 99:   MatDestroy(&B);
100:   MatDestroy(&C);
101:   VecDestroy(&x);
102:   VecDestroy(&y);
103:   PetscRandomDestroy(&r);
104:   PetscFinalize();
105:   return 0;
106: }

108: /*TEST

110:    build:
111:       requires: !complex

113:    test:
114:       args: -f ${DATAFILESPATH}/matrices/cfd.1.10 -mat_block_size 5
115:       requires: !complex double datafilespath !defined(PETSC_USE_64BIT_INDICES)

117: TEST*/