Actual source code: ex6.c

petsc-3.7.5 2017-01-01
Report Typos and Errors
  1: static char help[] = "Working out corner cases of the ASM preconditioner.\n\n";

  3: #include <petscksp.h>

  7: int main(int argc,char **args)
  8: {
  9:   KSP            ksp;
 10:   PC             pc;
 11:   Mat            A;
 12:   Vec            u, x, b;
 13:   PetscReal      error;
 14:   PetscMPIInt    rank, size, sized;
 15:   PetscInt       M = 8, N = 8, m, n, rstart, rend, r;
 16:   PetscBool      userSubdomains = PETSC_FALSE;

 19:   PetscInitialize(&argc, &args, NULL, help);
 20:   PetscOptionsGetInt(NULL,NULL, "-M", &M, NULL);
 21:   PetscOptionsGetInt(NULL,NULL, "-N", &N, NULL);
 22:   PetscOptionsGetBool(NULL,NULL, "-user_subdomains", &userSubdomains, NULL);
 23:   /* Do parallel decomposition */
 24:   MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
 25:   MPI_Comm_size(PETSC_COMM_WORLD, &size);
 26:   sized = (PetscMPIInt) PetscSqrtReal((PetscReal) size);
 27:   if (PetscSqr(sized) != size) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "This test may only be run on a nubmer of processes which is a perfect square, not %d", (int) size);
 28:   if (M % sized) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "The number of x-vertices %D does not divide the number of x-processes %d", M, (int) sized);
 29:   if (N % sized) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "The number of y-vertices %D does not divide the number of y-processes %d", N, (int) sized);
 30:   /* Assemble the matrix for the five point stencil, YET AGAIN
 31:        Every other process will be empty */
 32:   MatCreate(PETSC_COMM_WORLD, &A);
 33:   m    = (sized > 1) ? (rank % 2) ? 0 : 2*M/sized : M;
 34:   n    = N/sized;
 35:   MatSetSizes(A, m*n, m*n, M*N, M*N);
 36:   MatSetFromOptions(A);
 37:   MatSetUp(A);
 38:   MatGetOwnershipRange(A, &rstart, &rend);
 39:   for (r = rstart; r < rend; ++r) {
 40:     const PetscScalar diag = 4.0, offdiag = -1.0;
 41:     const PetscInt    i    = r/N;
 42:     const PetscInt    j    = r - i*N;
 43:     PetscInt          c;

 45:     if (i > 0)   {c = r - n; MatSetValues(A, 1, &r, 1, &c, &offdiag, INSERT_VALUES);}
 46:     if (i < M-1) {c = r + n; MatSetValues(A, 1, &r, 1, &c, &offdiag, INSERT_VALUES);}
 47:     if (j > 0)   {c = r - 1; MatSetValues(A, 1, &r, 1, &c, &offdiag, INSERT_VALUES);}
 48:     if (j < N-1) {c = r + 1; MatSetValues(A, 1, &r, 1, &c, &offdiag, INSERT_VALUES);}
 49:     MatSetValues(A, 1, &r, 1, &r, &diag, INSERT_VALUES);
 50:   }
 51:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 52:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 53:   /* Setup Solve */
 54:   VecCreate(PETSC_COMM_WORLD, &b);
 55:   VecSetSizes(b, m*n, PETSC_DETERMINE);
 56:   VecSetFromOptions(b);
 57:   VecDuplicate(b, &u);
 58:   VecDuplicate(b, &x);
 59:   VecSet(u, 1.0);
 60:   MatMult(A, u, b);
 61:   KSPCreate(PETSC_COMM_WORLD, &ksp);
 62:   KSPSetOperators(ksp, A, A);
 63:   KSPGetPC(ksp, &pc);
 64:   PCSetType(pc, PCASM);
 65:   /* Setup ASM by hand */
 66:   if (userSubdomains) {
 67:     IS        is;
 68:     PetscInt *rows;

 70:     /* Use no overlap for now */
 71:     PetscMalloc1(rend-rstart, &rows);
 72:     for (r = rstart; r < rend; ++r) rows[r-rstart] = r;
 73:     ISCreateGeneral(PETSC_COMM_SELF, rend-rstart, rows, PETSC_OWN_POINTER, &is);
 74:     PCASMSetLocalSubdomains(pc, 1, &is, &is);
 75:     ISDestroy(&is);
 76:   }
 77:   KSPSetFromOptions(ksp);
 78:   /* Solve and Compare */
 79:   KSPSolve(ksp, b, x);
 80:   VecAXPY(x, -1.0, u);
 81:   VecNorm(x, NORM_INFINITY, &error);
 82:   PetscPrintf(PETSC_COMM_WORLD, "Infinity norm of the error: %g\n", (double) error);
 83:   /* Cleanup */
 84:   KSPDestroy(&ksp);
 85:   MatDestroy(&A);
 86:   VecDestroy(&u);
 87:   VecDestroy(&x);
 88:   VecDestroy(&b);
 89:   PetscFinalize();
 90:   return 0;
 91: }