Actual source code: ex234.c

  1: static char help[] = "Basic test of various routines with SBAIJ matrices\n\n";

  3: #include <petscmat.h>

  5: int main(int argc, char **argv)
  6: {
  7:   PetscInt    ia[3] = {0, 2, 4};
  8:   PetscInt    ja[4] = {0, 1, 0, 1};
  9:   PetscScalar c[16] = {0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15};
 10:   PetscMPIInt size;
 11:   Mat         ssbaij, msbaij;
 12:   Vec         x, y;

 15:   PetscInitialize(&argc, &argv, NULL, help);
 16:   MPI_Comm_size(PETSC_COMM_WORLD, &size);
 18:   MatCreate(PETSC_COMM_SELF, &ssbaij);
 19:   MatSetType(ssbaij, MATSEQSBAIJ);
 20:   MatSetBlockSize(ssbaij, 2);
 21:   MatSetSizes(ssbaij, 4, 8, 4, 8);
 22:   MatSeqSBAIJSetPreallocationCSR(ssbaij, 2, ia, ja, c);
 23:   MatCreateMPIMatConcatenateSeqMat(PETSC_COMM_WORLD, ssbaij, PETSC_DECIDE, MAT_INITIAL_MATRIX, &msbaij);
 24:   MatView(msbaij, PETSC_VIEWER_STDOUT_WORLD);
 25:   MatDestroy(&msbaij);
 26:   MatCreateMPIMatConcatenateSeqMat(PETSC_COMM_WORLD, ssbaij, 4, MAT_INITIAL_MATRIX, &msbaij);
 27:   MatView(msbaij, PETSC_VIEWER_STDOUT_WORLD);
 28:   MatCreateVecs(msbaij, &x, &y);
 29:   VecSet(x, 1);
 30:   MatMult(msbaij, x, y);
 31:   VecView(y, PETSC_VIEWER_STDOUT_WORLD);
 32:   MatMultAdd(msbaij, x, x, y);
 33:   VecView(y, PETSC_VIEWER_STDOUT_WORLD);
 34:   MatGetDiagonal(msbaij, y);
 35:   VecView(y, PETSC_VIEWER_STDOUT_WORLD);
 36:   VecDestroy(&x);
 37:   VecDestroy(&y);
 38:   MatDestroy(&msbaij);
 39:   MatDestroy(&ssbaij);
 40:   PetscFinalize();
 41:   return 0;
 42: }

 44: /*TEST

 46:    test:
 47:      nsize: 2
 48:      filter: sed "s?\.??g"

 50: TEST*/