Actual source code: ex53.c

petsc-3.8.4 2018-03-24
Report Typos and Errors

  2: static char help[] = "Solves a tridiagonal linear system with KSP. \n\
  3:                       Modified from ex1.c to illustrate reuse of preconditioner \n\
  4:                       Written as requested by [petsc-maint #63875] \n\n";

  6:  #include <petscksp.h>

  8: int main(int argc,char **args)
  9: {
 10:   Vec            x,x2,b,u;     /* approx solution, RHS, exact solution */
 11:   Mat            A;            /* linear system matrix */
 12:   KSP            ksp;          /* linear solver context */
 13:   PC             pc;           /* preconditioner context */
 14:   PetscReal      norm,tol=100.*PETSC_MACHINE_EPSILON; /* norm of solution error */
 16:   PetscInt       i,n = 10,col[3],its;
 17:   PetscMPIInt    rank,size;
 18:   PetscScalar    one = 1.0,value[3];

 20:   PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
 21:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 22:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 23:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);

 25:   /* Create vectors.*/
 26:   VecCreate(PETSC_COMM_WORLD,&x);
 27:   PetscObjectSetName((PetscObject) x, "Solution");
 28:   VecSetSizes(x,PETSC_DECIDE,n);
 29:   VecSetFromOptions(x);
 30:   VecDuplicate(x,&b);
 31:   VecDuplicate(x,&u);
 32:   VecDuplicate(x,&x2);

 34:   /* Create matrix. Only proc[0] sets values - not efficient for parallel processing!
 35:      See ex23.c for efficient parallel assembly matrix */
 36:   MatCreate(PETSC_COMM_WORLD,&A);
 37:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
 38:   MatSetFromOptions(A);
 39:   MatSetUp(A);

 41:   if (!rank) {
 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;
 48:     MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
 49:     i    = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0;
 50:     MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);

 52:     i    = 0; col[0] = n-1; value[0] = 0.5; /* make A non-symmetric */
 53:     MatSetValues(A,1,&i,1,col,value,INSERT_VALUES);
 54:   }
 55:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 56:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 58:   /* Set exact solution */
 59:   VecSet(u,one);

 61:   /* Create linear solver context */
 62:   KSPCreate(PETSC_COMM_WORLD,&ksp);
 63:   KSPSetOperators(ksp,A,A);
 64:   KSPGetPC(ksp,&pc);
 65:   PCSetType(pc,PCLU);
 66: #if defined(PETSC_HAVE_MUMPS)
 67:   if (size > 1) {
 68:     PCFactorSetMatSolverPackage(pc,MATSOLVERMUMPS);
 69:   }
 70: #endif
 71:   KSPSetFromOptions(ksp);

 73:   /* 1. Solve linear system A x = b */
 74:   MatMult(A,u,b);
 75:   KSPSolve(ksp,b,x);

 77:   /* Check the error */
 78:   VecAXPY(x,-1.0,u);
 79:   VecNorm(x,NORM_2,&norm);
 80:   KSPGetIterationNumber(ksp,&its);
 81:   if (norm > tol) {
 82:     PetscPrintf(PETSC_COMM_WORLD,"1. Norm of error for Ax=b: %g, Iterations %D\n",(double)norm,its);
 83:   }

 85:   /* 2. Solve linear system A^T x = b*/
 86:   MatMultTranspose(A,u,b);
 87:   KSPSolveTranspose(ksp,b,x2);

 89:   /* Check the error */
 90:   VecAXPY(x2,-1.0,u);
 91:   VecNorm(x2,NORM_2,&norm);
 92:   KSPGetIterationNumber(ksp,&its);
 93:   if (norm > tol) {
 94:     PetscPrintf(PETSC_COMM_WORLD,"2. Norm of error for A^T x=b: %g, Iterations %D\n",(double)norm,its);
 95:   }

 97:   /* 3. Change A and solve A x = b with an iterative solver using A=LU as a preconditioner*/
 98:   if (!rank) {
 99:     i    = 0; col[0] = n-1; value[0] = 1.e-2;
100:     MatSetValues(A,1,&i,1,col,value,ADD_VALUES);
101:   }
102:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
103:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

105:   MatMult(A,u,b);
106:   KSPSolve(ksp,b,x);

108:   /* Check the error */
109:   VecAXPY(x,-1.0,u);
110:   VecNorm(x,NORM_2,&norm);
111:   KSPGetIterationNumber(ksp,&its);
112:   if (norm > tol) {
113:     PetscPrintf(PETSC_COMM_WORLD,"3. Norm of error for (A+Delta) x=b: %g, Iterations %D\n",(double)norm,its);
114:   }

116:   /* Free work space. */
117:   VecDestroy(&x);
118:   VecDestroy(&u);
119:   VecDestroy(&x2);
120:   VecDestroy(&b);
121:   MatDestroy(&A);
122:   KSPDestroy(&ksp);

124:   PetscFinalize();
125:   return ierr;
126: }