Actual source code: ex24.c
1: /*$Id: ex24.c,v 1.13 2001/08/07 21:30:50 bsmith Exp $*/
3: 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";
5: #include petscksp.h
10: int main(int argc,char **args)
11: {
12: Mat C;
13: PetscScalar v,none = -1.0;
14: int i,j,I,J,ierr,Istart,Iend,N,m = 4,n = 4,rank,size,its,k;
15: PetscReal err_norm,res_norm;
16: Vec x,b,u,u_tmp;
17: PetscRandom r;
18: PC pc;
19: KSP ksp;
21: PetscInitialize(&argc,&args,(char *)0,help);
22: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
23: MPI_Comm_size(PETSC_COMM_WORLD,&size);
24: PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
25: PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
26: N = m*n;
29: /* Generate matrix */
30: MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,N,N,&C);
31: MatSetFromOptions(C);
32: MatGetOwnershipRange(C,&Istart,&Iend);
33: for (I=Istart; I<Iend; I++) {
34: v = -1.0; i = I/n; j = I - i*n;
35: if (i>0) {J = I - n; MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES);}
36: if (i<m-1) {J = I + n; MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES);}
37: if (j>0) {J = I - 1; MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES);}
38: if (j<n-1) {J = I + 1; MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES);}
39: v = 4.0; MatSetValues(C,1,&I,1,&I,&v,ADD_VALUES);
40: }
41: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
42: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
44: /* a shift can make C indefinite. Preconditioners LU, ILU (for BAIJ format) and ICC may fail */
45: /* MatShift(&alpha, C); */
46: /* MatView(C,PETSC_VIEWER_STDOUT_WORLD); */
48: /* Setup and solve for system */
49:
50: /* Create vectors. */
51: VecCreate(PETSC_COMM_WORLD,&x);
52: VecSetSizes(x,PETSC_DECIDE,N);
53: VecSetFromOptions(x);
54: VecDuplicate(x,&b);
55: VecDuplicate(x,&u);
56: VecDuplicate(x,&u_tmp);
58: /* Set exact solution u; then compute right-hand-side vector b. */
59: PetscRandomCreate(PETSC_COMM_SELF,RANDOM_DEFAULT,&r);
60: VecSetRandom(r,u);
61: PetscRandomDestroy(r);
62:
63: MatMult(C,u,b);
65: for (k=0; k<3; k++){
66: if (k == 0){ /* CG */
67: KSPCreate(PETSC_COMM_WORLD,&ksp);
68: KSPSetOperators(ksp,C,C,DIFFERENT_NONZERO_PATTERN);
69: PetscPrintf(PETSC_COMM_WORLD,"\n CG: \n");
70: KSPSetType(ksp,KSPCG);
71: } else if (k == 1){ /* MINRES */
72: KSPCreate(PETSC_COMM_WORLD,&ksp);
73: KSPSetOperators(ksp,C,C,DIFFERENT_NONZERO_PATTERN);
74: PetscPrintf(PETSC_COMM_WORLD,"\n MINRES: \n");
75: KSPSetType(ksp,KSPMINRES);
76: } else { /* SYMMLQ */
77: KSPCreate(PETSC_COMM_WORLD,&ksp);
78: KSPSetOperators(ksp,C,C,DIFFERENT_NONZERO_PATTERN);
79: PetscPrintf(PETSC_COMM_WORLD,"\n SYMMLQ: \n");
80: KSPSetType(ksp,KSPSYMMLQ);
81: }
83: KSPGetPC(ksp,&pc);
84: /* PCSetType(pc,PCICC); */
85: PCSetType(pc,PCJACOBI);
86: KSPSetTolerances(ksp,1.e-7,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
88: /*
89: Set runtime options, e.g.,
90: -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
91: These options will override those specified above as long as
92: KSPSetFromOptions() is called _after_ any other customization
93: routines.
94: */
95: KSPSetFromOptions(ksp);
97: /* Solve linear system; */
98: KSPSetRhs(ksp,b);
99: KSPSetSolution(ksp,x);
100: KSPSolve(ksp);
102: KSPGetIterationNumber(ksp,&its);
103: /* Check error */
104: VecCopy(u,u_tmp);
105: VecAXPY(&none,x,u_tmp);
106: VecNorm(u_tmp,NORM_2,&err_norm);
107: MatMult(C,x,u_tmp);
108: VecAXPY(&none,b,u_tmp);
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: }