Actual source code: ex6.c
petsc-3.7.2 2016-06-05
2: static char help[] = "Reads a PETSc matrix and vector from a file and solves a linear system.\n\
3: Input arguments are:\n\
4: -f <input_file> : file to load. For example see $PETSC_DIR/share/petsc/datafiles/matrices\n\n";
6: #include <petscksp.h>
7: #include <petsclog.h>
11: int main(int argc,char **args)
12: {
13: #if !defined(PETSC_USE_COMPLEX)
15: PetscInt its;
16: PetscLogStage stage1,stage2;
17: PetscReal norm;
18: Vec x,b,u;
19: Mat A;
20: char file[PETSC_MAX_PATH_LEN];
21: PetscViewer fd;
22: PetscBool table = PETSC_FALSE,flg;
23: KSP ksp;
24: #endif
26: PetscInitialize(&argc,&args,(char*)0,help);
28: #if defined(PETSC_USE_COMPLEX)
29: SETERRQ(PETSC_COMM_WORLD,1,"This example does not work with complex numbers");
30: #else
31: PetscOptionsGetBool(NULL,NULL,"-table",&table,NULL);
34: /* Read matrix and RHS */
35: PetscOptionsGetString(NULL,NULL,"-f",file,PETSC_MAX_PATH_LEN,&flg);
36: if (!flg) SETERRQ(PETSC_COMM_WORLD,1,"Must indicate binary file with the -f option");
37: PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd);
38: MatCreate(PETSC_COMM_WORLD,&A);
39: MatLoad(A,fd);
40: VecCreate(PETSC_COMM_WORLD,&b);
41: VecLoad(b,fd);
42: PetscViewerDestroy(&fd);
44: /*
45: If the load matrix is larger then the vector, due to being padded
46: to match the blocksize then create a new padded vector
47: */
48: {
49: PetscInt m,n,j,mvec,start,end,indx;
50: Vec tmp;
51: PetscScalar *bold;
53: MatGetLocalSize(A,&m,&n);
54: VecCreate(PETSC_COMM_WORLD,&tmp);
55: VecSetSizes(tmp,m,PETSC_DECIDE);
56: VecSetFromOptions(tmp);
57: VecGetOwnershipRange(b,&start,&end);
58: VecGetLocalSize(b,&mvec);
59: VecGetArray(b,&bold);
60: for (j=0; j<mvec; j++) {
61: indx = start+j;
62: VecSetValues(tmp,1,&indx,bold+j,INSERT_VALUES);
63: }
64: VecRestoreArray(b,&bold);
65: VecDestroy(&b);
66: VecAssemblyBegin(tmp);
67: VecAssemblyEnd(tmp);
68: b = tmp;
69: }
70: VecDuplicate(b,&x);
71: VecDuplicate(b,&u);
73: VecSet(x,0.0);
74: PetscBarrier((PetscObject)A);
76: PetscLogStageRegister("mystage 1",&stage1);
77: PetscLogStagePush(stage1);
78: KSPCreate(PETSC_COMM_WORLD,&ksp);
79: KSPSetOperators(ksp,A,A);
80: KSPSetFromOptions(ksp);
81: KSPSetUp(ksp);
82: KSPSetUpOnBlocks(ksp);
83: PetscLogStagePop();
84: PetscBarrier((PetscObject)A);
86: PetscLogStageRegister("mystage 2",&stage2);
87: PetscLogStagePush(stage2);
88: KSPSolve(ksp,b,x);
89: PetscLogStagePop();
91: /* Show result */
92: MatMult(A,x,u);
93: VecAXPY(u,-1.0,b);
94: VecNorm(u,NORM_2,&norm);
95: KSPGetIterationNumber(ksp,&its);
96: /* matrix PC KSP Options its residual */
97: if (table) {
98: char *matrixname,kspinfo[120];
99: PetscViewer viewer;
100: PetscViewerStringOpen(PETSC_COMM_WORLD,kspinfo,120,&viewer);
101: KSPView(ksp,viewer);
102: PetscStrrchr(file,'/',&matrixname);
103: PetscPrintf(PETSC_COMM_WORLD,"%-8.8s %3D %2.0e %s \n",matrixname,its,norm,kspinfo);
104: PetscViewerDestroy(&viewer);
105: } else {
106: PetscPrintf(PETSC_COMM_WORLD,"Number of iterations = %3D\n",its);
107: PetscPrintf(PETSC_COMM_WORLD,"Residual norm = %g\n",(double)norm);
108: }
110: /* Cleanup */
111: KSPDestroy(&ksp);
112: VecDestroy(&x);
113: VecDestroy(&b);
114: VecDestroy(&u);
115: MatDestroy(&A);
117: PetscFinalize();
118: #endif
119: return 0;
120: }