Actual source code: ex30.c
2: static char help[] = "Reads a PETSc matrix and vector from a file and solves a linear system.\n\
3: It is copied and intended to move dirty codes from ksp/examples/tutorials/ex10.c and simplify ex10.c.\n\
4: Input parameters include\n\
5: -f0 <input_file> : first file to load (small system)\n\
6: -f1 <input_file> : second file to load (larger system)\n\n\
7: -trans : solve transpose system instead\n\n";
8: /*
9: This code can be used to test PETSc interface to other packages.\n\
10: Examples of command line options: \n\
11: ex30 -f0 <datafile> -ksp_type preonly \n\
12: -help -ksp_view \n\
13: -num_numfac <num_numfac> -num_rhs <num_rhs> \n\
14: -ksp_type preonly -pc_type lu -pc_factor_mat_solver_package spooles or superlu or superlu_dist or mumps \n\
15: -ksp_type preonly -pc_type cholesky -pc_factor_mat_solver_package spooles or mumps \n\
16: mpiexec -n <np> ex30 -f0 <datafile> -ksp_type cg -pc_type asm -pc_asm_type basic -sub_pc_type icc -mat_type sbaij
18: ./ex30 -f0 $D/small -mat_sigma -3.999999999999999 -ksp_type fgmres -pc_type lu -pc_factor_mat_solver_package superlu -mat_superlu_conditionnumber -ckerror -mat_superlu_diagpivotthresh 0
19: ./ex30 -f0 $D/small -mat_sigma -3.999999999999999 -ksp_type fgmres -pc_type hypre -pc_hypre_type boomeramg -ksp_type fgmres -ckError
20: ./ex30 -f0 $D/small -mat_sigma -3.999999999999999 -ksp_type fgmres -pc_type lu -pc_factor_mat_solver_package petsc -pc_factor_shift_type NONZERO -pc_factor_shift_amount 1.e-5 -ckerror
21: \n\n";
22: */
23: /*T
24: Concepts: KSP solving a linear system
25: Processors: n
26: T*/
28: #include <petscksp.h>
32: int main(int argc,char **args)
33: {
34: KSP ksp;
35: Mat A,B;
36: Vec x,b,u,b2; /* approx solution, RHS, exact solution */
37: PetscViewer fd; /* viewer */
38: char file[4][PETSC_MAX_PATH_LEN]; /* input file name */
39: PetscBool table = PETSC_FALSE,flg,flgB=PETSC_FALSE,trans=PETSC_FALSE,partition=PETSC_FALSE,initialguess = PETSC_FALSE;
40: PetscBool outputSoln=PETSC_FALSE;
42: PetscInt its,num_numfac,n,M;
43: PetscReal rnorm,enorm;
44: PetscLogDouble tsetup,tsetup1,tsetup2,tsolve,tsolve1,tsolve2;
45: PetscBool preload=PETSC_TRUE,diagonalscale,isSymmetric,ckrnorm=PETSC_TRUE,Test_MatDuplicate=PETSC_FALSE,ckerror=PETSC_FALSE;
46: PetscMPIInt rank;
47: PetscScalar sigma;
48: PetscInt m;
50: PetscInitialize(&argc,&args,(char *)0,help);
51: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
52: PetscOptionsGetBool(PETSC_NULL,"-table",&table,PETSC_NULL);
53: PetscOptionsGetBool(PETSC_NULL,"-trans",&trans,PETSC_NULL);
54: PetscOptionsGetBool(PETSC_NULL,"-partition",&partition,PETSC_NULL);
55: PetscOptionsGetBool(PETSC_NULL,"-initialguess",&initialguess,PETSC_NULL);
56: PetscOptionsGetBool(PETSC_NULL,"-output_solution",&outputSoln,PETSC_NULL);
57: PetscOptionsGetBool(PETSC_NULL,"-ckrnorm",&ckrnorm,PETSC_NULL);
58: PetscOptionsGetBool(PETSC_NULL,"-ckerror",&ckerror,PETSC_NULL);
60: /*
61: Determine files from which we read the two linear systems
62: (matrix and right-hand-side vector).
63: */
64: PetscOptionsGetString(PETSC_NULL,"-f",file[0],PETSC_MAX_PATH_LEN,&flg);
65: if (flg) {
66: PetscStrcpy(file[1],file[0]);
67: preload = PETSC_FALSE;
68: } else {
69: PetscOptionsGetString(PETSC_NULL,"-f0",file[0],PETSC_MAX_PATH_LEN,&flg);
70: if (!flg) SETERRQ(PETSC_COMM_WORLD,1,"Must indicate binary file with the -f0 or -f option");
71: PetscOptionsGetString(PETSC_NULL,"-f1",file[1],PETSC_MAX_PATH_LEN,&flg);
72: if (!flg) {preload = PETSC_FALSE;} /* don't bother with second system */
73: }
75: /* -----------------------------------------------------------
76: Beginning of linear solver loop
77: ----------------------------------------------------------- */
78: /*
79: Loop through the linear solve 2 times.
80: - The intention here is to preload and solve a small system;
81: then load another (larger) system and solve it as well.
82: This process preloads the instructions with the smaller
83: system so that more accurate performance monitoring (via
84: -log_summary) can be done with the larger one (that actually
85: is the system of interest).
86: */
87: PetscPreLoadBegin(preload,"Load system");
89: /* - - - - - - - - - - - New Stage - - - - - - - - - - - - -
90: Load system
91: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
93: /*
94: Open binary file. Note that we use FILE_MODE_READ to indicate
95: reading from this file.
96: */
97: PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[PetscPreLoadIt],FILE_MODE_READ,&fd);
98:
99: /*
100: Load the matrix and vector; then destroy the viewer.
101: */
102: MatCreate(PETSC_COMM_WORLD,&A);
103: MatLoad(A,fd);
104:
105: if (!preload){
106: flg = PETSC_FALSE;
107: PetscOptionsGetString(PETSC_NULL,"-rhs",file[2],PETSC_MAX_PATH_LEN,&flg);
108: if (flg){ /* rhs is stored in a separate file */
109: PetscViewerDestroy(&fd);
110: PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[2],FILE_MODE_READ,&fd);
111: } else {
112: /* if file contains no RHS, then use a vector of all ones */
113: PetscInfo(0,"Using vector of ones for RHS\n");
114: MatGetLocalSize(A,&m,PETSC_NULL);
115: VecCreate(PETSC_COMM_WORLD,&b);
116: VecSetSizes(b,m,PETSC_DECIDE);
117: VecSetFromOptions(b);
118: VecSet(b,1.0);
119: PetscObjectSetName((PetscObject)b, "Rhs vector");
120: }
121: }
122: PetscViewerDestroy(&fd);
124: /* Test MatDuplicate() */
125: if (Test_MatDuplicate){
126: MatDuplicate(A,MAT_COPY_VALUES,&B);
127: MatEqual(A,B,&flg);
128: if (!flg){
129: PetscPrintf(PETSC_COMM_WORLD," A != B \n");
130: }
131: MatDestroy(&B);
132: }
134: /* Add a shift to A */
135: PetscOptionsGetScalar(PETSC_NULL,"-mat_sigma",&sigma,&flg);
136: if (flg) {
137: PetscOptionsGetString(PETSC_NULL,"-fB",file[2],PETSC_MAX_PATH_LEN,&flgB);
138: if (flgB){
139: /* load B to get A = A + sigma*B */
140: PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[2],FILE_MODE_READ,&fd);
141: MatCreate(PETSC_COMM_WORLD,&B);
142: MatSetOptionsPrefix(B,"B_"); /* e.g., ./ex30 -f0 <A> -fB <B> -mat_sigma 1.0 -B_mat_view_draw */
143: MatLoad(B,fd);
144: PetscViewerDestroy(&fd);
145: MatAXPY(A,sigma,B,DIFFERENT_NONZERO_PATTERN); /* A <- sigma*B + A */
146: } else {
147: MatShift(A,sigma);
148: }
149: }
151: /* Make A singular for testing zero-pivot of ilu factorization */
152: /* Example: ./ex30 -f0 <datafile> -test_zeropivot -set_row_zero -pc_factor_shift_nonzero */
153: flg = PETSC_FALSE;
154: PetscOptionsGetBool(PETSC_NULL, "-test_zeropivot", &flg,PETSC_NULL);
155: if (flg) {
156: PetscInt row,ncols;
157: const PetscInt *cols;
158: const PetscScalar *vals;
159: PetscBool flg1=PETSC_FALSE;
160: PetscScalar *zeros;
161: row = 0;
162: MatGetRow(A,row,&ncols,&cols,&vals);
163: PetscMalloc(sizeof(PetscScalar)*(ncols+1),&zeros);
164: PetscMemzero(zeros,(ncols+1)*sizeof(PetscScalar));
165: flg1 = PETSC_FALSE;
166: PetscOptionsGetBool(PETSC_NULL, "-set_row_zero", &flg1,PETSC_NULL);
167: if (flg1){ /* set entire row as zero */
168: MatSetValues(A,1,&row,ncols,cols,zeros,INSERT_VALUES);
169: } else { /* only set (row,row) entry as zero */
170: MatSetValues(A,1,&row,1,&row,zeros,INSERT_VALUES);
171: }
172: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
173: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
174: }
176: /* Check whether A is symmetric */
177: flg = PETSC_FALSE;
178: PetscOptionsGetBool(PETSC_NULL, "-check_symmetry", &flg,PETSC_NULL);
179: if (flg) {
180: Mat Atrans;
181: MatTranspose(A, MAT_INITIAL_MATRIX,&Atrans);
182: MatEqual(A, Atrans, &isSymmetric);
183: if (isSymmetric) {
184: PetscPrintf(PETSC_COMM_WORLD,"A is symmetric \n");
185: } else {
186: PetscPrintf(PETSC_COMM_WORLD,"A is non-symmetric \n");
187: }
188: MatDestroy(&Atrans);
189: }
191: /*
192: If the loaded matrix is larger than the vector (due to being padded
193: to match the block size of the system), then create a new padded vector.
194: */
195:
196: MatGetLocalSize(A,&m,&n);
197: if (m != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ, "This example is not intended for rectangular matrices (%d, %d)", m, n);
198: MatGetSize(A,&M,PETSC_NULL);
199: VecGetSize(b,&m);
200: if (M != m) { /* Create a new vector b by padding the old one */
201: PetscInt j,mvec,start,end,indx;
202: Vec tmp;
203: PetscScalar *bold;
205: VecCreate(PETSC_COMM_WORLD,&tmp);
206: VecSetSizes(tmp,n,PETSC_DECIDE);
207: VecSetFromOptions(tmp);
208: VecGetOwnershipRange(b,&start,&end);
209: VecGetLocalSize(b,&mvec);
210: VecGetArray(b,&bold);
211: for (j=0; j<mvec; j++) {
212: indx = start+j;
213: VecSetValues(tmp,1,&indx,bold+j,INSERT_VALUES);
214: }
215: VecRestoreArray(b,&bold);
216: VecDestroy(&b);
217: VecAssemblyBegin(tmp);
218: VecAssemblyEnd(tmp);
219: b = tmp;
220: }
221: VecDuplicate(b,&b2);
222: VecDuplicate(b,&x);
223: PetscObjectSetName((PetscObject)x, "Solution vector");
224: VecDuplicate(b,&u);
225: PetscObjectSetName((PetscObject)u, "True Solution vector");
226: VecSet(x,0.0);
228: if (ckerror){ /* Set true solution */
229: VecSet(u,1.0);
230: MatMult(A,u,b);
231: }
233: /* - - - - - - - - - - - New Stage - - - - - - - - - - - - -
234: Setup solve for system
235: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
237: if (partition) {
238: MatPartitioning mpart;
239: IS mis,nis,is;
240: PetscInt *count;
241: PetscMPIInt size;
242: Mat BB;
243: MPI_Comm_size(PETSC_COMM_WORLD,&size);
244: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
245: PetscMalloc(size*sizeof(PetscInt),&count);
246: MatPartitioningCreate(PETSC_COMM_WORLD, &mpart);
247: MatPartitioningSetAdjacency(mpart, A);
248: /* MatPartitioningSetVertexWeights(mpart, weight); */
249: MatPartitioningSetFromOptions(mpart);
250: MatPartitioningApply(mpart, &mis);
251: MatPartitioningDestroy(&mpart);
252: ISPartitioningToNumbering(mis,&nis);
253: ISPartitioningCount(mis,size,count);
254: ISDestroy(&mis);
255: ISInvertPermutation(nis, count[rank], &is);
256: PetscFree(count);
257: ISDestroy(&nis);
258: ISSort(is);
259: MatGetSubMatrix(A,is,is,MAT_INITIAL_MATRIX,&BB);
261: /* need to move the vector also */
262: ISDestroy(&is);
263: MatDestroy(&A);
264: A = BB;
265: }
267: /*
268: We also explicitly time this stage via PetscGetTime()
269: */
270: PetscGetTime(&tsetup1);
272: /*
273: Create linear solver; set operators; set runtime options.
274: */
275: KSPCreate(PETSC_COMM_WORLD,&ksp);
276: KSPSetInitialGuessNonzero(ksp,initialguess);
277: num_numfac = 1;
278: PetscOptionsGetInt(PETSC_NULL,"-num_numfac",&num_numfac,PETSC_NULL);
279: while ( num_numfac-- ){
280:
282: KSPSetOperators(ksp,A,A,SAME_NONZERO_PATTERN);
283: KSPSetFromOptions(ksp);
285: /*
286: Here we explicitly call KSPSetUp() and KSPSetUpOnBlocks() to
287: enable more precise profiling of setting up the preconditioner.
288: These calls are optional, since both will be called within
289: KSPSolve() if they haven't been called already.
290: */
291: KSPSetUp(ksp);
292: KSPSetUpOnBlocks(ksp);
293: PetscGetTime(&tsetup2);
294: tsetup = tsetup2 - tsetup1;
296: /*
297: Tests "diagonal-scaling of preconditioned residual norm" as used
298: by many ODE integrator codes including SUNDIALS. Note this is different
299: than diagonally scaling the matrix before computing the preconditioner
300: */
301: diagonalscale = PETSC_FALSE;
302: PetscOptionsGetBool(PETSC_NULL,"-diagonal_scale",&diagonalscale,PETSC_NULL);
303: if (diagonalscale) {
304: PC pc;
305: PetscInt j,start,end,n;
306: Vec scale;
307:
308: KSPGetPC(ksp,&pc);
309: VecGetSize(x,&n);
310: VecDuplicate(x,&scale);
311: VecGetOwnershipRange(scale,&start,&end);
312: for (j=start; j<end; j++) {
313: VecSetValue(scale,j,((PetscReal)(j+1))/((PetscReal)n),INSERT_VALUES);
314: }
315: VecAssemblyBegin(scale);
316: VecAssemblyEnd(scale);
317: PCSetDiagonalScale(pc,scale);
318: VecDestroy(&scale);
319: }
321: /* - - - - - - - - - - - New Stage - - - - - - - - - - - - -
322: Solve system
323: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
324: /*
325: Solve linear system; we also explicitly time this stage.
326: */
327: PetscGetTime(&tsolve1);
328: if (trans) {
329: KSPSolveTranspose(ksp,b,x);
330: KSPGetIterationNumber(ksp,&its);
331: } else {
332: PetscInt num_rhs=1;
333: PetscOptionsGetInt(PETSC_NULL,"-num_rhs",&num_rhs,PETSC_NULL);
334:
335: while ( num_rhs-- ) {
336: KSPSolve(ksp,b,x);
337: }
338: KSPGetIterationNumber(ksp,&its);
339: if (ckrnorm){ /* Check residual for each rhs */
340: if (trans) {
341: MatMultTranspose(A,x,b2);
342: } else {
343: MatMult(A,x,b2);
344: }
345: VecAXPY(b2,-1.0,b);
346: VecNorm(b2,NORM_2,&rnorm);
347: PetscPrintf(PETSC_COMM_WORLD," Number of iterations = %3D\n",its);
348: PetscPrintf(PETSC_COMM_WORLD," Residual norm %A\n",rnorm);
349: }
350: if (ckerror && !trans){ /* Check error for each rhs */
351: /* VecView(x,PETSC_VIEWER_STDOUT_WORLD); */
352: VecAXPY(u,-1.0,x);
353: VecNorm(u,NORM_2,&enorm);
354: PetscPrintf(PETSC_COMM_WORLD," Error norm %A\n",enorm);
355: }
356:
357: } /* while ( num_rhs-- ) */
358: PetscGetTime(&tsolve2);
359: tsolve = tsolve2 - tsolve1;
361:
362: /*
363: Write output (optinally using table for solver details).
364: - PetscPrintf() handles output for multiprocessor jobs
365: by printing from only one processor in the communicator.
366: - KSPView() prints information about the linear solver.
367: */
368: if (table && ckrnorm) {
369: char *matrixname,kspinfo[120];
370: PetscViewer viewer;
372: /*
373: Open a string viewer; then write info to it.
374: */
375: PetscViewerStringOpen(PETSC_COMM_WORLD,kspinfo,120,&viewer);
376: KSPView(ksp,viewer);
377: PetscStrrchr(file[PetscPreLoadIt],'/',&matrixname);
378: PetscPrintf(PETSC_COMM_WORLD,"%-8.8s %3D %2.0e %2.1e %2.1e %2.1e %s \n",
379: matrixname,its,rnorm,tsetup+tsolve,tsetup,tsolve,kspinfo);
381: /*
382: Destroy the viewer
383: */
384: PetscViewerDestroy(&viewer);
385: }
387: PetscOptionsGetString(PETSC_NULL,"-solution",file[3],PETSC_MAX_PATH_LEN,&flg);
388: if (flg) {
389: PetscViewer viewer;
390: Vec xstar;
392: PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[3],FILE_MODE_READ,&viewer);
393: VecCreate(PETSC_COMM_WORLD,&xstar);
394: VecLoad(xstar,viewer);
395: VecAXPY(xstar, -1.0, x);
396: VecNorm(xstar, NORM_2, &enorm);
397: PetscPrintf(PETSC_COMM_WORLD, "Error norm %A\n", enorm);
398: VecDestroy(&xstar);
399: PetscViewerDestroy(&viewer);
400: }
401: if (outputSoln) {
402: PetscViewer viewer;
404: PetscViewerBinaryOpen(PETSC_COMM_WORLD,"solution.petsc",FILE_MODE_WRITE,&viewer);
405: VecView(x, viewer);
406: PetscViewerDestroy(&viewer);
407: }
409: flg = PETSC_FALSE;
410: PetscOptionsGetBool(PETSC_NULL, "-ksp_reason", &flg,PETSC_NULL);
411: if (flg){
412: KSPConvergedReason reason;
413: KSPGetConvergedReason(ksp,&reason);
414: PetscPrintf(PETSC_COMM_WORLD,"KSPConvergedReason: %D\n", reason);
415: }
416:
417: } /* while ( num_numfac-- ) */
419: /*
420: Free work space. All PETSc objects should be destroyed when they
421: are no longer needed.
422: */
423: MatDestroy(&A); VecDestroy(&b);
424: VecDestroy(&u); VecDestroy(&x);
425: VecDestroy(&b2);
426: KSPDestroy(&ksp);
427: if (flgB) { MatDestroy(&B); }
428: PetscPreLoadEnd();
429: /* -----------------------------------------------------------
430: End of linear solver loop
431: ----------------------------------------------------------- */
433: PetscFinalize();
434: return 0;
435: }