Actual source code: ex30.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: 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,NULL,"-table",&table,NULL);
52: PetscOptionsGetBool(NULL,NULL,"-trans",&trans,NULL);
53: PetscOptionsGetBool(NULL,NULL,"-partition",&partition,NULL);
54: PetscOptionsGetBool(NULL,NULL,"-initialguess",&initialguess,NULL);
55: PetscOptionsGetBool(NULL,NULL,"-output_solution",&outputSoln,NULL);
56: PetscOptionsGetBool(NULL,NULL,"-ckrnorm",&ckrnorm,NULL);
57: PetscOptionsGetBool(NULL,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,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,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,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,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,NULL,"-mat_sigma",&sigma,&flg);
136: if (flg) {
137: PetscOptionsGetString(NULL,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,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: PetscCalloc1(ncols+1,&zeros);
164: flg1 = PETSC_FALSE;
165: PetscOptionsGetBool(NULL,NULL, "-set_row_zero", &flg1,NULL);
166: if (flg1) { /* set entire row as zero */
167: MatSetValues(A,1,&row,ncols,cols,zeros,INSERT_VALUES);
168: } else { /* only set (row,row) entry as zero */
169: MatSetValues(A,1,&row,1,&row,zeros,INSERT_VALUES);
170: }
171: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
172: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
173: }
175: /* Check whether A is symmetric */
176: flg = PETSC_FALSE;
177: PetscOptionsGetBool(NULL,NULL, "-check_symmetry", &flg,NULL);
178: if (flg) {
179: Mat Atrans;
180: MatTranspose(A, MAT_INITIAL_MATRIX,&Atrans);
181: MatEqual(A, Atrans, &isSymmetric);
182: if (isSymmetric) {
183: PetscPrintf(PETSC_COMM_WORLD,"A is symmetric \n");
184: } else {
185: PetscPrintf(PETSC_COMM_WORLD,"A is non-symmetric \n");
186: }
187: MatDestroy(&Atrans);
188: }
190: /*
191: If the loaded matrix is larger than the vector (due to being padded
192: to match the block size of the system), then create a new padded vector.
193: */
195: MatGetLocalSize(A,&m,&n);
196: if (m != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ, "This example is not intended for rectangular matrices (%d, %d)", m, n);
197: MatGetSize(A,&M,NULL);
198: VecGetSize(b,&m);
199: if (M != m) { /* Create a new vector b by padding the old one */
200: PetscInt j,mvec,start,end,indx;
201: Vec tmp;
202: PetscScalar *bold;
204: VecCreate(PETSC_COMM_WORLD,&tmp);
205: VecSetSizes(tmp,n,PETSC_DECIDE);
206: VecSetFromOptions(tmp);
207: VecGetOwnershipRange(b,&start,&end);
208: VecGetLocalSize(b,&mvec);
209: VecGetArray(b,&bold);
210: for (j=0; j<mvec; j++) {
211: indx = start+j;
212: VecSetValues(tmp,1,&indx,bold+j,INSERT_VALUES);
213: }
214: VecRestoreArray(b,&bold);
215: VecDestroy(&b);
216: VecAssemblyBegin(tmp);
217: VecAssemblyEnd(tmp);
218: b = tmp;
219: }
220: VecDuplicate(b,&b2);
221: VecDuplicate(b,&x);
222: PetscObjectSetName((PetscObject)x, "Solution vector");
223: VecDuplicate(b,&u);
224: PetscObjectSetName((PetscObject)u, "True Solution vector");
225: VecSet(x,0.0);
227: if (ckerror) { /* Set true solution */
228: VecSet(u,1.0);
229: MatMult(A,u,b);
230: }
232: /* - - - - - - - - - - - New Stage - - - - - - - - - - - - -
233: Setup solve for system
234: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
236: if (partition) {
237: MatPartitioning mpart;
238: IS mis,nis,is;
239: PetscInt *count;
240: PetscMPIInt size;
241: Mat BB;
242: MPI_Comm_size(PETSC_COMM_WORLD,&size);
243: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
244: PetscMalloc1(size,&count);
245: MatPartitioningCreate(PETSC_COMM_WORLD, &mpart);
246: MatPartitioningSetAdjacency(mpart, A);
247: /* MatPartitioningSetVertexWeights(mpart, weight); */
248: MatPartitioningSetFromOptions(mpart);
249: MatPartitioningApply(mpart, &mis);
250: MatPartitioningDestroy(&mpart);
251: ISPartitioningToNumbering(mis,&nis);
252: ISPartitioningCount(mis,size,count);
253: ISDestroy(&mis);
254: ISInvertPermutation(nis, count[rank], &is);
255: PetscFree(count);
256: ISDestroy(&nis);
257: ISSort(is);
258: MatGetSubMatrix(A,is,is,MAT_INITIAL_MATRIX,&BB);
260: /* need to move the vector also */
261: ISDestroy(&is);
262: MatDestroy(&A);
263: A = BB;
264: }
266: /*
267: Create linear solver; set operators; set runtime options.
268: */
269: KSPCreate(PETSC_COMM_WORLD,&ksp);
270: KSPSetInitialGuessNonzero(ksp,initialguess);
271: num_numfac = 1;
272: PetscOptionsGetInt(NULL,NULL,"-num_numfac",&num_numfac,NULL);
273: while (num_numfac--) {
275: KSPSetOperators(ksp,A,A);
276: KSPSetFromOptions(ksp);
278: /*
279: Here we explicitly call KSPSetUp() and KSPSetUpOnBlocks() to
280: enable more precise profiling of setting up the preconditioner.
281: These calls are optional, since both will be called within
282: KSPSolve() if they haven't been called already.
283: */
284: KSPSetUp(ksp);
285: KSPSetUpOnBlocks(ksp);
287: /*
288: Tests "diagonal-scaling of preconditioned residual norm" as used
289: by many ODE integrator codes including SUNDIALS. Note this is different
290: than diagonally scaling the matrix before computing the preconditioner
291: */
292: diagonalscale = PETSC_FALSE;
293: PetscOptionsGetBool(NULL,NULL,"-diagonal_scale",&diagonalscale,NULL);
294: if (diagonalscale) {
295: PC pc;
296: PetscInt j,start,end,n;
297: Vec scale;
299: KSPGetPC(ksp,&pc);
300: VecGetSize(x,&n);
301: VecDuplicate(x,&scale);
302: VecGetOwnershipRange(scale,&start,&end);
303: for (j=start; j<end; j++) {
304: VecSetValue(scale,j,((PetscReal)(j+1))/((PetscReal)n),INSERT_VALUES);
305: }
306: VecAssemblyBegin(scale);
307: VecAssemblyEnd(scale);
308: PCSetDiagonalScale(pc,scale);
309: VecDestroy(&scale);
310: }
312: /* - - - - - - - - - - - New Stage - - - - - - - - - - - - -
313: Solve system
314: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
315: /*
316: Solve linear system;
317: */
318: if (trans) {
319: KSPSolveTranspose(ksp,b,x);
320: KSPGetIterationNumber(ksp,&its);
321: } else {
322: PetscInt num_rhs=1;
323: PetscOptionsGetInt(NULL,NULL,"-num_rhs",&num_rhs,NULL);
325: while (num_rhs--) {
326: KSPSolve(ksp,b,x);
327: }
328: KSPGetIterationNumber(ksp,&its);
329: if (ckrnorm) { /* Check residual for each rhs */
330: if (trans) {
331: MatMultTranspose(A,x,b2);
332: } else {
333: MatMult(A,x,b2);
334: }
335: VecAXPY(b2,-1.0,b);
336: VecNorm(b2,NORM_2,&rnorm);
337: PetscPrintf(PETSC_COMM_WORLD," Number of iterations = %3D\n",its);
338: PetscPrintf(PETSC_COMM_WORLD," Residual norm %g\n",(double)rnorm);
339: }
340: if (ckerror && !trans) { /* Check error for each rhs */
341: /* VecView(x,PETSC_VIEWER_STDOUT_WORLD); */
342: VecAXPY(u,-1.0,x);
343: VecNorm(u,NORM_2,&enorm);
344: PetscPrintf(PETSC_COMM_WORLD," Error norm %g\n",(double)enorm);
345: }
347: } /* while (num_rhs--) */
350: /*
351: Write output (optinally using table for solver details).
352: - PetscPrintf() handles output for multiprocessor jobs
353: by printing from only one processor in the communicator.
354: - KSPView() prints information about the linear solver.
355: */
356: if (table && ckrnorm) {
357: char *matrixname,kspinfo[120];
358: PetscViewer viewer;
360: /*
361: Open a string viewer; then write info to it.
362: */
363: PetscViewerStringOpen(PETSC_COMM_WORLD,kspinfo,120,&viewer);
364: KSPView(ksp,viewer);
365: PetscStrrchr(file[PetscPreLoadIt],'/',&matrixname);
366: PetscPrintf(PETSC_COMM_WORLD,"%-8.8s %3D %2.0e %s \n", matrixname,its,rnorm,kspinfo);
368: /*
369: Destroy the viewer
370: */
371: PetscViewerDestroy(&viewer);
372: }
374: PetscOptionsGetString(NULL,NULL,"-solution",file[3],PETSC_MAX_PATH_LEN,&flg);
375: if (flg) {
376: PetscViewer viewer;
377: Vec xstar;
379: PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[3],FILE_MODE_READ,&viewer);
380: VecCreate(PETSC_COMM_WORLD,&xstar);
381: VecLoad(xstar,viewer);
382: VecAXPY(xstar, -1.0, x);
383: VecNorm(xstar, NORM_2, &enorm);
384: PetscPrintf(PETSC_COMM_WORLD, "Error norm %g\n", (double)enorm);
385: VecDestroy(&xstar);
386: PetscViewerDestroy(&viewer);
387: }
388: if (outputSoln) {
389: PetscViewer viewer;
391: PetscViewerBinaryOpen(PETSC_COMM_WORLD,"solution.petsc",FILE_MODE_WRITE,&viewer);
392: VecView(x, viewer);
393: PetscViewerDestroy(&viewer);
394: }
396: flg = PETSC_FALSE;
397: PetscOptionsGetBool(NULL,NULL, "-ksp_reason", &flg,NULL);
398: if (flg) {
399: KSPConvergedReason reason;
400: KSPGetConvergedReason(ksp,&reason);
401: PetscPrintf(PETSC_COMM_WORLD,"KSPConvergedReason: %D\n", reason);
402: }
404: } /* while (num_numfac--) */
406: /*
407: Free work space. All PETSc objects should be destroyed when they
408: are no longer needed.
409: */
410: MatDestroy(&A); VecDestroy(&b);
411: VecDestroy(&u); VecDestroy(&x);
412: VecDestroy(&b2);
413: KSPDestroy(&ksp);
414: if (flgB) { MatDestroy(&B); }
415: PetscPreLoadEnd();
416: /* -----------------------------------------------------------
417: End of linear solver loop
418: ----------------------------------------------------------- */
420: PetscFinalize();
421: return 0;
422: }