Actual source code: ex10.c
2: static char help[] = "Reads a PETSc matrix and vector from a file and solves a linear system.\n\
3: This version first preloads and solves a small system, then loads \n\
4: another (larger) system and solves it as well. This example illustrates\n\
5: preloading of instructions with the smaller system so that more accurate\n\
6: performance monitoring can be done with the larger one (that actually\n\
7: is the system of interest). See the 'Performance Hints' chapter of the\n\
8: users manual for a discussion of preloading. Input parameters include\n\
9: -f0 <input_file> : first file to load (small system)\n\
10: -f1 <input_file> : second file to load (larger system)\n\n\
11: -trans : solve transpose system instead\n\n";
12: /*
13: This code can be used to test PETSc interface to other packages.\n\
14: Examples of command line options: \n\
15: ./ex10 -f0 <datafile> -ksp_type preonly \n\
16: -help -ksp_view \n\
17: -num_numfac <num_numfac> -num_rhs <num_rhs> \n\
18: -ksp_type preonly -pc_type lu -pc_factor_mat_solver_package spooles or superlu or superlu_dist or mumps \n\
19: -ksp_type preonly -pc_type cholesky -pc_factor_mat_solver_package spooles or dscpack or mumps \n\
20: mpiexec -n <np> ./ex10 -f0 <datafile> -ksp_type cg -pc_type asm -pc_asm_type basic -sub_pc_type icc -mat_type sbaij
21: \n\n";
22: */
23: /*T
24: Concepts: KSP^solving a linear system
25: Processors: n
26: T*/
28: /*
29: Include "petscksp.h" so that we can use KSP solvers. Note that this file
30: automatically includes:
31: petscsys.h - base PETSc routines petscvec.h - vectors
32: petscmat.h - matrices
33: petscis.h - index sets petscksp.h - Krylov subspace methods
34: petscviewer.h - viewers petscpc.h - preconditioners
35: */
36: #include petscksp.h
40: int main(int argc,char **args)
41: {
42: KSP ksp; /* linear solver context */
43: Mat A,B = 0; /* matrix */
44: Vec x,b,u; /* approx solution, RHS, exact solution */
45: PetscViewer fd; /* viewer */
46: char file[4][PETSC_MAX_PATH_LEN]; /* input file name */
47: PetscTruth table=PETSC_FALSE,flg,flgB=PETSC_FALSE,trans=PETSC_FALSE,initialguess = PETSC_FALSE;
48: PetscTruth outputSoln=PETSC_FALSE;
50: PetscInt its,num_numfac,m,n,M;
51: PetscReal norm;
52: PetscLogDouble tsetup,tsetup1,tsetup2,tsolve,tsolve1,tsolve2;
53: PetscTruth preload=PETSC_TRUE,isSymmetric,cknorm=PETSC_FALSE,initialguessfile = PETSC_FALSE;
54: PetscMPIInt rank;
55: char initialguessfilename[PETSC_MAX_PATH_LEN];
57: PetscInitialize(&argc,&args,(char *)0,help);
58: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
59: PetscOptionsGetTruth(PETSC_NULL,"-table",&table,PETSC_NULL);
60: PetscOptionsGetTruth(PETSC_NULL,"-trans",&trans,PETSC_NULL);
61: PetscOptionsGetTruth(PETSC_NULL,"-initialguess",&initialguess,PETSC_NULL);
62: PetscOptionsGetTruth(PETSC_NULL,"-output_solution",&outputSoln,PETSC_NULL);
63: PetscOptionsGetString(PETSC_NULL,"-initialguessfilename",initialguessfilename,PETSC_MAX_PATH_LEN-1,&initialguessfile);
65: /*
66: Determine files from which we read the two linear systems
67: (matrix and right-hand-side vector).
68: */
69: PetscOptionsGetString(PETSC_NULL,"-f",file[0],PETSC_MAX_PATH_LEN-1,&flg);
70: if (flg) {
71: PetscStrcpy(file[1],file[0]);
72: preload = PETSC_FALSE;
73: } else {
74: PetscOptionsGetString(PETSC_NULL,"-f0",file[0],PETSC_MAX_PATH_LEN-1,&flg);
75: if (!flg) SETERRQ(1,"Must indicate binary file with the -f0 or -f option");
76: PetscOptionsGetString(PETSC_NULL,"-f1",file[1],PETSC_MAX_PATH_LEN-1,&flg);
77: if (!flg) {preload = PETSC_FALSE;} /* don't bother with second system */
78: }
80: /* -----------------------------------------------------------
81: Beginning of linear solver loop
82: ----------------------------------------------------------- */
83: /*
84: Loop through the linear solve 2 times.
85: - The intention here is to preload and solve a small system;
86: then load another (larger) system and solve it as well.
87: This process preloads the instructions with the smaller
88: system so that more accurate performance monitoring (via
89: -log_summary) can be done with the larger one (that actually
90: is the system of interest).
91: */
92: PreLoadBegin(preload,"Load system");
94: /* - - - - - - - - - - - New Stage - - - - - - - - - - - - -
95: Load system
96: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
98: /*
99: Open binary file. Note that we use FILE_MODE_READ to indicate
100: reading from this file.
101: */
102: PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[PreLoadIt],FILE_MODE_READ,&fd);
103:
104: /*
105: Load the matrix and vector; then destroy the viewer.
106: */
107: MatLoad(fd,MATAIJ,&A);
108:
109: if (!preload){
110: flg = PETSC_FALSE;
111: PetscOptionsGetString(PETSC_NULL,"-rhs",file[2],PETSC_MAX_PATH_LEN-1,&flg);
112: if (flg){ /* rhs is stored in a separate file */
113: if (file[2][0] == '0') {
114: PetscInt m;
115: PetscScalar one = 1.0;
116: PetscInfo(0,"Using vector of ones for RHS\n");
117: MatGetLocalSize(A,&m,PETSC_NULL);
118: VecCreate(PETSC_COMM_WORLD,&b);
119: VecSetSizes(b,m,PETSC_DECIDE);
120: VecSetFromOptions(b);
121: VecSet(b,one);
122: } else {
123: PetscViewerDestroy(fd);
124: PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[2],FILE_MODE_READ,&fd);
125: VecLoad(fd,PETSC_NULL,&b);
126: }
127: } else {
128: VecLoad(fd,PETSC_NULL,&b);
129: }
130: }
131: PetscViewerDestroy(fd);
133: /* Make A singular for testing zero-pivot of ilu factorization */
134: /* Example: ./ex10 -f0 <datafile> -test_zeropivot -set_row_zero -pc_factor_shift_nonzero */
135: flg = PETSC_FALSE;
136: PetscOptionsGetTruth(PETSC_NULL, "-test_zeropivot", &flg,PETSC_NULL);
137: if (flg) {
138: PetscInt row,ncols;
139: const PetscInt *cols;
140: const PetscScalar *vals;
141: PetscTruth flg1=PETSC_FALSE;
142: PetscScalar *zeros;
143: row = 0;
144: MatGetRow(A,row,&ncols,&cols,&vals);
145: PetscMalloc(sizeof(PetscScalar)*(ncols+1),&zeros);
146: PetscMemzero(zeros,(ncols+1)*sizeof(PetscScalar));
147: PetscOptionsGetTruth(PETSC_NULL, "-set_row_zero", &flg1,PETSC_NULL);
148: if (flg1){ /* set entire row as zero */
149: MatSetValues(A,1,&row,ncols,cols,zeros,INSERT_VALUES);
150: } else { /* only set (row,row) entry as zero */
151: MatSetValues(A,1,&row,1,&row,zeros,INSERT_VALUES);
152: }
153: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
154: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
155: }
157: /* Check whether A is symmetric */
158: flg = PETSC_FALSE;
159: PetscOptionsGetTruth(PETSC_NULL, "-check_symmetry", &flg,PETSC_NULL);
160: if (flg) {
161: Mat Atrans;
162: MatTranspose(A, MAT_INITIAL_MATRIX,&Atrans);
163: MatEqual(A, Atrans, &isSymmetric);
164: if (isSymmetric) {
165: MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE);
166: } else {
167: PetscPrintf(PETSC_COMM_WORLD,"Warning: A is non-symmetric \n");
168: }
169: MatDestroy(Atrans);
170: }
172: /*
173: If the loaded matrix is larger than the vector (due to being padded
174: to match the block size of the system), then create a new padded vector.
175: */
176:
177: MatGetLocalSize(A,&m,&n);
178: if (m != n) {
179: SETERRQ2(PETSC_ERR_ARG_SIZ, "This example is not intended for rectangular matrices (%d, %d)", m, n);
180: }
181: MatGetSize(A,&M,PETSC_NULL);
182: VecGetSize(b,&m);
183: if (M != m) { /* Create a new vector b by padding the old one */
184: PetscInt j,mvec,start,end,indx;
185: Vec tmp;
186: PetscScalar *bold;
188: VecCreate(PETSC_COMM_WORLD,&tmp);
189: VecSetSizes(tmp,n,PETSC_DECIDE);
190: VecSetFromOptions(tmp);
191: VecGetOwnershipRange(b,&start,&end);
192: VecGetLocalSize(b,&mvec);
193: VecGetArray(b,&bold);
194: for (j=0; j<mvec; j++) {
195: indx = start+j;
196: VecSetValues(tmp,1,&indx,bold+j,INSERT_VALUES);
197: }
198: VecRestoreArray(b,&bold);
199: VecDestroy(b);
200: VecAssemblyBegin(tmp);
201: VecAssemblyEnd(tmp);
202: b = tmp;
203: }
204: VecDuplicate(b,&x);
205: VecDuplicate(b,&u);
206: if (initialguessfile) {
207: PetscViewer viewer2;
208: PetscViewerBinaryOpen(PETSC_COMM_WORLD,initialguessfilename,FILE_MODE_READ,&viewer2);
209: VecLoadIntoVector(viewer2,x);
210: PetscViewerDestroy(viewer2);
211: initialguess = PETSC_TRUE;
212: } else {
213: VecSet(x,0.0);
214: }
217: /* Check scaling in A */
218: flg = PETSC_FALSE;
219: PetscOptionsGetTruth(PETSC_NULL, "-check_scaling", &flg,PETSC_NULL);
220: if (flg) {
221: Vec max, min;
222: PetscInt idx;
223: PetscReal val;
225: VecDuplicate(x, &max);
226: VecDuplicate(x, &min);
227: MatGetRowMaxAbs(A, max, PETSC_NULL);
228: MatGetRowMinAbs(A, min, PETSC_NULL);
229: {
230: PetscViewer viewer;
232: PetscViewerASCIIOpen(PETSC_COMM_WORLD, "max.data", &viewer);
233: VecView(max, viewer);
234: PetscViewerDestroy(viewer);
235: PetscViewerASCIIOpen(PETSC_COMM_WORLD, "min.data", &viewer);
236: VecView(min, viewer);
237: PetscViewerDestroy(viewer);
238: }
239: VecView(max, PETSC_VIEWER_DRAW_WORLD);
240: VecMax(max, &idx, &val);
241: PetscPrintf(PETSC_COMM_WORLD, "Largest max row element %G at row %d\n", val, idx);
242: VecView(min, PETSC_VIEWER_DRAW_WORLD);
243: VecMin(min, &idx, &val);
244: PetscPrintf(PETSC_COMM_WORLD, "Smallest min row element %G at row %d\n", val, idx);
245: VecMin(max, &idx, &val);
246: PetscPrintf(PETSC_COMM_WORLD, "Smallest max row element %G at row %d\n", val, idx);
247: VecPointwiseDivide(max, max, min);
248: VecMax(max, &idx, &val);
249: PetscPrintf(PETSC_COMM_WORLD, "Largest row ratio %G at row %d\n", val, idx);
250: VecView(max, PETSC_VIEWER_DRAW_WORLD);
251: VecDestroy(max);
252: VecDestroy(min);
253: }
255: /* - - - - - - - - - - - New Stage - - - - - - - - - - - - -
256: Setup solve for system
257: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
258: /*
259: Conclude profiling last stage; begin profiling next stage.
260: */
261: PreLoadStage("KSPSetUpSolve");
263: /*
264: We also explicitly time this stage via PetscGetTime()
265: */
266: PetscGetTime(&tsetup1);
268: /*
269: Create linear solver; set operators; set runtime options.
270: */
271: KSPCreate(PETSC_COMM_WORLD,&ksp);
272: KSPSetInitialGuessNonzero(ksp,initialguess);
273: num_numfac = 1;
274: PetscOptionsGetInt(PETSC_NULL,"-num_numfac",&num_numfac,PETSC_NULL);
275: while ( num_numfac-- ){
276: PetscTruth lsqr;
277: char str[32];
278: PetscOptionsGetString(PETSC_NULL,"-ksp_type",str,32,&lsqr);
279: if (lsqr) {
280: PetscStrcmp("lsqr",str,&lsqr);
281: }
282: if (lsqr) {
283: Mat B;
284: MatMatMultTranspose(A,A,MAT_INITIAL_MATRIX,4,&B);
285: KSPSetOperators(ksp,A,B,SAME_NONZERO_PATTERN);
286: MatDestroy(B);
287: } else {
288: KSPSetOperators(ksp,A,A,SAME_NONZERO_PATTERN);
289: }
290: KSPSetFromOptions(ksp);
292: /*
293: Here we explicitly call KSPSetUp() and KSPSetUpOnBlocks() to
294: enable more precise profiling of setting up the preconditioner.
295: These calls are optional, since both will be called within
296: KSPSolve() if they haven't been called already.
297: */
298: KSPSetUp(ksp);
299: KSPSetUpOnBlocks(ksp);
300: PetscGetTime(&tsetup2);
301: tsetup = tsetup2 - tsetup1;
303: /* - - - - - - - - - - - New Stage - - - - - - - - - - - - -
304: Solve system
305: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
307: /*
308: Solve linear system; we also explicitly time this stage.
309: */
310: PetscGetTime(&tsolve1);
311: if (trans) {
312: KSPSolveTranspose(ksp,b,x);
313: KSPGetIterationNumber(ksp,&its);
314: } else {
315: PetscInt num_rhs=1;
316: PetscOptionsGetInt(PETSC_NULL,"-num_rhs",&num_rhs,PETSC_NULL);
317: cknorm = PETSC_FALSE;
318: PetscOptionsGetTruth(PETSC_NULL,"-cknorm",&cknorm,PETSC_NULL);
319: while ( num_rhs-- ) {
320: if (num_rhs == 1) VecSet(x,0.0);
321: KSPSolve(ksp,b,x);
322: }
323: KSPGetIterationNumber(ksp,&its);
324: if (cknorm){ /* Check error for each rhs */
325: if (trans) {
326: MatMultTranspose(A,x,u);
327: } else {
328: MatMult(A,x,u);
329: }
330: VecAXPY(u,-1.0,b);
331: VecNorm(u,NORM_2,&norm);
332: PetscPrintf(PETSC_COMM_WORLD," Number of iterations = %3D\n",its);
333: PetscPrintf(PETSC_COMM_WORLD," Residual norm %A\n",norm);
334: }
335: } /* while ( num_rhs-- ) */
336: PetscGetTime(&tsolve2);
337: tsolve = tsolve2 - tsolve1;
339: /* - - - - - - - - - - - New Stage - - - - - - - - - - - - -
340: Check error, print output, free data structures.
341: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
343: /*
344: Check error
345: */
346: if (trans) {
347: MatMultTranspose(A,x,u);
348: } else {
349: MatMult(A,x,u);
350: }
351: VecAXPY(u,-1.0,b);
352: VecNorm(u,NORM_2,&norm);
354: /*
355: Write output (optinally using table for solver details).
356: - PetscPrintf() handles output for multiprocessor jobs
357: by printing from only one processor in the communicator.
358: - KSPView() prints information about the linear solver.
359: */
360: if (table) {
361: char *matrixname,kspinfo[120];
362: PetscViewer viewer;
364: /*
365: Open a string viewer; then write info to it.
366: */
367: PetscViewerStringOpen(PETSC_COMM_WORLD,kspinfo,120,&viewer);
368: KSPView(ksp,viewer);
369: PetscStrrchr(file[PreLoadIt],'/',&matrixname);
370: PetscPrintf(PETSC_COMM_WORLD,"%-8.8s %3D %2.0e %2.1e %2.1e %2.1e %s \n",
371: matrixname,its,norm,tsetup+tsolve,tsetup,tsolve,kspinfo);
373: /*
374: Destroy the viewer
375: */
376: PetscViewerDestroy(viewer);
377: } else {
378: PetscPrintf(PETSC_COMM_WORLD,"Number of iterations = %3D\n",its);
379: PetscPrintf(PETSC_COMM_WORLD,"Residual norm %A\n",norm);
380: }
381: PetscOptionsGetString(PETSC_NULL,"-solution",file[3],PETSC_MAX_PATH_LEN-1,&flg);
382: if (flg) {
383: PetscViewer viewer;
384: Vec xstar;
385: PetscReal norm;
387: PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[3],FILE_MODE_READ,&viewer);
388: VecLoad(viewer, VECMPI, &xstar);
389: VecAXPY(xstar, -1.0, x);
390: VecNorm(xstar, NORM_2, &norm);
391: PetscPrintf(PETSC_COMM_WORLD, "Error norm %A\n", norm);
392: VecDestroy(xstar);
393: PetscViewerDestroy(viewer);
394: }
395: if (outputSoln) {
396: PetscViewer viewer;
398: PetscViewerBinaryOpen(PETSC_COMM_WORLD,"solution.petsc",FILE_MODE_WRITE,&viewer);
399: VecView(x, viewer);
400: PetscViewerDestroy(viewer);
401: }
403: flg = PETSC_FALSE;
404: PetscOptionsGetTruth(PETSC_NULL, "-ksp_reason", &flg,PETSC_NULL);
405: if (flg){
406: KSPConvergedReason reason;
407: KSPGetConvergedReason(ksp,&reason);
408: PetscPrintf(PETSC_COMM_WORLD,"KSPConvergedReason: %D\n", reason);
409: }
410:
411: } /* while ( num_numfac-- ) */
413: /*
414: Free work space. All PETSc objects should be destroyed when they
415: are no longer needed.
416: */
417: MatDestroy(A); VecDestroy(b);
418: VecDestroy(u); VecDestroy(x);
419: KSPDestroy(ksp);
420: if (flgB) { MatDestroy(B); }
421: PreLoadEnd();
422: /* -----------------------------------------------------------
423: End of linear solver loop
424: ----------------------------------------------------------- */
426: PetscFinalize();
427: return 0;
428: }