Actual source code: ex24.c

  1: /*$Id: ex24.c,v 1.13 2001/08/07 21:30:50 bsmith Exp $*/

  3: static char help[] = "Tests CG, MINRES and SYMMLQ on symmetric matrices with SBAIJ format. The preconditioner ICC only works on sequential SBAIJ format. nn";

 5:  #include petscsles.h


  8: int main(int argc,char **args)
  9: {
 10:   Mat         C;
 11:   PetscScalar v,none = -1.0;
 12:   int         i,j,I,J,ierr,Istart,Iend,N,m = 4,n = 4,rank,size,its,k;
 13:   PetscReal   err_norm,res_norm;
 14:   Vec         x,b,u,u_tmp;
 15:   PetscRandom r;
 16:   SLES        sles;
 17:   PC          pc;
 18:   KSP         ksp;

 20:   PetscInitialize(&argc,&args,(char *)0,help);
 21:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 22:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 23:   PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
 24:   PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
 25:   N = m*n;


 28:   /* Generate matrix */
 29:   MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,N,N,&C);
 30:   MatSetFromOptions(C);
 31:   MatGetOwnershipRange(C,&Istart,&Iend);
 32:   for (I=Istart; I<Iend; I++) {
 33:     v = -1.0; i = I/n; j = I - i*n;
 34:     if (i>0)   {J = I - n; MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES);}
 35:     if (i<m-1) {J = I + n; MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES);}
 36:     if (j>0)   {J = I - 1; MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES);}
 37:     if (j<n-1) {J = I + 1; MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES);}
 38:     v = 4.0; MatSetValues(C,1,&I,1,&I,&v,ADD_VALUES);
 39:   }
 40:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 41:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

 43:   /* a shift can make C indefinite. Preconditioners LU, ILU (for BAIJ format) and ICC may fail */
 44:   /* MatShift(&alpha, C); */
 45:   /* MatView(C,PETSC_VIEWER_STDOUT_WORLD); */

 47:   /* Setup and solve for system */
 48: 
 49:   /* Create vectors.  */
 50:   VecCreate(PETSC_COMM_WORLD,&x);
 51:   VecSetSizes(x,PETSC_DECIDE,N);
 52:   VecSetFromOptions(x);
 53:   VecDuplicate(x,&b);
 54:   VecDuplicate(x,&u);
 55:   VecDuplicate(x,&u_tmp);

 57:   /* Set exact solution u; then compute right-hand-side vector b. */
 58:   PetscRandomCreate(PETSC_COMM_SELF,RANDOM_DEFAULT,&r);
 59:   VecSetRandom(r,u);
 60:   PetscRandomDestroy(r);
 61: 
 62:   MatMult(C,u,b);

 64:   for (k=0; k<3; k++){
 65:     if (k == 0){                              /* CG  */
 66:       SLESCreate(PETSC_COMM_WORLD,&sles);
 67:       SLESSetOperators(sles,C,C,DIFFERENT_NONZERO_PATTERN);
 68:       SLESGetKSP(sles,&ksp);
 69:       PetscPrintf(PETSC_COMM_WORLD,"n CG: n");
 70:       KSPSetType(ksp,KSPCG);
 71:     } else if (k == 1){                       /* MINRES */
 72:       SLESCreate(PETSC_COMM_WORLD,&sles);
 73:       SLESSetOperators(sles,C,C,DIFFERENT_NONZERO_PATTERN);
 74:       SLESGetKSP(sles,&ksp);
 75:       PetscPrintf(PETSC_COMM_WORLD,"n MINRES: n");
 76:       KSPSetType(ksp,KSPMINRES);
 77:     } else {                                 /* SYMMLQ */
 78:       SLESCreate(PETSC_COMM_WORLD,&sles);
 79:       SLESSetOperators(sles,C,C,DIFFERENT_NONZERO_PATTERN);
 80:       SLESGetKSP(sles,&ksp);
 81:       PetscPrintf(PETSC_COMM_WORLD,"n SYMMLQ: n");
 82:       KSPSetType(ksp,KSPSYMMLQ);
 83:     }

 85:     SLESGetPC(sles,&pc);
 86:     /* PCSetType(pc,PCICC); */
 87:     PCSetType(pc,PCJACOBI);
 88:     KSPSetTolerances(ksp,1.e-7,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);

 90:     /*
 91:     Set runtime options, e.g.,
 92:         -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
 93:     These options will override those specified above as long as
 94:     SLESSetFromOptions() is called _after_ any other customization
 95:     routines.
 96:     */
 97:     SLESSetFromOptions(sles);

 99:     /* Solve linear system; */
100:     SLESSolve(sles,b,x,&its);
101: 
102:   /* Check error */
103:     VecCopy(u,u_tmp);
104:     VecAXPY(&none,x,u_tmp);
105:     VecNorm(u_tmp,NORM_2,&err_norm);
106:     MatMult(C,x,u_tmp);
107:     VecAXPY(&none,b,u_tmp);
108:     VecNorm(u_tmp,NORM_2,&res_norm);
109: 
110:     PetscPrintf(PETSC_COMM_WORLD,"Number of iterations = %3dn",its);
111:     PetscPrintf(PETSC_COMM_WORLD,"Residual norm %A;",res_norm);
112:     PetscPrintf(PETSC_COMM_WORLD,"  Error norm %A.n",err_norm);

114:     SLESDestroy(sles);
115:   }
116: 
117:   /* 
118:        Free work space.  All PETSc objects should be destroyed when they
119:        are no longer needed.
120:   */
121:   VecDestroy(b);
122:   VecDestroy(u);
123:   VecDestroy(x);
124:   VecDestroy(u_tmp);
125:   MatDestroy(C);

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