Actual source code: ex60.c


  2: static char help[] = "Tests MatGetColumnVector().";

  4: #include <petscmat.h>

  6: int main(int argc, char **args)
  7: {
  8:   Mat         C;
  9:   PetscInt    i, j, m = 3, n = 2, Ii, J, col = 0;
 10:   PetscMPIInt size, rank;
 11:   PetscScalar v;
 12:   Vec         yy;

 14:   PetscFunctionBeginUser;
 15:   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
 16:   PetscCall(PetscOptionsGetInt(NULL, NULL, "-col", &col, NULL));

 18:   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
 19:   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
 20:   n = 2 * size;

 22:   /* create the matrix for the five point stencil, YET AGAIN*/
 23:   PetscCall(MatCreate(PETSC_COMM_WORLD, &C));
 24:   PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, m * n, m * n));
 25:   PetscCall(MatSetFromOptions(C));
 26:   PetscCall(MatSeqAIJSetPreallocation(C, 5, NULL));
 27:   PetscCall(MatMPIAIJSetPreallocation(C, 5, NULL, 5, NULL));
 28:   PetscCall(MatSetUp(C));

 30:   for (i = 0; i < m; i++) {
 31:     for (j = 2 * rank; j < 2 * rank + 2; j++) {
 32:       v  = -1.0;
 33:       Ii = j + n * i;
 34:       if (i > 0) {
 35:         J = Ii - n;
 36:         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
 37:       }
 38:       if (i < m - 1) {
 39:         J = Ii + n;
 40:         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
 41:       }
 42:       if (j > 0) {
 43:         J = Ii - 1;
 44:         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
 45:       }
 46:       if (j < n - 1) {
 47:         J = Ii + 1;
 48:         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
 49:       }
 50:       v = 4.0;
 51:       PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, INSERT_VALUES));
 52:     }
 53:   }
 54:   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
 55:   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
 56:   PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD));

 58:   PetscCall(MatCreateVecs(C, NULL, &yy));
 59:   PetscCall(VecSetFromOptions(yy));

 61:   PetscCall(MatGetColumnVector(C, yy, col));

 63:   PetscCall(VecView(yy, PETSC_VIEWER_STDOUT_WORLD));

 65:   PetscCall(VecDestroy(&yy));
 66:   PetscCall(MatDestroy(&C));
 67:   PetscCall(PetscFinalize());
 68:   return 0;
 69: }

 71: /*TEST

 73:    test:
 74:       nsize: 3
 75:       args: -col 7

 77:    test:
 78:       suffix: dense
 79:       nsize: 3
 80:       args: -col 7 -mat_type dense -vec_type {{mpi standard}}
 81:       filter: grep -v type

 83:    test:
 84:       requires: cuda
 85:       suffix: dense_cuda
 86:       nsize: 3
 87:       output_file: output/ex60_dense.out
 88:       args: -col 7 -mat_type {{mpidense mpidensecuda}} -vec_type {{mpi standard cuda mpicuda}}
 89:       filter: grep -v type

 91: TEST*/