Actual source code: ex24.c

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

 4:  #include petscksp.h


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

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


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

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

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

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

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

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

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

 98:     /* Solve linear system; */
 99:     KSPSolve(ksp,b,x);

101:     KSPGetIterationNumber(ksp,&its);
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 = %3D\n",its);
111:     PetscPrintf(PETSC_COMM_WORLD,"Residual norm %A;",res_norm);
112:     PetscPrintf(PETSC_COMM_WORLD,"  Error norm %A.\n",err_norm);

114:     KSPDestroy(ksp);
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: }