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;
12: PetscMPIInt size,rank;
13: PetscReal err_norm,res_norm;
14: Vec x,b,u,u_tmp;
15: PC pc;
16: KSP ksp;
18: PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
19: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
20: MPI_Comm_size(PETSC_COMM_WORLD,&size);
21: PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
22: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
23: N = m*n;
25: /* Generate matrix */
26: MatCreate(PETSC_COMM_WORLD,&C);
27: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N);
28: MatSetFromOptions(C);
29: MatSetUp(C);
30: MatGetOwnershipRange(C,&Istart,&Iend);
31: for (Ii=Istart; Ii<Iend; Ii++) {
32: v = -1.0; i = Ii/n; j = Ii - i*n;
33: if (i>0) {J = Ii - n; MatSetValues(C,1,&Ii,1,&J,&v,ADD_VALUES);}
34: if (i<m-1) {J = Ii + n; MatSetValues(C,1,&Ii,1,&J,&v,ADD_VALUES);}
35: if (j>0) {J = Ii - 1; MatSetValues(C,1,&Ii,1,&J,&v,ADD_VALUES);}
36: if (j<n-1) {J = Ii + 1; MatSetValues(C,1,&Ii,1,&J,&v,ADD_VALUES);}
37: v = 4.0; MatSetValues(C,1,&Ii,1,&Ii,&v,ADD_VALUES);
38: }
39: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
40: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
42: /* a shift can make C indefinite. Preconditioners LU, ILU (for BAIJ format) and ICC may fail */
43: /* MatShift(C,alpha); */
44: /* MatView(C,PETSC_VIEWER_STDOUT_WORLD); */
46: /* Setup and solve for system */
47: /* Create vectors. */
48: VecCreate(PETSC_COMM_WORLD,&x);
49: VecSetSizes(x,PETSC_DECIDE,N);
50: VecSetFromOptions(x);
51: VecDuplicate(x,&b);
52: VecDuplicate(x,&u);
53: VecDuplicate(x,&u_tmp);
54: /* Set exact solution u; then compute right-hand-side vector b. */
55: VecSet(u,1.0);
56: MatMult(C,u,b);
58: for (k=0; k<3; k++) {
59: if (k == 0) { /* CG */
60: KSPCreate(PETSC_COMM_WORLD,&ksp);
61: KSPSetOperators(ksp,C,C);
62: PetscPrintf(PETSC_COMM_WORLD,"\n CG: \n");
63: KSPSetType(ksp,KSPCG);
64: } else if (k == 1) { /* MINRES */
65: KSPCreate(PETSC_COMM_WORLD,&ksp);
66: KSPSetOperators(ksp,C,C);
67: PetscPrintf(PETSC_COMM_WORLD,"\n MINRES: \n");
68: KSPSetType(ksp,KSPMINRES);
69: } else { /* SYMMLQ */
70: KSPCreate(PETSC_COMM_WORLD,&ksp);
71: KSPSetOperators(ksp,C,C);
72: PetscPrintf(PETSC_COMM_WORLD,"\n SYMMLQ: \n");
73: KSPSetType(ksp,KSPSYMMLQ);
74: }
75: KSPGetPC(ksp,&pc);
76: /* PCSetType(pc,PCICC); */
77: PCSetType(pc,PCJACOBI);
78: KSPSetTolerances(ksp,1.e-7,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
80: /*
81: Set runtime options, e.g.,
82: -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
83: These options will override those specified above as long as
84: KSPSetFromOptions() is called _after_ any other customization
85: routines.
86: */
87: KSPSetFromOptions(ksp);
89: /* Solve linear system; */
90: KSPSetUp(ksp);
91: KSPSolve(ksp,b,x);
93: KSPGetIterationNumber(ksp,&its);
94: /* Check error */
95: VecCopy(u,u_tmp);
96: VecAXPY(u_tmp,none,x);
97: VecNorm(u_tmp,NORM_2,&err_norm);
98: MatMult(C,x,u_tmp);
99: VecAXPY(u_tmp,none,b);
100: VecNorm(u_tmp,NORM_2,&res_norm);
102: PetscPrintf(PETSC_COMM_WORLD,"Number of iterations = %3D\n",its);
103: PetscPrintf(PETSC_COMM_WORLD,"Residual norm %g;",(double)res_norm);
104: PetscPrintf(PETSC_COMM_WORLD," Error norm %g.\n",(double)err_norm);
105: KSPDestroy(&ksp);
106: }
108: /*
109: Free work space. All PETSc objects should be destroyed when they
110: are no longer needed.
111: */
112: VecDestroy(&b);
113: VecDestroy(&u);
114: VecDestroy(&x);
115: VecDestroy(&u_tmp);
116: MatDestroy(&C);
118: PetscFinalize();
119: return ierr;
120: }
122: /*TEST
124: test:
125: args: -pc_type icc -mat_type seqsbaij -mat_ignore_lower_triangular
127: test:
128: suffix: 2
129: args: -pc_type icc -pc_factor_levels 2 -mat_type seqsbaij -mat_ignore_lower_triangular
131: test:
132: suffix: 3
133: nsize: 2
134: args: -pc_type bjacobi -sub_pc_type icc -mat_type mpisbaij -mat_ignore_lower_triangular -ksp_max_it 8
136: test:
137: suffix: 4
138: nsize: 2
139: args: -pc_type bjacobi -sub_pc_type icc -sub_pc_factor_levels 1 -mat_type mpisbaij -mat_ignore_lower_triangular
141: TEST*/