Actual source code: ex24.c
petsc-3.7.5 2017-01-01
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: 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(NULL,NULL,"-m",&m,NULL);
25: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
26: N = m*n;
29: /* Generate matrix */
30: MatCreate(PETSC_COMM_WORLD,&C);
31: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N);
32: MatSetFromOptions(C);
33: MatSetUp(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: 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: /* Create vectors. */
52: VecCreate(PETSC_COMM_WORLD,&x);
53: VecSetSizes(x,PETSC_DECIDE,N);
54: VecSetFromOptions(x);
55: VecDuplicate(x,&b);
56: VecDuplicate(x,&u);
57: VecDuplicate(x,&u_tmp);
58: /* Set exact solution u; then compute right-hand-side vector b. */
59: VecSet(u,1.0);
60: MatMult(C,u,b);
62: for (k=0; k<3; k++) {
63: if (k == 0) { /* CG */
64: KSPCreate(PETSC_COMM_WORLD,&ksp);
65: KSPSetOperators(ksp,C,C);
66: PetscPrintf(PETSC_COMM_WORLD,"\n CG: \n");
67: KSPSetType(ksp,KSPCG);
68: } else if (k == 1) { /* MINRES */
69: KSPCreate(PETSC_COMM_WORLD,&ksp);
70: KSPSetOperators(ksp,C,C);
71: PetscPrintf(PETSC_COMM_WORLD,"\n MINRES: \n");
72: KSPSetType(ksp,KSPMINRES);
73: } else { /* SYMMLQ */
74: KSPCreate(PETSC_COMM_WORLD,&ksp);
75: KSPSetOperators(ksp,C,C);
76: PetscPrintf(PETSC_COMM_WORLD,"\n SYMMLQ: \n");
77: KSPSetType(ksp,KSPSYMMLQ);
78: }
79: KSPGetPC(ksp,&pc);
80: /* PCSetType(pc,PCICC); */
81: PCSetType(pc,PCJACOBI);
82: KSPSetTolerances(ksp,1.e-7,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
84: /*
85: Set runtime options, e.g.,
86: -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
87: These options will override those specified above as long as
88: KSPSetFromOptions() is called _after_ any other customization
89: routines.
90: */
91: KSPSetFromOptions(ksp);
93: /* Solve linear system; */
94: KSPSetUp(ksp);
95: KSPSolve(ksp,b,x);
97: KSPGetIterationNumber(ksp,&its);
98: /* Check error */
99: VecCopy(u,u_tmp);
100: VecAXPY(u_tmp,none,x);
101: VecNorm(u_tmp,NORM_2,&err_norm);
102: MatMult(C,x,u_tmp);
103: VecAXPY(u_tmp,none,b);
104: VecNorm(u_tmp,NORM_2,&res_norm);
106: PetscPrintf(PETSC_COMM_WORLD,"Number of iterations = %3D\n",its);
107: PetscPrintf(PETSC_COMM_WORLD,"Residual norm %g;",(double)res_norm);
108: PetscPrintf(PETSC_COMM_WORLD," Error norm %g.\n",(double)err_norm);
109: KSPDestroy(&ksp);
110: }
112: /*
113: Free work space. All PETSc objects should be destroyed when they
114: are no longer needed.
115: */
116: VecDestroy(&b);
117: VecDestroy(&u);
118: VecDestroy(&x);
119: VecDestroy(&u_tmp);
120: MatDestroy(&C);
122: PetscFinalize();
123: return 0;
124: }