Actual source code: ex30.c
petsc-3.5.4 2015-05-23
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 superlu or superlu_dist or mumps \n\
15: -ksp_type preonly -pc_type cholesky -pc_factor_mat_solver_package 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: PetscBool preload=PETSC_TRUE,diagonalscale,isSymmetric,ckrnorm=PETSC_TRUE,Test_MatDuplicate=PETSC_FALSE,ckerror=PETSC_FALSE;
45: PetscMPIInt rank;
46: PetscScalar sigma;
47: PetscInt m;
49: PetscInitialize(&argc,&args,(char*)0,help);
50: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
51: PetscOptionsGetBool(NULL,"-table",&table,NULL);
52: PetscOptionsGetBool(NULL,"-trans",&trans,NULL);
53: PetscOptionsGetBool(NULL,"-partition",&partition,NULL);
54: PetscOptionsGetBool(NULL,"-initialguess",&initialguess,NULL);
55: PetscOptionsGetBool(NULL,"-output_solution",&outputSoln,NULL);
56: PetscOptionsGetBool(NULL,"-ckrnorm",&ckrnorm,NULL);
57: PetscOptionsGetBool(NULL,"-ckerror",&ckerror,NULL);
59: /*
60: Determine files from which we read the two linear systems
61: (matrix and right-hand-side vector).
62: */
63: PetscOptionsGetString(NULL,"-f",file[0],PETSC_MAX_PATH_LEN,&flg);
64: if (flg) {
65: PetscStrcpy(file[1],file[0]);
66: preload = PETSC_FALSE;
67: } else {
68: PetscOptionsGetString(NULL,"-f0",file[0],PETSC_MAX_PATH_LEN,&flg);
69: if (!flg) SETERRQ(PETSC_COMM_WORLD,1,"Must indicate binary file with the -f0 or -f option");
70: PetscOptionsGetString(NULL,"-f1",file[1],PETSC_MAX_PATH_LEN,&flg);
71: if (!flg) preload = PETSC_FALSE; /* don't bother with second system */
72: }
74: /* -----------------------------------------------------------
75: Beginning of linear solver loop
76: ----------------------------------------------------------- */
77: /*
78: Loop through the linear solve 2 times.
79: - The intention here is to preload and solve a small system;
80: then load another (larger) system and solve it as well.
81: This process preloads the instructions with the smaller
82: system so that more accurate performance monitoring (via
83: -log_summary) can be done with the larger one (that actually
84: is the system of interest).
85: */
86: PetscPreLoadBegin(preload,"Load system");
88: /* - - - - - - - - - - - New Stage - - - - - - - - - - - - -
89: Load system
90: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
92: /*
93: Open binary file. Note that we use FILE_MODE_READ to indicate
94: reading from this file.
95: */
96: PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[PetscPreLoadIt],FILE_MODE_READ,&fd);
98: /*
99: Load the matrix and vector; then destroy the viewer.
100: */
101: MatCreate(PETSC_COMM_WORLD,&A);
102: MatSetFromOptions(A);
103: MatLoad(A,fd);
105: if (!preload) {
106: flg = PETSC_FALSE;
107: PetscOptionsGetString(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,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(NULL,"-mat_sigma",&sigma,&flg);
136: if (flg) {
137: PetscOptionsGetString(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_");
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(NULL, "-test_zeropivot", &flg,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(NULL, "-set_row_zero", &flg1,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(NULL, "-check_symmetry", &flg,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: */
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,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: PetscMalloc1(size,&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: Create linear solver; set operators; set runtime options.
269: */
270: KSPCreate(PETSC_COMM_WORLD,&ksp);
271: KSPSetInitialGuessNonzero(ksp,initialguess);
272: num_numfac = 1;
273: PetscOptionsGetInt(NULL,"-num_numfac",&num_numfac,NULL);
274: while (num_numfac--) {
276: KSPSetOperators(ksp,A,A);
277: KSPSetFromOptions(ksp);
279: /*
280: Here we explicitly call KSPSetUp() and KSPSetUpOnBlocks() to
281: enable more precise profiling of setting up the preconditioner.
282: These calls are optional, since both will be called within
283: KSPSolve() if they haven't been called already.
284: */
285: KSPSetUp(ksp);
286: KSPSetUpOnBlocks(ksp);
288: /*
289: Tests "diagonal-scaling of preconditioned residual norm" as used
290: by many ODE integrator codes including SUNDIALS. Note this is different
291: than diagonally scaling the matrix before computing the preconditioner
292: */
293: diagonalscale = PETSC_FALSE;
294: PetscOptionsGetBool(NULL,"-diagonal_scale",&diagonalscale,NULL);
295: if (diagonalscale) {
296: PC pc;
297: PetscInt j,start,end,n;
298: Vec scale;
300: KSPGetPC(ksp,&pc);
301: VecGetSize(x,&n);
302: VecDuplicate(x,&scale);
303: VecGetOwnershipRange(scale,&start,&end);
304: for (j=start; j<end; j++) {
305: VecSetValue(scale,j,((PetscReal)(j+1))/((PetscReal)n),INSERT_VALUES);
306: }
307: VecAssemblyBegin(scale);
308: VecAssemblyEnd(scale);
309: PCSetDiagonalScale(pc,scale);
310: VecDestroy(&scale);
311: }
313: /* - - - - - - - - - - - New Stage - - - - - - - - - - - - -
314: Solve system
315: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
316: /*
317: Solve linear system;
318: */
319: if (trans) {
320: KSPSolveTranspose(ksp,b,x);
321: KSPGetIterationNumber(ksp,&its);
322: } else {
323: PetscInt num_rhs=1;
324: PetscOptionsGetInt(NULL,"-num_rhs",&num_rhs,NULL);
326: while (num_rhs--) {
327: KSPSolve(ksp,b,x);
328: }
329: KSPGetIterationNumber(ksp,&its);
330: if (ckrnorm) { /* Check residual for each rhs */
331: if (trans) {
332: MatMultTranspose(A,x,b2);
333: } else {
334: MatMult(A,x,b2);
335: }
336: VecAXPY(b2,-1.0,b);
337: VecNorm(b2,NORM_2,&rnorm);
338: PetscPrintf(PETSC_COMM_WORLD," Number of iterations = %3D\n",its);
339: PetscPrintf(PETSC_COMM_WORLD," Residual norm %g\n",(double)rnorm);
340: }
341: if (ckerror && !trans) { /* Check error for each rhs */
342: /* VecView(x,PETSC_VIEWER_STDOUT_WORLD); */
343: VecAXPY(u,-1.0,x);
344: VecNorm(u,NORM_2,&enorm);
345: PetscPrintf(PETSC_COMM_WORLD," Error norm %g\n",(double)enorm);
346: }
348: } /* while (num_rhs--) */
351: /*
352: Write output (optinally using table for solver details).
353: - PetscPrintf() handles output for multiprocessor jobs
354: by printing from only one processor in the communicator.
355: - KSPView() prints information about the linear solver.
356: */
357: if (table && ckrnorm) {
358: char *matrixname,kspinfo[120];
359: PetscViewer viewer;
361: /*
362: Open a string viewer; then write info to it.
363: */
364: PetscViewerStringOpen(PETSC_COMM_WORLD,kspinfo,120,&viewer);
365: KSPView(ksp,viewer);
366: PetscStrrchr(file[PetscPreLoadIt],'/',&matrixname);
367: PetscPrintf(PETSC_COMM_WORLD,"%-8.8s %3D %2.0e %s \n", matrixname,its,rnorm,kspinfo);
369: /*
370: Destroy the viewer
371: */
372: PetscViewerDestroy(&viewer);
373: }
375: PetscOptionsGetString(NULL,"-solution",file[3],PETSC_MAX_PATH_LEN,&flg);
376: if (flg) {
377: PetscViewer viewer;
378: Vec xstar;
380: PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[3],FILE_MODE_READ,&viewer);
381: VecCreate(PETSC_COMM_WORLD,&xstar);
382: VecLoad(xstar,viewer);
383: VecAXPY(xstar, -1.0, x);
384: VecNorm(xstar, NORM_2, &enorm);
385: PetscPrintf(PETSC_COMM_WORLD, "Error norm %g\n", (double)enorm);
386: VecDestroy(&xstar);
387: PetscViewerDestroy(&viewer);
388: }
389: if (outputSoln) {
390: PetscViewer viewer;
392: PetscViewerBinaryOpen(PETSC_COMM_WORLD,"solution.petsc",FILE_MODE_WRITE,&viewer);
393: VecView(x, viewer);
394: PetscViewerDestroy(&viewer);
395: }
397: flg = PETSC_FALSE;
398: PetscOptionsGetBool(NULL, "-ksp_reason", &flg,NULL);
399: if (flg) {
400: KSPConvergedReason reason;
401: KSPGetConvergedReason(ksp,&reason);
402: PetscPrintf(PETSC_COMM_WORLD,"KSPConvergedReason: %D\n", reason);
403: }
405: } /* while (num_numfac--) */
407: /*
408: Free work space. All PETSc objects should be destroyed when they
409: are no longer needed.
410: */
411: MatDestroy(&A); VecDestroy(&b);
412: VecDestroy(&u); VecDestroy(&x);
413: VecDestroy(&b2);
414: KSPDestroy(&ksp);
415: if (flgB) { MatDestroy(&B); }
416: PetscPreLoadEnd();
417: /* -----------------------------------------------------------
418: End of linear solver loop
419: ----------------------------------------------------------- */
421: PetscFinalize();
422: return 0;
423: }