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,&C);
 32:   MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N);
 33:   MatSetFromOptions(C);
 34:   MatGetOwnershipRange(C,&Istart,&Iend);
 35:   for (I=Istart; I<Iend; I++) {
 36:     v = -1.0; i = I/n; j = I - i*n;
 37:     if (i>0)   {J = I - n; MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES);}
 38:     if (i<m-1) {J = I + n; MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES);}
 39:     if (j>0)   {J = I - 1; MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES);}
 40:     if (j<n-1) {J = I + 1; MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES);}
 41:     v = 4.0; MatSetValues(C,1,&I,1,&I,&v,ADD_VALUES);
 42:   }
 43:     CHKMEMQ;
 44:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 45:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

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

 51:   /* Setup and solve for system */
 52:     CHKMEMQ;
 53:   /* Create vectors.  */
 54:   VecCreate(PETSC_COMM_WORLD,&x);
 55:   VecSetSizes(x,PETSC_DECIDE,N);
 56:   VecSetFromOptions(x);
 57:   VecDuplicate(x,&b);
 58:   VecDuplicate(x,&u);
 59:   VecDuplicate(x,&u_tmp);
 60:     CHKMEMQ;
 61:   /* Set exact solution u; then compute right-hand-side vector b. */
 62:   PetscRandomCreate(PETSC_COMM_SELF,RANDOM_DEFAULT,&r);
 63:   VecSetRandom(u,r);
 64:   PetscRandomDestroy(r);
 65:       CHKMEMQ;
 66:   MatMult(C,u,b);
 67:     CHKMEMQ;

 69:   for (k=0; k<3; k++){
 70:     CHKMEMQ;
 71:     if (k == 0){                              /* CG  */
 72:       KSPCreate(PETSC_COMM_WORLD,&ksp);
 73:       KSPSetOperators(ksp,C,C,DIFFERENT_NONZERO_PATTERN);
 74:       PetscPrintf(PETSC_COMM_WORLD,"\n CG: \n");
 75:       KSPSetType(ksp,KSPCG);
 76:     } else if (k == 1){                       /* MINRES */
 77:       KSPCreate(PETSC_COMM_WORLD,&ksp);
 78:       KSPSetOperators(ksp,C,C,DIFFERENT_NONZERO_PATTERN);
 79:       PetscPrintf(PETSC_COMM_WORLD,"\n MINRES: \n");
 80:       KSPSetType(ksp,KSPMINRES);
 81:     } else {                                 /* SYMMLQ */
 82:       KSPCreate(PETSC_COMM_WORLD,&ksp);
 83:       KSPSetOperators(ksp,C,C,DIFFERENT_NONZERO_PATTERN);
 84:       PetscPrintf(PETSC_COMM_WORLD,"\n SYMMLQ: \n");
 85:       KSPSetType(ksp,KSPSYMMLQ);
 86:     }
 87:     CHKMEMQ;
 88:     KSPGetPC(ksp,&pc);
 89:     /* PCSetType(pc,PCICC); */
 90:     PCSetType(pc,PCJACOBI);
 91:     KSPSetTolerances(ksp,1.e-7,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);

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

103:     /* Solve linear system; */
104:     CHKMEMQ;
105:     KSPSetUp(ksp);
106:     CHKMEMQ;
107:     KSPSolve(ksp,b,x);

109:     KSPGetIterationNumber(ksp,&its);
110:   /* Check error */
111:     VecCopy(u,u_tmp);
112:     VecAXPY(u_tmp,none,x);
113:     VecNorm(u_tmp,NORM_2,&err_norm);
114:     MatMult(C,x,u_tmp);
115:     VecAXPY(u_tmp,none,b);
116:     VecNorm(u_tmp,NORM_2,&res_norm);
117: 
118:     PetscPrintf(PETSC_COMM_WORLD,"Number of iterations = %3D\n",its);
119:     PetscPrintf(PETSC_COMM_WORLD,"Residual norm %A;",res_norm);
120:     PetscPrintf(PETSC_COMM_WORLD,"  Error norm %A.\n",err_norm);

122:     KSPDestroy(ksp);
123:   }
124: 
125:   /* 
126:        Free work space.  All PETSc objects should be destroyed when they
127:        are no longer needed.
128:   */
129:   VecDestroy(b);
130:   VecDestroy(u);
131:   VecDestroy(x);
132:   VecDestroy(u_tmp);
133:   MatDestroy(C);

135:   PetscFinalize();
136:   return 0;
137: }