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: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
44: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
46: /* a shift can make C indefinite. Preconditioners LU, ILU (for BAIJ format) and ICC may fail */
47: /* MatShift(C,alpha); */
48: /* MatView(C,PETSC_VIEWER_STDOUT_WORLD); */
50: /* Setup and solve for system */
51:
52: /* Create vectors. */
53: VecCreate(PETSC_COMM_WORLD,&x);
54: VecSetSizes(x,PETSC_DECIDE,N);
55: VecSetFromOptions(x);
56: VecDuplicate(x,&b);
57: VecDuplicate(x,&u);
58: VecDuplicate(x,&u_tmp);
60: /* Set exact solution u; then compute right-hand-side vector b. */
61: PetscRandomCreate(PETSC_COMM_SELF,RANDOM_DEFAULT,&r);
62: VecSetRandom(u,r);
63: PetscRandomDestroy(r);
64:
65: MatMult(C,u,b);
67: for (k=0; k<3; k++){
68: if (k == 0){ /* CG */
69: KSPCreate(PETSC_COMM_WORLD,&ksp);
70: KSPSetOperators(ksp,C,C,DIFFERENT_NONZERO_PATTERN);
71: PetscPrintf(PETSC_COMM_WORLD,"\n CG: \n");
72: KSPSetType(ksp,KSPCG);
73: } else if (k == 1){ /* MINRES */
74: KSPCreate(PETSC_COMM_WORLD,&ksp);
75: KSPSetOperators(ksp,C,C,DIFFERENT_NONZERO_PATTERN);
76: PetscPrintf(PETSC_COMM_WORLD,"\n MINRES: \n");
77: KSPSetType(ksp,KSPMINRES);
78: } else { /* SYMMLQ */
79: KSPCreate(PETSC_COMM_WORLD,&ksp);
80: KSPSetOperators(ksp,C,C,DIFFERENT_NONZERO_PATTERN);
81: PetscPrintf(PETSC_COMM_WORLD,"\n SYMMLQ: \n");
82: KSPSetType(ksp,KSPSYMMLQ);
83: }
85: KSPGetPC(ksp,&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: KSPSetFromOptions() is called _after_ any other customization
95: routines.
96: */
97: KSPSetFromOptions(ksp);
99: /* Solve linear system; */
100: KSPSolve(ksp,b,x);
102: KSPGetIterationNumber(ksp,&its);
103: /* Check error */
104: VecCopy(u,u_tmp);
105: VecAXPY(u_tmp,none,x);
106: VecNorm(u_tmp,NORM_2,&err_norm);
107: MatMult(C,x,u_tmp);
108: VecAXPY(u_tmp,none,b);
109: VecNorm(u_tmp,NORM_2,&res_norm);
110:
111: PetscPrintf(PETSC_COMM_WORLD,"Number of iterations = %3D\n",its);
112: PetscPrintf(PETSC_COMM_WORLD,"Residual norm %A;",res_norm);
113: PetscPrintf(PETSC_COMM_WORLD," Error norm %A.\n",err_norm);
115: KSPDestroy(ksp);
116: }
117:
118: /*
119: Free work space. All PETSc objects should be destroyed when they
120: are no longer needed.
121: */
122: VecDestroy(b);
123: VecDestroy(u);
124: VecDestroy(x);
125: VecDestroy(u_tmp);
126: MatDestroy(C);
128: PetscFinalize();
129: return 0;
130: }