Actual source code: ex6.c
1: /*$Id: ex6.c,v 1.74 2001/09/11 16:33:24 bsmith Exp $*/
3: static char help[] = "Reads a PETSc matrix and vector from a file and solves a linear system.n
4: Input arguments are:n
5: -f <input_file> : file to load. For a 5X5 example of the 5-pt. stencil,n
6: use the file petsc/src/mat/examples/matbinary.exnn";
8: #include petscsles.h
9: #include petsclog.h
11: int main(int argc,char **args)
12: {
13: #if !defined(PETSC_USE_COMPLEX)
14: int ierr,its,stage1,stage2;
15: PetscReal norm;
16: PetscLogDouble tsetup1,tsetup2,tsetup,tsolve1,tsolve2,tsolve;
17: PetscScalar zero = 0.0,none = -1.0;
18: Vec x,b,u;
19: Mat A;
20: SLES sles;
21: char file[128];
22: PetscViewer fd;
23: PetscTruth table,flg;
24: #endif
26: PetscInitialize(&argc,&args,(char *)0,help);
28: #if defined(PETSC_USE_COMPLEX)
29: SETERRQ(1,"This example does not work with complex numbers");
30: #else
31: PetscOptionsHasName(PETSC_NULL,"-table",&table);
34: /* Read matrix and RHS */
35: PetscOptionsGetString(PETSC_NULL,"-f",file,127,&flg);
36: if (!flg) SETERRQ(1,"Must indicate binary file with the -f option");
37: PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,PETSC_BINARY_RDONLY,&fd);
39: MatLoad(fd,MATSEQAIJ,&A);
40: VecLoad(fd,&b);
41: PetscViewerDestroy(fd);
43: /*
44: If the load matrix is larger then the vector, due to being padded
45: to match the blocksize then create a new padded vector
46: */
47: {
48: int m,n,j,mvec,start,end,index;
49: Vec tmp;
50: PetscScalar *bold;
52: MatGetLocalSize(A,&m,&n);
53: VecCreate(PETSC_COMM_WORLD,&tmp);
54: VecSetSizes(tmp,m,PETSC_DECIDE);
55: VecSetFromOptions(tmp);
56: VecGetOwnershipRange(b,&start,&end);
57: VecGetLocalSize(b,&mvec);
58: VecGetArray(b,&bold);
59: for (j=0; j<mvec; j++) {
60: index = start+j;
61: ierr = VecSetValues(tmp,1,&index,bold+j,INSERT_VALUES);
62: }
63: VecRestoreArray(b,&bold);
64: VecDestroy(b);
65: VecAssemblyBegin(tmp);
66: VecAssemblyEnd(tmp);
67: b = tmp;
68: }
69: VecDuplicate(b,&x);
70: VecDuplicate(b,&u);
72: VecSet(&zero,x);
73: PetscBarrier((PetscObject)A);
75: PetscLogStageRegister(&stage1,"mystage 1");
76: PetscLogStagePush(stage1);
77: PetscGetTime(&tsetup1);
78: SLESCreate(PETSC_COMM_WORLD,&sles);
79: SLESSetOperators(sles,A,A,DIFFERENT_NONZERO_PATTERN);
80: SLESSetFromOptions(sles);
81: SLESSetUp(sles,b,x);
82: SLESSetUpOnBlocks(sles);
83: PetscGetTime(&tsetup2);
84: tsetup = tsetup2 -tsetup1;
85: PetscLogStagePop();
86: PetscBarrier((PetscObject)A);
88: PetscLogStageRegister(&stage2,"mystage 2");
89: PetscLogStagePush(stage2);
90: PetscGetTime(&tsolve1);
91: SLESSolve(sles,b,x,&its);
92: PetscGetTime(&tsolve2);
93: tsolve = tsolve2 - tsolve1;
94: PetscLogStagePop();
96: /* Show result */
97: MatMult(A,x,u);
98: VecAXPY(&none,b,u);
99: VecNorm(u,NORM_2,&norm);
100: /* matrix PC KSP Options its residual setuptime solvetime */
101: if (table) {
102: char *matrixname,slesinfo[120];
103: PetscViewer viewer;
104: PetscViewerStringOpen(PETSC_COMM_WORLD,slesinfo,120,&viewer);
105: SLESView(sles,viewer);
106: PetscStrrchr(file,'/',&matrixname);
107: PetscPrintf(PETSC_COMM_WORLD,"%-8.8s %3d %2.0e %2.1e %2.1e %2.1e %s n",
108: matrixname,its,norm,tsetup+tsolve,tsetup,tsolve,slesinfo);
109: PetscViewerDestroy(viewer);
110: } else {
111: PetscPrintf(PETSC_COMM_WORLD,"Number of iterations = %3dn",its);
112: PetscPrintf(PETSC_COMM_WORLD,"Residual norm = %An",norm);
113: }
115: /* Cleanup */
116: SLESDestroy(sles);
117: VecDestroy(x);
118: VecDestroy(b);
119: VecDestroy(u);
120: MatDestroy(A);
122: PetscFinalize();
123: #endif
124: return 0;
125: }