Actual source code: ex15.c
2: static char help[] = "KSP linear solver on an operator with a null space.\n\n";
4: #include <petscksp.h>
6: int main(int argc, char **args)
7: {
8: Vec x, b, u; /* approx solution, RHS, exact solution */
9: Mat A; /* linear system matrix */
10: KSP ksp; /* KSP context */
11: PetscInt i, n = 10, col[3], its, i1, i2;
12: PetscScalar none = -1.0, value[3], avalue;
13: PetscReal norm;
14: PC pc;
17: PetscInitialize(&argc, &args, (char *)0, help);
18: PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL);
20: /* Create vectors */
21: VecCreate(PETSC_COMM_WORLD, &x);
22: VecSetSizes(x, PETSC_DECIDE, n);
23: VecSetFromOptions(x);
24: VecDuplicate(x, &b);
25: VecDuplicate(x, &u);
27: /* create a solution that is orthogonal to the constants */
28: VecGetOwnershipRange(u, &i1, &i2);
29: for (i = i1; i < i2; i++) {
30: avalue = i;
31: VecSetValues(u, 1, &i, &avalue, INSERT_VALUES);
32: }
33: VecAssemblyBegin(u);
34: VecAssemblyEnd(u);
35: VecSum(u, &avalue);
36: avalue = -avalue / (PetscReal)n;
37: VecShift(u, avalue);
39: /* Create and assemble matrix */
40: MatCreate(PETSC_COMM_WORLD, &A);
41: MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, n, n);
42: MatSetFromOptions(A);
43: value[0] = -1.0;
44: value[1] = 2.0;
45: value[2] = -1.0;
46: for (i = 1; i < n - 1; i++) {
47: col[0] = i - 1;
48: col[1] = i;
49: col[2] = i + 1;
50: MatSetValues(A, 1, &i, 3, col, value, INSERT_VALUES);
51: }
52: i = n - 1;
53: col[0] = n - 2;
54: col[1] = n - 1;
55: value[1] = 1.0;
56: MatSetValues(A, 1, &i, 2, col, value, INSERT_VALUES);
57: i = 0;
58: col[0] = 0;
59: col[1] = 1;
60: value[0] = 1.0;
61: value[1] = -1.0;
62: MatSetValues(A, 1, &i, 2, col, value, INSERT_VALUES);
63: MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);
64: MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);
65: MatMult(A, u, b);
67: /* Create KSP context; set operators and options; solve linear system */
68: KSPCreate(PETSC_COMM_WORLD, &ksp);
69: KSPSetOperators(ksp, A, A);
71: /* Insure that preconditioner has same null space as matrix */
72: /* currently does not do anything */
73: KSPGetPC(ksp, &pc);
75: KSPSetFromOptions(ksp);
76: KSPSolve(ksp, b, x);
77: /* KSPView(ksp,PETSC_VIEWER_STDOUT_WORLD); */
79: /* Check error */
80: VecAXPY(x, none, u);
81: VecNorm(x, NORM_2, &norm);
82: KSPGetIterationNumber(ksp, &its);
83: PetscPrintf(PETSC_COMM_WORLD, "Norm of error %g, Iterations %" PetscInt_FMT "\n", (double)norm, its);
85: /* Free work space */
86: VecDestroy(&x);
87: VecDestroy(&u);
88: VecDestroy(&b);
89: MatDestroy(&A);
90: KSPDestroy(&ksp);
91: PetscFinalize();
92: return 0;
93: }