Actual source code: ex15.c

  2: static char help[] = "KSP on an operator with a null space.\n\n";

 4:  #include petscksp.h

  8: int main(int argc,char **args)
  9: {
 10:   Vec            x,b,u;      /* approx solution, RHS, exact solution */
 11:   Mat            A;            /* linear system matrix */
 12:   KSP            ksp;         /* KSP context */
 14:   PetscInt       i,n = 10,col[3],its,i1,i2;
 15:   PetscScalar    none = -1.0,value[3],avalue;
 16:   PetscReal      norm;
 17:   PC             pc;

 19:   PetscInitialize(&argc,&args,(char *)0,help);
 20:   PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);

 22:   /* Create vectors */
 23:   VecCreate(PETSC_COMM_WORLD,&x);
 24:   VecSetSizes(x,PETSC_DECIDE,n);
 25:   VecSetFromOptions(x);
 26:   VecDuplicate(x,&b);
 27:   VecDuplicate(x,&u);

 29:   /* create a solution that is orthogonal to the constants */
 30:   VecGetOwnershipRange(u,&i1,&i2);
 31:   for (i=i1; i<i2; i++) {
 32:     avalue = i;
 33:     VecSetValues(u,1,&i,&avalue,INSERT_VALUES);
 34:   }
 35:   VecAssemblyBegin(u);
 36:   VecAssemblyEnd(u);
 37:   VecSum(u,&avalue);
 38:   avalue = -avalue/(PetscReal)n;
 39:   VecShift(&avalue,u);

 41:   /* Create and assemble matrix */
 42:   MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,&A);
 43:   MatSetFromOptions(A);
 44:   value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
 45:   for (i=1; i<n-1; i++) {
 46:     col[0] = i-1; col[1] = i; col[2] = i+1;
 47:     MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
 48:   }
 49:   i = n - 1; col[0] = n - 2; col[1] = n - 1; value[1] = 1.0;
 50:   MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
 51:   i = 0; col[0] = 0; col[1] = 1; value[0] = 1.0; value[1] = -1.0;
 52:   MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
 53:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 54:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 55:   MatMult(A,u,b);

 57:   /* Create KSP context; set operators and options; solve linear system */
 58:   KSPCreate(PETSC_COMM_WORLD,&ksp);
 59:   KSPSetOperators(ksp,A,A,DIFFERENT_NONZERO_PATTERN);

 61:   /* Insure that preconditioner has same null space as matrix */
 62:   /* currently does not do anything */
 63:   KSPGetPC(ksp,&pc);

 65:   KSPSetFromOptions(ksp);
 66:   KSPSolve(ksp,b,x);
 67:   /* KSPView(ksp,PETSC_VIEWER_STDOUT_WORLD); */

 69:   /* Check error */
 70:   VecAXPY(&none,u,x);
 71:   VecNorm(x,NORM_2,&norm);
 72:   KSPGetIterationNumber(ksp,&its);
 73:   PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A, Iterations %D\n",norm,its);

 75:   /* Free work space */
 76:   VecDestroy(x);VecDestroy(u);
 77:   VecDestroy(b);MatDestroy(A);
 78:   KSPDestroy(ksp);
 79:   PetscFinalize();
 80:   return 0;
 81: }