Actual source code: ex25.c
petsc-3.7.3 2016-07-24
1: static char help[] = "Tests CG, MINRES and SYMMLQ on the symmetric indefinite matrices: afiro \n\n";
3: #include <petscksp.h>
7: int main(int argc,char **args)
8: {
9: Mat C;
10: PetscScalar none = -1.0;
11: PetscMPIInt rank,size;
13: PetscInt its,k;
14: PetscReal err_norm,res_norm;
15: Vec x,b,u,u_tmp;
16: PC pc;
17: KSP ksp;
18: PetscViewer view;
19: char filein[128]; /* input file name */
21: PetscInitialize(&argc,&args,(char*)0,help);
22: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
23: MPI_Comm_size(PETSC_COMM_WORLD,&size);
25: /* Load the binary data file "filein". Set runtime option: -f filein */
26: PetscPrintf(PETSC_COMM_WORLD,"\n Load dataset ...\n");
27: PetscOptionsGetString(NULL,NULL,"-f",filein,128,NULL);
28: PetscViewerBinaryOpen(PETSC_COMM_WORLD,filein,FILE_MODE_READ,&view);
29: MatCreate(PETSC_COMM_WORLD,&C);
30: MatSetType(C,MATMPISBAIJ);
31: MatLoad(C,view);
32: VecCreate(PETSC_COMM_WORLD,&b);
33: VecCreate(PETSC_COMM_WORLD,&u);
34: VecLoad(b,view);
35: VecLoad(u,view);
36: PetscViewerDestroy(&view);
37: /* VecView(b,VIEWER_STDOUT_WORLD); */
38: /* MatView(C,VIEWER_STDOUT_WORLD); */
40: VecDuplicate(u,&u_tmp);
42: /* Check accuracy of the data */
43: /*
44: MatMult(C,u,u_tmp);
45: VecAXPY(u_tmp,none,b);
46: VecNorm(u_tmp,NORM_2,&res_norm);
47: PetscPrintf(PETSC_COMM_WORLD,"Accuracy of the loading data: | b - A*u |_2 : %g \n",(double)res_norm);
48: */
50: /* Setup and solve for system */
51: VecDuplicate(b,&x);
52: for (k=0; k<3; k++) {
53: if (k == 0) { /* CG */
54: KSPCreate(PETSC_COMM_WORLD,&ksp);
55: KSPSetType(ksp,KSPCG);
56: KSPSetOperators(ksp,C,C);
57: PetscPrintf(PETSC_COMM_WORLD,"\n CG: \n");
58: } else if (k == 1) { /* MINRES */
59: KSPCreate(PETSC_COMM_WORLD,&ksp);
60: KSPSetType(ksp,KSPMINRES);
61: KSPSetOperators(ksp,C,C);
62: PetscPrintf(PETSC_COMM_WORLD,"\n MINRES: \n");
63: } else { /* SYMMLQ */
64: KSPCreate(PETSC_COMM_WORLD,&ksp);
65: KSPSetOperators(ksp,C,C);
66: KSPSetType(ksp,KSPSYMMLQ);
67: PetscPrintf(PETSC_COMM_WORLD,"\n SYMMLQ: \n");
68: }
70: KSPGetPC(ksp,&pc);
71: PCSetType(pc,PCNONE);
73: /*
74: Set runtime options, e.g.,
75: -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
76: -pc_type jacobi -pc_jacobi_type rowmax
77: These options will override those specified above as long as
78: KSPSetFromOptions() is called _after_ any other customization routines.
79: */
80: KSPSetFromOptions(ksp);
82: /* Solve linear system; */
83: KSPSolve(ksp,b,x);
84: KSPGetIterationNumber(ksp,&its);
86: /* Check error */
87: VecCopy(u,u_tmp);
88: VecAXPY(u_tmp,none,x);
89: VecNorm(u_tmp,NORM_2,&err_norm);
90: MatMult(C,x,u_tmp);
91: VecAXPY(u_tmp,none,b);
92: VecNorm(u_tmp,NORM_2,&res_norm);
94: PetscPrintf(PETSC_COMM_WORLD,"Number of iterations = %3d\n",its);
95: PetscPrintf(PETSC_COMM_WORLD,"Residual norm: %g;",(double)res_norm);
96: PetscPrintf(PETSC_COMM_WORLD," Error norm: %g.\n",(double)err_norm);
98: KSPDestroy(&ksp);
99: }
101: /*
102: Free work space. All PETSc objects should be destroyed when they
103: are no longer needed.
104: */
105: VecDestroy(&b);
106: VecDestroy(&u);
107: VecDestroy(&x);
108: VecDestroy(&u_tmp);
109: MatDestroy(&C);
111: PetscFinalize();
112: return 0;
113: }