Actual source code: ex72.c
petsc-3.10.5 2019-03-28
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: -nearnulldim <0> : number of vectors in the near-null space immediately following matrix\n\n\
12: -trans : solve transpose system instead\n\n";
13: /*
14: This code can be used to test PETSc interface to other packages.\n\
15: Examples of command line options: \n\
16: ./ex72 -f0 <datafile> -ksp_type preonly \n\
17: -help -ksp_view \n\
18: -num_numfac <num_numfac> -num_rhs <num_rhs> \n\
19: -ksp_type preonly -pc_type lu -pc_factor_mat_solver_type superlu or superlu_dist or mumps \n\
20: -ksp_type preonly -pc_type cholesky -pc_factor_mat_solver_type mumps \n\
21: mpiexec -n <np> ./ex72 -f0 <datafile> -ksp_type cg -pc_type asm -pc_asm_type basic -sub_pc_type icc -mat_type sbaij
22: \n\n";
23: */
24: /*T
25: Concepts: KSP^solving a linear system
26: Processors: n
27: T*/
33: /*
34: Include "petscksp.h" so that we can use KSP solvers. Note that this file
35: automatically includes:
36: petscsys.h - base PETSc routines petscvec.h - vectors
37: petscmat.h - matrices
38: petscis.h - index sets petscksp.h - Krylov subspace methods
39: petscviewer.h - viewers petscpc.h - preconditioners
40: */
41: #include <petscksp.h>
43: int main(int argc,char **args)
44: {
45: KSP ksp; /* linear solver context */
46: Mat A; /* matrix */
47: Vec x,b,u; /* approx solution, RHS, exact solution */
48: PetscViewer viewer; /* viewer */
49: char file[4][PETSC_MAX_PATH_LEN]; /* input file name */
50: PetscBool table =PETSC_FALSE,flg,trans=PETSC_FALSE,initialguess = PETSC_FALSE;
51: PetscBool outputSoln=PETSC_FALSE,constantnullspace = PETSC_FALSE;
53: PetscInt its,num_numfac,m,n,M,nearnulldim = 0;
54: PetscReal norm;
55: PetscBool preload=PETSC_TRUE,isSymmetric,cknorm=PETSC_FALSE,initialguessfile = PETSC_FALSE;
56: PetscMPIInt rank;
57: char initialguessfilename[PETSC_MAX_PATH_LEN];
59: PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
60: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
61: PetscOptionsGetBool(NULL,NULL,"-table",&table,NULL);
62: PetscOptionsGetBool(NULL,NULL,"-constantnullspace",&constantnullspace,NULL);
63: PetscOptionsGetBool(NULL,NULL,"-trans",&trans,NULL);
64: PetscOptionsGetBool(NULL,NULL,"-initialguess",&initialguess,NULL);
65: PetscOptionsGetBool(NULL,NULL,"-output_solution",&outputSoln,NULL);
66: PetscOptionsGetString(NULL,NULL,"-initialguessfilename",initialguessfilename,PETSC_MAX_PATH_LEN,&initialguessfile);
67: PetscOptionsGetInt(NULL,NULL,"-nearnulldim",&nearnulldim,NULL);
69: /*
70: Determine files from which we read the two linear systems
71: (matrix and right-hand-side vector).
72: */
73: PetscOptionsGetString(NULL,NULL,"-f",file[0],PETSC_MAX_PATH_LEN,&flg);
74: if (flg) {
75: PetscStrcpy(file[1],file[0]);
76: preload = PETSC_FALSE;
77: } else {
78: PetscOptionsGetString(NULL,NULL,"-f0",file[0],PETSC_MAX_PATH_LEN,&flg);
79: if (!flg) SETERRQ(PETSC_COMM_WORLD,1,"Must indicate binary file with the -f0 or -f option");
80: PetscOptionsGetString(NULL,NULL,"-f1",file[1],PETSC_MAX_PATH_LEN,&flg);
81: if (!flg) preload = PETSC_FALSE; /* don't bother with second system */
82: }
84: /* -----------------------------------------------------------
85: Beginning of linear solver loop
86: ----------------------------------------------------------- */
87: /*
88: Loop through the linear solve 2 times.
89: - The intention here is to preload and solve a small system;
90: then load another (larger) system and solve it as well.
91: This process preloads the instructions with the smaller
92: system so that more accurate performance monitoring (via
93: -log_view) can be done with the larger one (that actually
94: is the system of interest).
95: */
96: PetscPreLoadBegin(preload,"Load system");
98: /* - - - - - - - - - - - New Stage - - - - - - - - - - - - -
99: Load system
100: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
102: /*
103: Open binary file. Note that we use FILE_MODE_READ to indicate
104: reading from this file.
105: */
106: PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[PetscPreLoadIt],FILE_MODE_READ,&viewer);
108: /*
109: Load the matrix and vector; then destroy the viewer.
110: */
111: MatCreate(PETSC_COMM_WORLD,&A);
112: MatSetFromOptions(A);
113: MatLoad(A,viewer);
114: if (nearnulldim) {
115: MatNullSpace nullsp;
116: Vec *nullvecs;
117: PetscInt i;
118: PetscMalloc1(nearnulldim,&nullvecs);
119: for (i=0; i<nearnulldim; i++) {
120: VecCreate(PETSC_COMM_WORLD,&nullvecs[i]);
121: VecLoad(nullvecs[i],viewer);
122: }
123: MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_FALSE,nearnulldim,nullvecs,&nullsp);
124: MatSetNearNullSpace(A,nullsp);
125: for (i=0; i<nearnulldim; i++) {VecDestroy(&nullvecs[i]);}
126: PetscFree(nullvecs);
127: MatNullSpaceDestroy(&nullsp);
128: }
129: if (constantnullspace) {
130: MatNullSpace constant;
131: MatNullSpaceCreate(PETSC_COMM_WORLD,PETSC_TRUE,0,NULL,&constant);
132: MatSetNullSpace(A,constant);
133: MatNullSpaceDestroy(&constant);
134: }
135: flg = PETSC_FALSE;
136: PetscOptionsGetString(NULL,NULL,"-rhs",file[2],PETSC_MAX_PATH_LEN,&flg);
137: VecCreate(PETSC_COMM_WORLD,&b);
138: if (flg) { /* rhs is stored in a separate file */
139: if (file[2][0] == '0' || file[2][0] == 0) {
140: PetscInt m;
141: PetscScalar one = 1.0;
142: PetscInfo(0,"Using vector of ones for RHS\n");
143: MatGetLocalSize(A,&m,NULL);
144: VecSetSizes(b,m,PETSC_DECIDE);
145: VecSetFromOptions(b);
146: VecSet(b,one);
147: } else {
148: PetscViewerDestroy(&viewer);
149: PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[2],FILE_MODE_READ,&viewer);
150: VecSetFromOptions(b);
151: VecLoad(b,viewer);
152: }
153: } else { /* rhs is stored in the same file as matrix */
154: VecSetFromOptions(b);
155: VecLoad(b,viewer);
156: }
157: PetscViewerDestroy(&viewer);
159: /* Make A singular for testing zero-pivot of ilu factorization */
160: /* Example: ./ex72 -f0 <datafile> -test_zeropivot -pc_factor_shift_type <shift_type> */
161: flg = PETSC_FALSE;
162: PetscOptionsGetBool(NULL,NULL, "-test_zeropivot", &flg,NULL);
163: if (flg) { /* set a row as zeros */
164: PetscInt row=0;
165: MatSetOption(A,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE);
166: MatZeroRows(A,1,&row,0.0,NULL,NULL);
167: }
169: /* Check whether A is symmetric, then set A->symmetric option */
170: flg = PETSC_FALSE;
171: PetscOptionsGetBool(NULL,NULL, "-check_symmetry", &flg,NULL);
172: if (flg) {
173: MatIsSymmetric(A,0.0,&isSymmetric);
174: if (!isSymmetric) {
175: PetscPrintf(PETSC_COMM_WORLD,"Warning: A is non-symmetric \n");
176: }
177: }
179: /*
180: If the loaded matrix is larger than the vector (due to being padded
181: to match the block size of the system), then create a new padded vector.
182: */
184: MatGetLocalSize(A,&m,&n);
185: /* if (m != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ, "This example is not intended for rectangular matrices (%d, %d)", m, n);*/
186: MatGetSize(A,&M,NULL);
187: VecGetSize(b,&m);
188: if (M != m) { /* Create a new vector b by padding the old one */
189: PetscInt j,mvec,start,end,indx;
190: Vec tmp;
191: PetscScalar *bold;
193: VecCreate(PETSC_COMM_WORLD,&tmp);
194: VecSetSizes(tmp,n,PETSC_DECIDE);
195: VecSetFromOptions(tmp);
196: VecGetOwnershipRange(b,&start,&end);
197: VecGetLocalSize(b,&mvec);
198: VecGetArray(b,&bold);
199: for (j=0; j<mvec; j++) {
200: indx = start+j;
201: VecSetValues(tmp,1,&indx,bold+j,INSERT_VALUES);
202: }
203: VecRestoreArray(b,&bold);
204: VecDestroy(&b);
205: VecAssemblyBegin(tmp);
206: VecAssemblyEnd(tmp);
207: b = tmp;
208: }
210: MatCreateVecs(A,&x,NULL);
211: VecDuplicate(b,&u);
212: if (initialguessfile) {
213: PetscViewerBinaryOpen(PETSC_COMM_WORLD,initialguessfilename,FILE_MODE_READ,&viewer);
214: VecLoad(x,viewer);
215: PetscViewerDestroy(&viewer);
216: initialguess = PETSC_TRUE;
217: } else if (initialguess) {
218: VecSet(x,1.0);
219: } else {
220: VecSet(x,0.0);
221: }
224: /* Check scaling in A */
225: flg = PETSC_FALSE;
226: PetscOptionsGetBool(NULL,NULL, "-check_scaling", &flg,NULL);
227: if (flg) {
228: Vec max, min;
229: PetscInt idx;
230: PetscReal val;
232: VecDuplicate(x, &max);
233: VecDuplicate(x, &min);
234: MatGetRowMaxAbs(A, max, NULL);
235: MatGetRowMinAbs(A, min, NULL);
236: {
237: PetscViewerASCIIOpen(PETSC_COMM_WORLD, "max.data", &viewer);
238: VecView(max, viewer);
239: PetscViewerDestroy(&viewer);
240: PetscViewerASCIIOpen(PETSC_COMM_WORLD, "min.data", &viewer);
241: VecView(min, viewer);
242: PetscViewerDestroy(&viewer);
243: }
244: VecView(max, PETSC_VIEWER_DRAW_WORLD);
245: VecMax(max, &idx, &val);
246: PetscPrintf(PETSC_COMM_WORLD, "Largest max row element %g at row %D\n", (double)val, idx);
247: VecView(min, PETSC_VIEWER_DRAW_WORLD);
248: VecMin(min, &idx, &val);
249: PetscPrintf(PETSC_COMM_WORLD, "Smallest min row element %g at row %D\n", (double)val, idx);
250: VecMin(max, &idx, &val);
251: PetscPrintf(PETSC_COMM_WORLD, "Smallest max row element %g at row %D\n", (double)val, idx);
252: VecPointwiseDivide(max, max, min);
253: VecMax(max, &idx, &val);
254: PetscPrintf(PETSC_COMM_WORLD, "Largest row ratio %g at row %D\n", (double)val, idx);
255: VecView(max, PETSC_VIEWER_DRAW_WORLD);
256: VecDestroy(&max);
257: VecDestroy(&min);
258: }
260: /* MatView(A,PETSC_VIEWER_STDOUT_WORLD); */
261: /* - - - - - - - - - - - New Stage - - - - - - - - - - - - -
262: Setup solve for system
263: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
264: /*
265: Conclude profiling last stage; begin profiling next stage.
266: */
267: PetscPreLoadStage("KSPSetUpSolve");
269: /*
270: Create linear solver; set operators; set runtime options.
271: */
272: KSPCreate(PETSC_COMM_WORLD,&ksp);
273: KSPSetInitialGuessNonzero(ksp,initialguess);
274: num_numfac = 1;
275: PetscOptionsGetInt(NULL,NULL,"-num_numfac",&num_numfac,NULL);
276: while (num_numfac--) {
277: PC pc;
278: PetscBool lsqr,isbddc,ismatis;
279: char str[32];
281: PetscOptionsGetString(NULL,NULL,"-ksp_type",str,32,&lsqr);
282: if (lsqr) {
283: PetscStrcmp("lsqr",str,&lsqr);
284: }
285: if (lsqr) {
286: Mat BtB;
287: MatTransposeMatMult(A,A,MAT_INITIAL_MATRIX,4,&BtB);
288: KSPSetOperators(ksp,A,BtB);
289: MatDestroy(&BtB);
290: } else {
291: KSPSetOperators(ksp,A,A);
292: }
293: KSPSetFromOptions(ksp);
295: /* if we test BDDC, make sure pmat is of type MATIS */
296: KSPGetPC(ksp,&pc);
297: PetscObjectTypeCompare((PetscObject)pc,PCBDDC,&isbddc);
298: PetscObjectTypeCompare((PetscObject)A,MATIS,&ismatis);
299: if (isbddc && !ismatis) {
300: Mat J;
302: MatConvert(A,MATIS,MAT_INITIAL_MATRIX,&J);
303: KSPSetOperators(ksp,A,J);
304: MatDestroy(&J);
305: }
307: /*
308: Here we explicitly call KSPSetUp() and KSPSetUpOnBlocks() to
309: enable more precise profiling of setting up the preconditioner.
310: These calls are optional, since both will be called within
311: KSPSolve() if they haven't been called already.
312: */
313: KSPSetUp(ksp);
314: KSPSetUpOnBlocks(ksp);
316: /* - - - - - - - - - - - New Stage - - - - - - - - - - - - -
317: Solve system
318: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
320: /*
321: Solve linear system;
322: */
323: if (trans) {
324: KSPSolveTranspose(ksp,b,x);
325: KSPGetIterationNumber(ksp,&its);
326: } else {
327: PetscInt num_rhs=1;
328: PetscOptionsGetInt(NULL,NULL,"-num_rhs",&num_rhs,NULL);
329: cknorm = PETSC_FALSE;
330: PetscOptionsGetBool(NULL,NULL,"-cknorm",&cknorm,NULL);
331: while (num_rhs--) {
332: if (num_rhs == 1) VecSet(x,0.0);
333: KSPSolve(ksp,b,x);
334: }
335: KSPGetIterationNumber(ksp,&its);
336: if (cknorm) { /* Check error for each rhs */
337: if (trans) {
338: MatMultTranspose(A,x,u);
339: } else {
340: MatMult(A,x,u);
341: }
342: VecAXPY(u,-1.0,b);
343: VecNorm(u,NORM_2,&norm);
344: PetscPrintf(PETSC_COMM_WORLD," Number of iterations = %3D\n",its);
345: if (!PetscIsNanScalar(norm)) {
346: if (norm < 1.e-12) {
347: PetscPrintf(PETSC_COMM_WORLD," Residual norm < 1.e-12\n");
348: } else {
349: PetscPrintf(PETSC_COMM_WORLD," Residual norm %g\n",(double)norm);
350: }
351: }
352: }
353: } /* while (num_rhs--) */
355: /* - - - - - - - - - - - New Stage - - - - - - - - - - - - -
356: Check error, print output, free data structures.
357: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
359: /*
360: Check error
361: */
362: if (trans) {
363: MatMultTranspose(A,x,u);
364: } else {
365: MatMult(A,x,u);
366: }
367: VecAXPY(u,-1.0,b);
368: VecNorm(u,NORM_2,&norm);
369: /*
370: Write output (optinally using table for solver details).
371: - PetscPrintf() handles output for multiprocessor jobs
372: by printing from only one processor in the communicator.
373: - KSPView() prints information about the linear solver.
374: */
375: if (table) {
376: char *matrixname,kspinfo[120];
378: /*
379: Open a string viewer; then write info to it.
380: */
381: PetscViewerStringOpen(PETSC_COMM_WORLD,kspinfo,120,&viewer);
382: KSPView(ksp,viewer);
383: PetscStrrchr(file[PetscPreLoadIt],'/',&matrixname);
384: PetscPrintf(PETSC_COMM_WORLD,"%-8.8s %3D %2.0e %s \n",matrixname,its,norm,kspinfo);
386: /*
387: Destroy the viewer
388: */
389: PetscViewerDestroy(&viewer);
390: } else {
391: PetscPrintf(PETSC_COMM_WORLD,"Number of iterations = %3D\n",its);
392: if (!PetscIsNanScalar(norm)) {
393: if (norm < 1.e-12 && !PetscIsNanScalar((PetscScalar)norm)) {
394: PetscPrintf(PETSC_COMM_WORLD," Residual norm < 1.e-12\n");
395: } else {
396: PetscPrintf(PETSC_COMM_WORLD,"Residual norm %g\n",(double)norm);
397: }
398: }
399: }
400: PetscOptionsGetString(NULL,NULL,"-solution",file[3],PETSC_MAX_PATH_LEN,&flg);
401: if (flg) {
402: Vec xstar;
403: PetscReal norm;
405: PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[3],FILE_MODE_READ,&viewer);
406: VecCreate(PETSC_COMM_WORLD,&xstar);
407: VecLoad(xstar,viewer);
408: VecAXPY(xstar, -1.0, x);
409: VecNorm(xstar, NORM_2, &norm);
410: PetscPrintf(PETSC_COMM_WORLD, "Error norm %g\n", (double)norm);
411: VecDestroy(&xstar);
412: PetscViewerDestroy(&viewer);
413: }
414: if (outputSoln) {
415: PetscViewerBinaryOpen(PETSC_COMM_WORLD,"solution.petsc",FILE_MODE_WRITE,&viewer);
416: VecView(x, viewer);
417: PetscViewerDestroy(&viewer);
418: }
420: flg = PETSC_FALSE;
421: PetscOptionsGetBool(NULL,NULL, "-ksp_reason", &flg,NULL);
422: if (flg) {
423: KSPConvergedReason reason;
424: KSPGetConvergedReason(ksp,&reason);
425: PetscPrintf(PETSC_COMM_WORLD,"KSPConvergedReason: %D\n", reason);
426: }
428: } /* while (num_numfac--) */
430: /*
431: Free work space. All PETSc objects should be destroyed when they
432: are no longer needed.
433: */
434: MatDestroy(&A); VecDestroy(&b);
435: VecDestroy(&u); VecDestroy(&x);
436: KSPDestroy(&ksp);
437: PetscPreLoadEnd();
438: /* -----------------------------------------------------------
439: End of linear solver loop
440: ----------------------------------------------------------- */
442: PetscFinalize();
443: return ierr;
444: }
447: /*TEST
449: build:
450: requires: !complex
452: testset:
453: suffix: 1
454: nsize: 2
455: args: -f0 ${wPETSC_DIR}/share/petsc/datafiles/matrices/spd-real-int@PETSC_INDEX_SIZE@-float@PETSC_SCALAR_SIZE@
456: requires: !__float128
458: testset:
459: suffix: 1a
460: args: -f0 ${wPETSC_DIR}/share/petsc/datafiles/matrices/spd-real-int@PETSC_INDEX_SIZE@-float@PETSC_SCALAR_SIZE@
461: requires: !__float128
463: testset:
464: nsize: 2
465: requires: datafilespath double !define(PETSC_USE_64BIT_INDICES)
466: args: -f0 ${DATAFILESPATH}/matrices/medium
467: args: -ksp_type bicg
468: test:
469: suffix: 2
471: testset:
472: requires: datafilespath double !define(PETSC_USE_64BIT_INDICES)
473: args: -f0 ${DATAFILESPATH}/matrices/medium
474: args: -ksp_type bicg
475: test:
476: suffix: 4
477: args: -pc_type lu
478: test:
479: suffix: 5
481: testset:
482: suffix: 6
483: requires: datafilespath double !define(PETSC_USE_64BIT_INDICES)
484: args: -f0 ${DATAFILESPATH}/matrices/fem1
485: args: -pc_factor_levels 2 -pc_factor_fill 1.73 -ksp_gmres_cgs_refinement_type refine_always
487: testset:
488: suffix: 7
489: requires: datafilespath double !define(PETSC_USE_64BIT_INDICES)
490: args: -f0 ${DATAFILESPATH}/matrices/medium
491: args: -viewer_binary_skip_info -mat_type seqbaij
492: args: -matload_block_size {{2 3 4 5 6 7 8}separate output}
493: args: -ksp_max_it 100 -ksp_gmres_cgs_refinement_type refine_always
494: args: -ksp_rtol 1.0e-15 -ksp_monitor_short
495: test:
496: suffix: a
497: test:
498: suffix: b
499: args: -pc_factor_mat_ordering_type nd
500: test:
501: suffix: c
502: args: -pc_factor_levels 1
504: testset:
505: suffix: 7_d
506: requires: datafilespath double !define(PETSC_USE_64BIT_INDICES)
507: args: -f0 ${DATAFILESPATH}/matrices/medium
508: args: -viewer_binary_skip_info -mat_type seqbaij
509: args: -matload_block_size {{2 3 4 5 6 7 8}shared output}
510: args: -ksp_type preonly -pc_type lu
512: testset:
513: suffix: 8
514: requires: datafilespath double !define(PETSC_USE_64BIT_INDICES)
515: args: -f0 ${DATAFILESPATH}/matrices/medium
516: args: -ksp_diagonal_scale -pc_type eisenstat -ksp_monitor_short -ksp_diagonal_scale_fix -ksp_gmres_cgs_refinement_type refine_always -mat_no_inode
518: testset:
519: suffix: 9
520: requires: datafilespath double !define(PETSC_USE_64BIT_INDICES)
521: args: -f0 ${DATAFILESPATH}/matrices/medium
522: args: -viewer_binary_skip_info -matload_block_size {{1 2 3 4 5 6 7}separate output} -ksp_max_it 100 -ksp_gmres_cgs_refinement_type refine_always -ksp_rtol 1.0e-15 -ksp_monitor_short
523: test:
524: suffix: a
525: args: -mat_type seqbaij
526: test:
527: suffix: b
528: args: -mat_type seqbaij -trans
529: test:
530: suffix: c
531: nsize: 2
532: args: -mat_type mpibaij
533: test:
534: suffix: d
535: nsize: 2
536: args: -mat_type mpibaij -trans
537: test:
538: suffix: e
539: nsize: 3
540: args: -mat_type mpibaij
541: test:
542: suffix: f
543: nsize: 3
544: args: -mat_type mpibaij -trans
547: testset:
548: suffix: 10
549: nsize: 2
550: requires: datafilespath double !define(PETSC_USE_64BIT_INDICES)
551: args: -ksp_type fgmres -pc_type ksp -f0 ${DATAFILESPATH}/matrices/medium -ksp_fgmres_modifypcksp -ksp_monitor_short
553: testset:
554: TODO: Need to determine goal of this test
555: suffix: 11
556: nsize: 2
557: args: -f0 http://ftp.mcs.anl.gov/pub/petsc/Datafiles/matrices/testmatrix.gz
559: testset:
560: suffix: 12
561: requires: matlab
562: args: -pc_type lu -pc_factor_mat_solver_type matlab -f0 ${DATAFILESPATH}/matrices/arco1
564: testset:
565: suffix: 13
566: requires: lusol
567: args: -f0 ${DATAFILESPATH}/matrices/arco1
568: args: -mat_type lusol -pc_type lu
570: testset:
571: nsize: 3
572: args: -f0 ${DATAFILESPATH}/matrices/medium
573: requires: datafilespath double !define(PETSC_USE_64BIT_INDICES)
574: test:
575: suffix: 14
576: requires: spai
577: args: -pc_type spai
578: test:
579: suffix: 15
580: requires: hypre
581: args: -pc_type hypre -pc_hypre_type pilut
582: test:
583: suffix: 16
584: requires: hypre
585: args: -pc_type hypre -pc_hypre_type parasails
586: test:
587: suffix: 17
588: requires: hypre
589: args: -pc_type hypre -pc_hypre_type boomeramg
591: testset:
592: suffix: 19
593: requires: datafilespath double !define(PETSC_USE_64BIT_INDICES)
594: args: -f0 ${DATAFILESPATH}/matrices/poisson1
595: args: -ksp_type cg -pc_type icc
596: args: -pc_factor_levels {{0 2 4}separate output}
597: test:
598: test:
599: args: -mat_type seqsbaij
602: testset:
603: suffix: ILU
604: requires: datafilespath double !define(PETSC_USE_64BIT_INDICES)
605: args: -f0 ${DATAFILESPATH}/matrices/small
606: args: -pc_factor_levels 1
607: test:
608: test:
609: # This is tested against regular ILU (used to be denoted ILUBAIJ)
610: args: -mat_type baij
612: testset:
613: suffix: aijcusparse
614: requires: datafilespath double !define(PETSC_USE_64BIT_INDICES) veccuda
615: args: -f0 ${DATAFILESPATH}/matrices/medium -ksp_monitor_short -ksp_view -mat_view ascii::ascii_info -mat_type aijcusparse -pc_factor_mat_solver_type cusparse -pc_type ilu -vec_type cuda
617: testset:
618: TODO: No output file. Need to determine if deprecated
619: suffix: asm_viennacl
620: nsize: 2
621: requires: viennacl
622: args: -pc_type asm -pc_asm_sub_mat_type aijviennacl -f0 ${wPETSC_DIR}/share/petsc/datafiles/matrices/spd-real-int${PETSC_INDEX_SIZE}-float${PETSC_SCALAR_SIZE}
624: testset:
625: nsize: 2
626: requires: datafilespath double !define(PETSC_USE_64BIT_INDICES) hypre
627: args: -f0 ${DATAFILESPATH}/matrices/poisson2.gz -ksp_monitor_short -ksp_rtol 1.E-9 -pc_type hypre -pc_hypre_type boomeramg
628: test:
629: suffix: boomeramg_euclid
630: args: -pc_hypre_boomeramg_smooth_type Euclid -pc_hypre_boomeramg_smooth_num_levels 2 -pc_hypre_boomeramg_eu_level 1 -pc_hypre_boomeramg_eu_droptolerance 0.01
631: TODO: Need to determine if deprecated
632: test:
633: suffix: boomeramg_euclid_bj
634: args: -pc_hypre_boomeramg_smooth_type Euclid -pc_hypre_boomeramg_smooth_num_levels 2 -pc_hypre_boomeramg_eu_level 1 -pc_hypre_boomeramg_eu_droptolerance 0.01 -pc_hypre_boomeramg_eu_bj
635: TODO: Need to determine if deprecated
636: test:
637: suffix: boomeramg_parasails
638: args: -pc_hypre_boomeramg_smooth_type ParaSails -pc_hypre_boomeramg_smooth_num_levels 2
639: test:
640: suffix: boomeramg_pilut
641: args: -pc_hypre_boomeramg_smooth_type Pilut -pc_hypre_boomeramg_smooth_num_levels 2
642: test:
643: suffix: boomeramg_schwarz
644: args: -pc_hypre_boomeramg_smooth_type Schwarz-smoothers
646: testset:
647: suffix: cg_singlereduction
648: requires: datafilespath double !define(PETSC_USE_64BIT_INDICES)
649: args: -f0 ${DATAFILESPATH}/matrices/small
650: args: -mat_type mpisbaij -ksp_type cg -pc_type eisenstat -ksp_monitor_short -ksp_converged_reason
651: test:
652: test:
653: args: -ksp_cg_single_reduction
655: testset:
656: requires: datafilespath double !define(PETSC_USE_64BIT_INDICES)
657: args: -f0 ${DATAFILESPATH}/matrices/poisson2.gz
658: args: -ksp_monitor_short -pc_type icc
659: test:
660: suffix: cr
661: args: -ksp_type cr
662: test:
663: suffix: lcd
664: args: -ksp_type lcd
666: testset:
667: requires: datafilespath double !define(PETSC_USE_64BIT_INDICES)
668: args: -f0 ${DATAFILESPATH}/matrices/small
669: args: -ksp_monitor_short -ksp_view -mat_view ascii::ascii_info
670: test:
671: suffix: seqaijcrl
672: args: -mat_type seqaijcrl
673: test:
674: suffix: seqaijperm
675: args: -mat_type seqaijperm
677: testset:
678: nsize: 2
679: requires: datafilespath double !define(PETSC_USE_64BIT_INDICES)
680: args: -f0 ${DATAFILESPATH}/matrices/small
681: args: -ksp_monitor_short -ksp_view
682: # Different output files
683: test:
684: suffix: mpiaijcrl
685: args: -mat_type mpiaijcrl
686: test:
687: suffix: mpiaijperm
688: args: -mat_type mpiaijperm
690: testset:
691: nsize: 4
692: requires: datafilespath double !define(PETSC_USE_64BIT_INDICES)
693: args: -ksp_monitor_short -ksp_view
694: test:
695: suffix: xxt
696: args: -f0 ${DATAFILESPATH}/matrices/poisson1 -check_symmetry -ksp_type cg -pc_type tfs
697: test:
698: suffix: xyt
699: args: -f0 ${DATAFILESPATH}/matrices/medium -ksp_type gmres -pc_type tfs
701: testset:
702: # The output file here is the same as mumps
703: suffix: mumps_cholesky
704: output_file: output/ex72_mumps.out
705: requires: datafilespath double !define(PETSC_USE_64BIT_INDICES) mumps
706: args: -f0 ${DATAFILESPATH}/matrices/small -ksp_type preonly -pc_type cholesky -pc_factor_mat_solver_type mumps -num_numfac 2 -num_rhs 2
707: nsize: {{1 2}}
708: test:
709: args: -mat_type sbaij -mat_ignore_lower_triangular
710: test:
711: args: -mat_type aij
712: test:
713: args: -mat_type aij -matload_spd
715: testset:
716: # The output file here is the same as mumps
717: suffix: mumps_lu
718: output_file: output/ex72_mumps.out
719: requires: datafilespath double !define(PETSC_USE_64BIT_INDICES) mumps
720: args: -f0 ${DATAFILESPATH}/matrices/small -ksp_type preonly -pc_type lu -pc_factor_mat_solver_type mumps -num_numfac 2 -num_rhs 2
721: test:
722: args: -mat_type seqaij
723: test:
724: nsize: 2
725: args: -mat_type mpiaij
726: test:
727: args: -mat_type seqbaij -matload_block_size 2
728: test:
729: nsize: 2
730: args: -mat_type mpibaij -matload_block_size 2
731: test:
732: args: -mat_type aij -mat_mumps_icntl_7 5
733: TODO: Need to determine if deprecated
734: test:
735: nsize: 2
736: args: -mat_type mpiaij -mat_mumps_icntl_28 2 -mat_mumps_icntl_29 2
737: TODO: Need to determine if deprecated
739: testset:
740: # The output file here is the same as mumps
741: suffix: mumps_redundant
742: output_file: output/ex72_mumps_redundant.out
743: nsize: 8
744: requires: datafilespath double !define(PETSC_USE_64BIT_INDICES) mumps
745: args: -f0 ${DATAFILESPATH}/matrices/medium -ksp_type preonly -pc_type redundant -pc_redundant_number {{8 7 6 5 4 3 2 1}} -redundant_pc_factor_mat_solver_type mumps -num_numfac 2 -num_rhs 2
747: testset:
748: suffix: pastix_cholesky
749: requires: datafilespath double !define(PETSC_USE_64BIT_INDICES) pastix
750: output_file: output/ex72_mumps.out
751: nsize: {{1 2}}
752: args: -f0 ${DATAFILESPATH}/matrices/small -ksp_type preonly -pc_factor_mat_solver_type pastix -num_numfac 2 -num_rhs 2 -pc_type cholesky -mat_type sbaij -mat_ignore_lower_triangular
754: testset:
755: suffix: pastix_lu
756: requires: datafilespath double !define(PETSC_USE_64BIT_INDICES) pastix
757: args: -f0 ${DATAFILESPATH}/matrices/small -ksp_type preonly -pc_type lu -pc_factor_mat_solver_type pastix -num_numfac 2 -num_rhs 2
758: output_file: output/ex72_mumps.out
759: test:
760: args: -mat_type seqaij
761: test:
762: nsize: 2
763: args: -mat_type mpiaij
765: testset:
766: suffix: pastix_redundant
767: output_file: output/ex72_mumps_redundant.out
768: nsize: 8
769: requires: datafilespath double !define(PETSC_USE_64BIT_INDICES) pastix
770: args: -f0 ${DATAFILESPATH}/matrices/medium -ksp_type preonly -pc_type redundant -pc_redundant_number {{8 7 6 5 4 3 2 1}} -redundant_pc_factor_mat_solver_type pastix -num_numfac 2 -num_rhs 2
773: testset:
774: suffix: superlu_dist_lu
775: requires: datafilespath double !define(PETSC_USE_64BIT_INDICES) superlu_dist
776: output_file: output/ex72_mumps.out
777: args: -f0 ${DATAFILESPATH}/matrices/small -ksp_type preonly -pc_type lu -pc_factor_mat_solver_type superlu_dist -num_numfac 2 -num_rhs 2
778: nsize: {{1 2}}
780: testset:
781: suffix: superlu_dist_redundant
782: nsize: 8
783: output_file: output/ex72_mumps_redundant.out
784: requires: datafilespath double !define(PETSC_USE_64BIT_INDICES) superlu_dist
785: args: -f0 ${DATAFILESPATH}/matrices/medium -ksp_type preonly -pc_type redundant -pc_redundant_number {{8 7 6 5 4 3 2 1}} -redundant_pc_factor_mat_solver_type superlu_dist -num_numfac 2 -num_rhs 2
787: testset:
788: suffix: superlu_lu
789: output_file: output/ex72_mumps.out
790: requires: datafilespath double !define(PETSC_USE_64BIT_INDICES) superlu
791: args: -f0 ${DATAFILESPATH}/matrices/small -ksp_type preonly -pc_type lu -pc_factor_mat_solver_type superlu -num_numfac 2 -num_rhs 2
793: testset:
794: suffix: umfpack
795: requires: datafilespath double !define(PETSC_USE_64BIT_INDICES) suitesparse
796: args: -f0 ${DATAFILESPATH}/matrices/small -ksp_type preonly -pc_type lu -mat_type seqaij -pc_factor_mat_solver_type umfpack -num_numfac 2 -num_rhs 2
799: testset:
800: suffix: zeropivot
801: requires: datafilespath double !define(PETSC_USE_64BIT_INDICES) mumps
802: args: -f0 ${DATAFILESPATH}/matrices/small -test_zeropivot -ksp_converged_reason -ksp_type fgmres -pc_type ksp
803: test:
804: nsize: 3
805: args: -ksp_pc_type bjacobi
806: test:
807: nsize: 2
808: args: -ksp_ksp_type cg -ksp_pc_type bjacobi -ksp_pc_bjacobi_blocks 1
809: #test:
810: #nsize: 3
811: #args: -ksp_ksp_converged_reason -ksp_pc_type bjacobi -ksp_sub_ksp_converged_reason
812: #TODO: Need to determine if deprecated
814: testset:
815: requires: datafilespath double !define(PETSC_USE_64BIT_INDICES)
816: args: -f0 ${DATAFILESPATH}/matrices/medium -ksp_type fgmres
817: test:
818: suffix: bddc_seq
819: nsize: 1
820: args: -pc_type bddc
821: test:
822: suffix: bddc_par
823: nsize: 2
824: args: -pc_type bddc
825: test:
826: requires: parmetis
827: suffix: bddc_par_nd_parmetis
828: filter: sed -e "s/Number of iterations = [0-9]/Number of iterations = 9/g"
829: nsize: 4
830: args: -ksp_error_if_not_converged -pc_type bddc -mat_is_disassemble_l2g_type nd -mat_partitioning_type parmetis
831: test:
832: requires: ptscotch define(PETSC_HAVE_SCOTCH_PARMETIS_V3_NODEND)
833: suffix: bddc_par_nd_ptscotch
834: filter: sed -e "s/Number of iterations = [0-9]/Number of iterations = 9/g"
835: nsize: 4
836: args: -ksp_error_if_not_converged -pc_type bddc -mat_is_disassemble_l2g_type nd -mat_partitioning_type ptscotch
837: TEST*/