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;
16: PetscInitialize(&argc,&args,(char*)0,help);
17: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
19: /* Create vectors */
20: VecCreate(PETSC_COMM_WORLD,&x);
21: VecSetSizes(x,PETSC_DECIDE,n);
22: VecSetFromOptions(x);
23: VecDuplicate(x,&b);
24: VecDuplicate(x,&u);
26: /* create a solution that is orthogonal to the constants */
27: VecGetOwnershipRange(u,&i1,&i2);
28: for (i=i1; i<i2; i++) {
29: avalue = i;
30: VecSetValues(u,1,&i,&avalue,INSERT_VALUES);
31: }
32: VecAssemblyBegin(u);
33: VecAssemblyEnd(u);
34: VecSum(u,&avalue);
35: avalue = -avalue/(PetscReal)n;
36: VecShift(u,avalue);
38: /* Create and assemble matrix */
39: MatCreate(PETSC_COMM_WORLD,&A);
40: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
41: MatSetFromOptions(A);
42: value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
43: for (i=1; i<n-1; i++) {
44: col[0] = i-1; col[1] = i; col[2] = i+1;
45: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
46: }
47: i = n - 1; col[0] = n - 2; col[1] = n - 1; value[1] = 1.0;
48: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
49: i = 0; col[0] = 0; col[1] = 1; value[0] = 1.0; value[1] = -1.0;
50: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
51: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
52: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
53: MatMult(A,u,b);
55: /* Create KSP context; set operators and options; solve linear system */
56: KSPCreate(PETSC_COMM_WORLD,&ksp);
57: KSPSetOperators(ksp,A,A);
59: /* Insure that preconditioner has same null space as matrix */
60: /* currently does not do anything */
61: KSPGetPC(ksp,&pc);
63: KSPSetFromOptions(ksp);
64: KSPSolve(ksp,b,x);
65: /* KSPView(ksp,PETSC_VIEWER_STDOUT_WORLD); */
67: /* Check error */
68: VecAXPY(x,none,u);
69: VecNorm(x,NORM_2,&norm);
70: KSPGetIterationNumber(ksp,&its);
71: PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g, Iterations %D\n",(double)norm,its);
73: /* Free work space */
74: VecDestroy(&x));PetscCall(VecDestroy(&u);
75: VecDestroy(&b));PetscCall(MatDestroy(&A);
76: KSPDestroy(&ksp);
77: PetscFinalize();
78: return 0;
79: }