Actual source code: ex15.c

petsc-3.12.1 2019-10-22
Report Typos and Errors

  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 */
 12:   PetscInt       i,n = 10,col[3],its,i1,i2;
 13:   PetscScalar    none = -1.0,value[3],avalue;
 14:   PetscReal      norm;
 15:   PC             pc;

 17:   PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
 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; value[1] = 2.0; value[2] = -1.0;
 44:   for (i=1; i<n-1; i++) {
 45:     col[0] = i-1; col[1] = i; col[2] = i+1;
 46:     MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
 47:   }
 48:   i    = n - 1; col[0] = n - 2; col[1] = n - 1; value[1] = 1.0;
 49:   MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
 50:   i    = 0; col[0] = 0; col[1] = 1; value[0] = 1.0; value[1] = -1.0;
 51:   MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
 52:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 53:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 54:   MatMult(A,u,b);

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

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

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

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

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