Actual source code: ex50.c


  2: static char help[] = "Tests SeqBAIJ point block Jacobi for different block sizes\n\n";

  4: #include <petscksp.h>

  6: int main(int argc, char **args)
  7: {
  8:   Vec         x, b, u;
  9:   Mat         A;    /* linear system matrix */
 10:   KSP         ksp;  /* linear solver context */
 11:   PetscRandom rctx; /* random number generator context */
 12:   PetscReal   norm; /* norm of solution error */
 13:   PetscInt    i, j, k, l, n = 27, its, bs = 2, Ii, J;
 14:   PetscScalar v;

 17:   PetscInitialize(&argc, &args, (char *)0, help);
 18:   PetscOptionsGetInt(NULL, NULL, "-bs", &bs, NULL);
 19:   PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL);

 21:   MatCreate(PETSC_COMM_WORLD, &A);
 22:   MatSetSizes(A, n * bs, n * bs, PETSC_DETERMINE, PETSC_DETERMINE);
 23:   MatSetBlockSize(A, bs);
 24:   MatSetType(A, MATSEQBAIJ);
 25:   MatSetFromOptions(A);
 26:   MatSetUp(A);

 28:   /*
 29:      Don't bother to preallocate matrix
 30:   */
 31:   PetscRandomCreate(PETSC_COMM_SELF, &rctx);
 32:   for (i = 0; i < n; i++) {
 33:     for (j = 0; j < n; j++) {
 34:       PetscRandomGetValue(rctx, &v);
 35:       if (PetscRealPart(v) < .25 || i == j) {
 36:         for (k = 0; k < bs; k++) {
 37:           for (l = 0; l < bs; l++) {
 38:             PetscRandomGetValue(rctx, &v);
 39:             Ii = i * bs + k;
 40:             J  = j * bs + l;
 41:             if (Ii == J) v += 10.;
 42:             MatSetValue(A, Ii, J, v, INSERT_VALUES);
 43:           }
 44:         }
 45:       }
 46:     }
 47:   }

 49:   MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);
 50:   MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);

 52:   VecCreate(PETSC_COMM_WORLD, &u);
 53:   VecSetSizes(u, PETSC_DECIDE, n * bs);
 54:   VecSetFromOptions(u);
 55:   VecDuplicate(u, &b);
 56:   VecDuplicate(b, &x);

 58:   VecSet(u, 1.0);
 59:   MatMult(A, u, b);

 61:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 62:                 Create the linear solver and set various options
 63:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 65:   /*
 66:      Create linear solver context
 67:   */
 68:   KSPCreate(PETSC_COMM_WORLD, &ksp);

 70:   /*
 71:      Set operators. Here the matrix that defines the linear system
 72:      also serves as the preconditioning matrix.
 73:   */
 74:   KSPSetOperators(ksp, A, A);

 76:   KSPSetFromOptions(ksp);

 78:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 79:                       Solve the linear system
 80:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 82:   KSPSolve(ksp, b, x);

 84:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 85:                       Check solution and clean up
 86:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 88:   /*
 89:      Check the error
 90:   */
 91:   VecAXPY(x, -1.0, u);
 92:   VecNorm(x, NORM_2, &norm);
 93:   KSPGetIterationNumber(ksp, &its);

 95:   /*
 96:      Print convergence information.  PetscPrintf() produces a single
 97:      print statement from all processes that share a communicator.
 98:      An alternative is PetscFPrintf(), which prints to a file.
 99:   */
100:   if (norm > .1) PetscPrintf(PETSC_COMM_WORLD, "Norm of residual %g iterations %" PetscInt_FMT " bs %" PetscInt_FMT "\n", (double)norm, its, bs);

102:   /*
103:      Free work space.  All PETSc objects should be destroyed when they
104:      are no longer needed.
105:   */
106:   KSPDestroy(&ksp);
107:   VecDestroy(&u);
108:   VecDestroy(&x);
109:   VecDestroy(&b);
110:   MatDestroy(&A);
111:   PetscRandomDestroy(&rctx);

113:   /*
114:      Always call PetscFinalize() before exiting a program.  This routine
115:        - finalizes the PETSc libraries as well as MPI
116:        - provides summary and diagnostic information if certain runtime
117:          options are chosen (e.g., -log_view).
118:   */
119:   PetscFinalize();
120:   return 0;
121: }

123: /*TEST

125:    test:
126:       args: -bs {{1 2 3 4 5 6 7}} -pc_type pbjacobi

128:    test:
129:       suffix: 2
130:       args: -bs {{8 9 10 11 12 13 14 15}} -pc_type ilu

132: TEST*/