Actual source code: ex100.c


  2: static char help[] = "Tests various routines in MatMAIJ format.\n";

  4: #include <petscmat.h>
  5: #define IMAX 15
  6: int main(int argc, char **args)
  7: {
  8:   Mat         A, B, MA;
  9:   PetscViewer fd;
 10:   char        file[PETSC_MAX_PATH_LEN];
 11:   PetscInt    m, n, M, N, dof = 1;
 12:   PetscMPIInt rank, size;
 13:   PetscBool   flg;

 16:   PetscInitialize(&argc, &args, (char *)0, help);
 17:   MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
 18:   MPI_Comm_size(PETSC_COMM_WORLD, &size);

 20:   /* Load aij matrix A */
 21:   PetscOptionsGetString(NULL, NULL, "-f", file, sizeof(file), &flg);
 23:   PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &fd);
 24:   MatCreate(PETSC_COMM_WORLD, &A);
 25:   MatLoad(A, fd);
 26:   PetscViewerDestroy(&fd);

 28:   /* Get dof, then create maij matrix MA */
 29:   PetscOptionsGetInt(NULL, NULL, "-dof", &dof, NULL);
 30:   MatCreateMAIJ(A, dof, &MA);
 31:   MatGetLocalSize(MA, &m, &n);
 32:   MatGetSize(MA, &M, &N);

 34:   if (size == 1) {
 35:     MatConvert(MA, MATSEQAIJ, MAT_INITIAL_MATRIX, &B);
 36:   } else {
 37:     MatConvert(MA, MATMPIAIJ, MAT_INITIAL_MATRIX, &B);
 38:   }

 40:   /* Test MatMult() */
 41:   MatMultEqual(MA, B, 10, &flg);
 43:   /* Test MatMultAdd() */
 44:   MatMultAddEqual(MA, B, 10, &flg);

 47:   /* Test MatMultTranspose() */
 48:   MatMultTransposeEqual(MA, B, 10, &flg);

 51:   /* Test MatMultTransposeAdd() */
 52:   MatMultTransposeAddEqual(MA, B, 10, &flg);

 55:   MatDestroy(&MA);
 56:   MatDestroy(&A);
 57:   MatDestroy(&B);
 58:   PetscFinalize();
 59:   return 0;
 60: }

 62: /*TEST

 64:    build:
 65:       requires: !complex

 67:    test:
 68:       nsize: {{1 3}}
 69:       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
 70:       args: -f ${DATAFILESPATH}/matrices/arco1 -dof {{1 2 3 4 5 6 8 9 16}} -viewer_binary_skip_info

 72: TEST*/