Actual source code: ex53.c

petsc-3.3-p7 2013-05-11
  2: /* Program usage:  ./ex53 [-help] [all PETSc options] */

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

  8: #include <petscksp.h>

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

 24:   PetscInitialize(&argc,&args,(char *)0,help);
 25:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 26:   PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);

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

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

 44:   if (!rank){
 45:     value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
 46:     for (i=1; i<n-1; i++) {
 47:       col[0] = i-1; col[1] = i; col[2] = i+1;
 48:       MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
 49:     }
 50:     i = n - 1; col[0] = n - 2; col[1] = n - 1;
 51:     MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
 52:     i = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0;
 53:     MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);

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

 61:   /* Set exact solution */
 62:   VecSet(u,one);

 64:   /* Create linear solver context */
 65:   KSPCreate(PETSC_COMM_WORLD,&ksp);
 66:   KSPSetOperators(ksp,A,A,SAME_PRECONDITIONER);
 67:   KSPGetPC(ksp,&pc);
 68:   PCSetType(pc,PCLU);
 69: #ifdef PETSC_HAVE_MUMPS 
 70:   PCFactorSetMatSolverPackage(pc,MATSOLVERMUMPS);
 71: #endif
 72:   KSPSetFromOptions(ksp);
 73: 
 74:   /* 1. Solve linear system A x = b */
 75:   MatMult(A,u,b);
 76:   KSPSolve(ksp,b,x);

 78:   /* Check the error */
 79:   VecAXPY(x,neg_one,u);
 80:   VecNorm(x,NORM_2,&norm);
 81:   KSPGetIterationNumber(ksp,&its);
 82:   if (norm > tol){
 83:     PetscPrintf(PETSC_COMM_WORLD,"1. Norm of error for Ax=b: %G, Iterations %D\n",
 84:                      norm,its);
 85:   }
 86: 
 87:   /* 2. Solve linear system A^T x = b*/
 88:   MatMultTranspose(A,u,b);
 89:   KSPSolveTranspose(ksp,b,x2);
 90: 
 91:   /* Check the error */
 92:   VecAXPY(x2,neg_one,u);
 93:   VecNorm(x2,NORM_2,&norm);
 94:   KSPGetIterationNumber(ksp,&its);
 95:   if (norm > tol){
 96:     PetscPrintf(PETSC_COMM_WORLD,"2. Norm of error for A^T x=b: %G, Iterations %D\n",
 97:                      norm,its);
 98:   }

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

108:   MatMult(A,u,b);
109:   KSPSolve(ksp,b,x);

111:   /* Check the error */
112:   VecAXPY(x,neg_one,u);
113:   VecNorm(x,NORM_2,&norm);
114:   KSPGetIterationNumber(ksp,&its);
115:   if (norm > tol){
116:     PetscPrintf(PETSC_COMM_WORLD,"3. Norm of error for (A+Delta) x=b: %G, Iterations %D\n",
117:                      norm,its);
118:   }

120:   /* Free work space. */
121:   VecDestroy(&x);
122:   VecDestroy(&u);
123:   VecDestroy(&x2);
124:   VecDestroy(&b);
125:   MatDestroy(&A);
126:   KSPDestroy(&ksp);

128:   PetscFinalize();
129:   return 0;
130: }