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 nonzero 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:   PetscReal   tol = PetscDefined(USE_REAL___FLOAT128) ? 1e-14 : 100. * PETSC_MACHINE_EPSILON;
 13:   PetscInt    i, col[3], its, rstart, rend, N = 10, num_numfac;
 14:   PetscScalar value[3];

 16:   PetscFunctionBeginUser;
 17:   PetscCall(PetscInitialize(&argc, &args, NULL, help));
 18:   PetscCall(PetscOptionsGetInt(NULL, NULL, "-N", &N, NULL));

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

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

 42:   /* Create vectors */
 43:   PetscCall(MatCreateVecs(A, &x, &b));
 44:   PetscCall(VecDuplicate(x, &u));

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

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

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

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

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

 80:   /* Free work space. */
 81:   PetscCall(VecDestroy(&x));
 82:   PetscCall(VecDestroy(&u));
 83:   PetscCall(VecDestroy(&b));
 84:   PetscCall(MatDestroy(&A));
 85:   PetscCall(KSPDestroy(&ksp));

 87:   PetscCall(PetscFinalize());
 88:   return 0;
 89: }

 91: /*TEST

 93:     test:
 94:       args: -num_numfac 2 -pc_type lu
 95:       output_file: output/empty.out

 97:     test:
 98:       suffix: 2
 99:       args: -num_numfac 2 -pc_type lu -pc_factor_mat_solver_type mumps
100:       requires: mumps
101:       output_file: output/empty.out

103:     test:
104:       suffix: 3
105:       nsize: 3
106:       args: -num_numfac 2 -pc_type lu -pc_factor_mat_solver_type mumps
107:       requires: mumps
108:       output_file: output/empty.out

110: TEST*/