Actual source code: ex6.c

  1: static char help[] = "Solves a tridiagonal linear system with KSP. \n\
  2: It illustrates how to do one symbolic factorization and multiple numeric factorizations using same matrix structure. \n\n";

  4: #include <petscksp.h>
  5: int main(int argc, char **args)
  6: {
  7:   Vec         x, b, u; /* approx solution, RHS, exact solution */
  8:   Mat         A;       /* linear system matrix */
  9:   KSP         ksp;     /* linear solver context */
 10:   PC          pc;      /* preconditioner context */
 11:   PetscReal   norm;    /* norm of solution error */
 12:   PetscInt    i, col[3], its, rstart, rend, N = 10, num_numfac;
 13:   PetscScalar value[3];

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

 19:   /* Create and assemble matrix. */
 20:   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
 21:   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, N, N));
 22:   PetscCall(MatSetFromOptions(A));
 23:   PetscCall(MatSetUp(A));
 24:   PetscCall(MatGetOwnershipRange(A, &rstart, &rend));

 26:   value[0] = -1.0;
 27:   value[1] = 2.0;
 28:   value[2] = -1.0;
 29:   for (i = rstart; i < rend; i++) {
 30:     col[0] = i - 1;
 31:     col[1] = i;
 32:     col[2] = i + 1;
 33:     if (i == 0) {
 34:       PetscCall(MatSetValues(A, 1, &i, 2, col + 1, value + 1, INSERT_VALUES));
 35:     } else if (i == N - 1) {
 36:       PetscCall(MatSetValues(A, 1, &i, 2, col, value, INSERT_VALUES));
 37:     } else {
 38:       PetscCall(MatSetValues(A, 1, &i, 3, col, value, INSERT_VALUES));
 39:     }
 40:   }
 41:   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
 42:   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
 43:   PetscCall(MatSetOption(A, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));

 45:   /* Create vectors */
 46:   PetscCall(MatCreateVecs(A, &x, &b));
 47:   PetscCall(VecDuplicate(x, &u));

 49:   /* Set exact solution; then compute right-hand-side vector. */
 50:   PetscCall(VecSet(u, 1.0));
 51:   PetscCall(MatMult(A, u, b));

 53:   /* Create the linear solver and set various options. */
 54:   PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp));
 55:   PetscCall(KSPGetPC(ksp, &pc));
 56:   PetscCall(PCSetType(pc, PCJACOBI));
 57:   PetscCall(KSPSetTolerances(ksp, 1.e-5, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
 58:   PetscCall(KSPSetOperators(ksp, A, A));
 59:   PetscCall(KSPSetFromOptions(ksp));

 61:   num_numfac = 1;
 62:   PetscCall(PetscOptionsGetInt(NULL, NULL, "-num_numfac", &num_numfac, NULL));
 63:   while (num_numfac--) {
 64:     /* An example on how to update matrix A for repeated numerical factorization and solve. */
 65:     PetscScalar one = 1.0;
 66:     PetscInt    i   = 0;
 67:     PetscCall(MatSetValues(A, 1, &i, 1, &i, &one, ADD_VALUES));
 68:     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
 69:     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
 70:     /* Update b */
 71:     PetscCall(MatMult(A, u, b));

 73:     /* Solve the linear system */
 74:     PetscCall(KSPSolve(ksp, b, x));

 76:     /* Check the solution and clean up */
 77:     PetscCall(VecAXPY(x, -1.0, u));
 78:     PetscCall(VecNorm(x, NORM_2, &norm));
 79:     PetscCall(KSPGetIterationNumber(ksp, &its));
 80:     if (norm > 100 * PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Norm of error %g, Iterations %" PetscInt_FMT "\n", (double)norm, its));
 81:   }

 83:   /* Free work space. */
 84:   PetscCall(VecDestroy(&x));
 85:   PetscCall(VecDestroy(&u));
 86:   PetscCall(VecDestroy(&b));
 87:   PetscCall(MatDestroy(&A));
 88:   PetscCall(KSPDestroy(&ksp));

 90:   PetscCall(PetscFinalize());
 91:   return 0;
 92: }

 94: /*TEST

 96:     test:
 97:       args: -num_numfac 2 -pc_type lu

 99:     test:
100:       suffix: 2
101:       args: -num_numfac 2 -pc_type lu -pc_factor_mat_solver_type mumps
102:       requires: mumps

104:     test:
105:       suffix: 3
106:       nsize: 3
107:       args: -num_numfac 2 -pc_type lu -pc_factor_mat_solver_type mumps
108:       requires: mumps

110: TEST*/