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>

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

 17:   PetscInitialize(&argc,&args,(char*)0,help);
 18:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 19:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 20:   PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
 21:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 22:   N    = m*n;

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

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

 45:   /* Setup and solve for system */
 46:   /* Create vectors.  */
 47:   VecCreate(PETSC_COMM_WORLD,&x);
 48:   VecSetSizes(x,PETSC_DECIDE,N);
 49:   VecSetFromOptions(x);
 50:   VecDuplicate(x,&b);
 51:   VecDuplicate(x,&u);
 52:   VecDuplicate(x,&u_tmp);
 53:   /* Set exact solution u; then compute right-hand-side vector b. */
 54:   VecSet(u,1.0);
 55:   MatMult(C,u,b);

 57:   for (k=0; k<3; k++) {
 58:     if (k == 0) {                              /* CG  */
 59:       KSPCreate(PETSC_COMM_WORLD,&ksp);
 60:       KSPSetOperators(ksp,C,C);
 61:       PetscPrintf(PETSC_COMM_WORLD,"\n CG: \n");
 62:       KSPSetType(ksp,KSPCG);
 63:     } else if (k == 1) {                       /* MINRES */
 64:       KSPCreate(PETSC_COMM_WORLD,&ksp);
 65:       KSPSetOperators(ksp,C,C);
 66:       PetscPrintf(PETSC_COMM_WORLD,"\n MINRES: \n");
 67:       KSPSetType(ksp,KSPMINRES);
 68:     } else {                                 /* SYMMLQ */
 69:       KSPCreate(PETSC_COMM_WORLD,&ksp);
 70:       KSPSetOperators(ksp,C,C);
 71:       PetscPrintf(PETSC_COMM_WORLD,"\n SYMMLQ: \n");
 72:       KSPSetType(ksp,KSPSYMMLQ);
 73:     }
 74:     KSPGetPC(ksp,&pc);
 75:     /* PCSetType(pc,PCICC); */
 76:     PCSetType(pc,PCJACOBI);
 77:     KSPSetTolerances(ksp,1.e-7,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);

 79:     /*
 80:     Set runtime options, e.g.,
 81:         -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
 82:     These options will override those specified above as long as
 83:     KSPSetFromOptions() is called _after_ any other customization
 84:     routines.
 85:     */
 86:     KSPSetFromOptions(ksp);

 88:     /* Solve linear system; */
 89:     KSPSetUp(ksp);
 90:     KSPSolve(ksp,b,x);

 92:     KSPGetIterationNumber(ksp,&its);
 93:     /* Check error */
 94:     VecCopy(u,u_tmp);
 95:     VecAXPY(u_tmp,none,x);
 96:     VecNorm(u_tmp,NORM_2,&err_norm);
 97:     MatMult(C,x,u_tmp);
 98:     VecAXPY(u_tmp,none,b);
 99:     VecNorm(u_tmp,NORM_2,&res_norm);

101:     PetscPrintf(PETSC_COMM_WORLD,"Number of iterations = %3D\n",its);
102:     PetscPrintf(PETSC_COMM_WORLD,"Residual norm %g;",(double)res_norm);
103:     PetscPrintf(PETSC_COMM_WORLD,"  Error norm %g.\n",(double)err_norm);
104:     KSPDestroy(&ksp);
105:   }

107:   /*
108:        Free work space.  All PETSc objects should be destroyed when they
109:        are no longer needed.
110:   */
111:   VecDestroy(&b);
112:   VecDestroy(&u);
113:   VecDestroy(&x);
114:   VecDestroy(&u_tmp);
115:   MatDestroy(&C);

117:   PetscFinalize();
118:   return 0;
119: }

121: /*TEST

123:     test:
124:       args: -pc_type icc -mat_type seqsbaij -mat_ignore_lower_triangular

126:     test:
127:       suffix: 2
128:       args: -pc_type icc -pc_factor_levels 2  -mat_type seqsbaij -mat_ignore_lower_triangular

130:     test:
131:       suffix: 3
132:       nsize: 2
133:       args: -pc_type bjacobi -sub_pc_type icc  -mat_type mpisbaij -mat_ignore_lower_triangular -ksp_max_it 8

135:     test:
136:       suffix: 4
137:       nsize: 2
138:       args: -pc_type bjacobi -sub_pc_type icc -sub_pc_factor_levels 1 -mat_type mpisbaij -mat_ignore_lower_triangular

140: TEST*/