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,Ii,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 (Ii=Istart; Ii<Iend; Ii++) {
36: v = -1.0; i = Ii/n; j = Ii - i*n;
37: if (i>0) {J = Ii - n; MatSetValues(C,1,&Ii,1,&J,&v,ADD_VALUES);}
38: if (i<m-1) {J = Ii + n; MatSetValues(C,1,&Ii,1,&J,&v,ADD_VALUES);}
39: if (j>0) {J = Ii - 1; MatSetValues(C,1,&Ii,1,&J,&v,ADD_VALUES);}
40: if (j<n-1) {J = Ii + 1; MatSetValues(C,1,&Ii,1,&J,&v,ADD_VALUES);}
41: v = 4.0; MatSetValues(C,1,&Ii,1,&Ii,&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,&r);
63: PetscRandomSetFromOptions(r);
64: VecSetRandom(u,r);
65: PetscRandomDestroy(r);
66: CHKMEMQ;
67: MatMult(C,u,b);
68: CHKMEMQ;
70: for (k=0; k<3; k++){
71: CHKMEMQ;
72: if (k == 0){ /* CG */
73: KSPCreate(PETSC_COMM_WORLD,&ksp);
74: KSPSetOperators(ksp,C,C,DIFFERENT_NONZERO_PATTERN);
75: PetscPrintf(PETSC_COMM_WORLD,"\n CG: \n");
76: KSPSetType(ksp,KSPCG);
77: } else if (k == 1){ /* MINRES */
78: KSPCreate(PETSC_COMM_WORLD,&ksp);
79: KSPSetOperators(ksp,C,C,DIFFERENT_NONZERO_PATTERN);
80: PetscPrintf(PETSC_COMM_WORLD,"\n MINRES: \n");
81: KSPSetType(ksp,KSPMINRES);
82: } else { /* SYMMLQ */
83: KSPCreate(PETSC_COMM_WORLD,&ksp);
84: KSPSetOperators(ksp,C,C,DIFFERENT_NONZERO_PATTERN);
85: PetscPrintf(PETSC_COMM_WORLD,"\n SYMMLQ: \n");
86: KSPSetType(ksp,KSPSYMMLQ);
87: }
88: CHKMEMQ;
89: KSPGetPC(ksp,&pc);
90: /* PCSetType(pc,PCICC); */
91: PCSetType(pc,PCJACOBI);
92: KSPSetTolerances(ksp,1.e-7,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
94: /*
95: Set runtime options, e.g.,
96: -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
97: These options will override those specified above as long as
98: KSPSetFromOptions() is called _after_ any other customization
99: routines.
100: */
101: CHKMEMQ;
102: KSPSetFromOptions(ksp);
104: /* Solve linear system; */
105: CHKMEMQ;
106: KSPSetUp(ksp);
107: CHKMEMQ;
108: KSPSolve(ksp,b,x);
110: KSPGetIterationNumber(ksp,&its);
111: /* Check error */
112: VecCopy(u,u_tmp);
113: VecAXPY(u_tmp,none,x);
114: VecNorm(u_tmp,NORM_2,&err_norm);
115: MatMult(C,x,u_tmp);
116: VecAXPY(u_tmp,none,b);
117: VecNorm(u_tmp,NORM_2,&res_norm);
118:
119: PetscPrintf(PETSC_COMM_WORLD,"Number of iterations = %3D\n",its);
120: PetscPrintf(PETSC_COMM_WORLD,"Residual norm %A;",res_norm);
121: PetscPrintf(PETSC_COMM_WORLD," Error norm %A.\n",err_norm);
123: KSPDestroy(ksp);
124: }
125:
126: /*
127: Free work space. All PETSc objects should be destroyed when they
128: are no longer needed.
129: */
130: VecDestroy(b);
131: VecDestroy(u);
132: VecDestroy(x);
133: VecDestroy(u_tmp);
134: MatDestroy(C);
136: PetscFinalize();
137: return 0;
138: }