Actual source code: ex10.c

petsc-3.8.4 2018-03-24
Report Typos and Errors

  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:    ./ex10 -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_package superlu or superlu_dist or mumps \n\
 20:         -ksp_type preonly -pc_type cholesky -pc_factor_mat_solver_package mumps \n\
 21:    mpiexec -n <np> ./ex10 -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    fd;              /* 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,&fd);

108:   /*
109:      Load the matrix and vector; then destroy the viewer.
110:   */
111:   MatCreate(PETSC_COMM_WORLD,&A);
112:   MatSetFromOptions(A);
113:   MatLoad(A,fd);
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],fd);
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(&fd);
149:       PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[2],FILE_MODE_READ,&fd);
150:       VecSetFromOptions(b);
151:       VecLoad(b,fd);
152:     }
153:   } else {   /* rhs is stored in the same file as matrix */
154:     VecSetFromOptions(b);
155:     VecLoad(b,fd);
156:   }
157:   PetscViewerDestroy(&fd);

159:   /* Make A singular for testing zero-pivot of ilu factorization        */
160:   /* Example: ./ex10 -f0 <datafile> -test_zeropivot -set_row_zero -pc_factor_shift_type <shift_type> */
161:   flg  = PETSC_FALSE;
162:   PetscOptionsGetBool(NULL,NULL, "-test_zeropivot", &flg,NULL);
163:   if (flg) {
164:     PetscInt          row=0;
165:     PetscBool         flg1=PETSC_FALSE;
166:     PetscScalar       zero=0.0;

168:     PetscOptionsGetBool(NULL,NULL, "-set_row_zero", &flg1,NULL);
169:     if (flg1) {   /* set a row as zeros */
170:       MatSetOption(A,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE);
171:       MatZeroRows(A,1,&row,0.0,NULL,NULL);
172:     } else if (!rank) {
173:       /* set (row,row) entry as zero */
174:       MatSetValues(A,1,&row,1,&row,&zero,INSERT_VALUES);
175:     }
176:     MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
177:     MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
178:   }

180:   /* Check whether A is symmetric */
181:   flg  = PETSC_FALSE;
182:   PetscOptionsGetBool(NULL,NULL, "-check_symmetry", &flg,NULL);
183:   if (flg) {
184:     Mat Atrans;
185:     MatTranspose(A, MAT_INITIAL_MATRIX,&Atrans);
186:     MatEqual(A, Atrans, &isSymmetric);
187:     if (isSymmetric) {
188:       MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE);
189:     } else {
190:       PetscPrintf(PETSC_COMM_WORLD,"Warning: A is non-symmetric \n");
191:     }
192:     MatDestroy(&Atrans);
193:   }

195:   /*
196:      If the loaded matrix is larger than the vector (due to being padded
197:      to match the block size of the system), then create a new padded vector.
198:   */

200:   MatGetLocalSize(A,&m,&n);
201:   /*  if (m != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ, "This example is not intended for rectangular matrices (%d, %d)", m, n);*/
202:   MatGetSize(A,&M,NULL);
203:   VecGetSize(b,&m);
204:   if (M != m) {   /* Create a new vector b by padding the old one */
205:     PetscInt    j,mvec,start,end,indx;
206:     Vec         tmp;
207:     PetscScalar *bold;

209:     VecCreate(PETSC_COMM_WORLD,&tmp);
210:     VecSetSizes(tmp,n,PETSC_DECIDE);
211:     VecSetFromOptions(tmp);
212:     VecGetOwnershipRange(b,&start,&end);
213:     VecGetLocalSize(b,&mvec);
214:     VecGetArray(b,&bold);
215:     for (j=0; j<mvec; j++) {
216:       indx = start+j;
217:       VecSetValues(tmp,1,&indx,bold+j,INSERT_VALUES);
218:     }
219:     VecRestoreArray(b,&bold);
220:     VecDestroy(&b);
221:     VecAssemblyBegin(tmp);
222:     VecAssemblyEnd(tmp);
223:     b    = tmp;
224:   }

226:   MatCreateVecs(A,&x,NULL);
227:   VecDuplicate(b,&u);
228:   if (initialguessfile) {
229:     PetscViewer viewer2;
230:     PetscViewerBinaryOpen(PETSC_COMM_WORLD,initialguessfilename,FILE_MODE_READ,&viewer2);
231:     VecLoad(x,viewer2);
232:     PetscViewerDestroy(&viewer2);
233:     initialguess = PETSC_TRUE;
234:   } else if (initialguess) {
235:     VecSet(x,1.0);
236:   } else {
237:     VecSet(x,0.0);
238:   }


241:   /* Check scaling in A */
242:   flg  = PETSC_FALSE;
243:   PetscOptionsGetBool(NULL,NULL, "-check_scaling", &flg,NULL);
244:   if (flg) {
245:     Vec       max, min;
246:     PetscInt  idx;
247:     PetscReal val;

249:     VecDuplicate(x, &max);
250:     VecDuplicate(x, &min);
251:     MatGetRowMaxAbs(A, max, NULL);
252:     MatGetRowMinAbs(A, min, NULL);
253:     {
254:       PetscViewer viewer;

256:       PetscViewerASCIIOpen(PETSC_COMM_WORLD, "max.data", &viewer);
257:       VecView(max, viewer);
258:       PetscViewerDestroy(&viewer);
259:       PetscViewerASCIIOpen(PETSC_COMM_WORLD, "min.data", &viewer);
260:       VecView(min, viewer);
261:       PetscViewerDestroy(&viewer);
262:     }
263:     VecView(max, PETSC_VIEWER_DRAW_WORLD);
264:     VecMax(max, &idx, &val);
265:     PetscPrintf(PETSC_COMM_WORLD, "Largest max row element %g at row %D\n", (double)val, idx);
266:     VecView(min, PETSC_VIEWER_DRAW_WORLD);
267:     VecMin(min, &idx, &val);
268:     PetscPrintf(PETSC_COMM_WORLD, "Smallest min row element %g at row %D\n", (double)val, idx);
269:     VecMin(max, &idx, &val);
270:     PetscPrintf(PETSC_COMM_WORLD, "Smallest max row element %g at row %D\n", (double)val, idx);
271:     VecPointwiseDivide(max, max, min);
272:     VecMax(max, &idx, &val);
273:     PetscPrintf(PETSC_COMM_WORLD, "Largest row ratio %g at row %D\n", (double)val, idx);
274:     VecView(max, PETSC_VIEWER_DRAW_WORLD);
275:     VecDestroy(&max);
276:     VecDestroy(&min);
277:   }

279:   /*  MatView(A,PETSC_VIEWER_STDOUT_WORLD); */
280:   /* - - - - - - - - - - - New Stage - - - - - - - - - - - - -
281:                     Setup solve for system
282:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
283:   /*
284:      Conclude profiling last stage; begin profiling next stage.
285:   */
286:   PetscPreLoadStage("KSPSetUpSolve");

288:   /*
289:      Create linear solver; set operators; set runtime options.
290:   */
291:   KSPCreate(PETSC_COMM_WORLD,&ksp);
292:   KSPSetInitialGuessNonzero(ksp,initialguess);
293:   num_numfac = 1;
294:   PetscOptionsGetInt(NULL,NULL,"-num_numfac",&num_numfac,NULL);
295:   while (num_numfac--) {
296:     PetscBool lsqr;
297:     char      str[32];
298:     PetscOptionsGetString(NULL,NULL,"-ksp_type",str,32,&lsqr);
299:     if (lsqr) {
300:       PetscStrcmp("lsqr",str,&lsqr);
301:     }
302:     if (lsqr) {
303:       Mat BtB;
304:       MatTransposeMatMult(A,A,MAT_INITIAL_MATRIX,4,&BtB);
305:       KSPSetOperators(ksp,A,BtB);
306:       MatDestroy(&BtB);
307:     } else {
308:       KSPSetOperators(ksp,A,A);
309:     }
310:     KSPSetFromOptions(ksp);

312:     /*
313:      Here we explicitly call KSPSetUp() and KSPSetUpOnBlocks() to
314:      enable more precise profiling of setting up the preconditioner.
315:      These calls are optional, since both will be called within
316:      KSPSolve() if they haven't been called already.
317:     */
318:     KSPSetUp(ksp);
319:     KSPSetUpOnBlocks(ksp);

321:     /* - - - - - - - - - - - New Stage - - - - - - - - - - - - -
322:                          Solve system
323:       - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

325:     /*
326:      Solve linear system; 
327:     */
328:     if (trans) {
329:       KSPSolveTranspose(ksp,b,x);
330:       KSPGetIterationNumber(ksp,&its);
331:     } else {
332:       PetscInt num_rhs=1;
333:       PetscOptionsGetInt(NULL,NULL,"-num_rhs",&num_rhs,NULL);
334:       cknorm = PETSC_FALSE;
335:       PetscOptionsGetBool(NULL,NULL,"-cknorm",&cknorm,NULL);
336:       while (num_rhs--) {
337:         if (num_rhs == 1) VecSet(x,0.0);
338:         KSPSolve(ksp,b,x);
339:       }
340:       KSPGetIterationNumber(ksp,&its);
341:       if (cknorm) {     /* Check error for each rhs */
342:         if (trans) {
343:           MatMultTranspose(A,x,u);
344:         } else {
345:           MatMult(A,x,u);
346:         }
347:         VecAXPY(u,-1.0,b);
348:         VecNorm(u,NORM_2,&norm);
349:         PetscPrintf(PETSC_COMM_WORLD,"  Number of iterations = %3D\n",its);
350:         if (!PetscIsNanScalar(norm)) {
351:           if (norm < 1.e-12) {
352:             PetscPrintf(PETSC_COMM_WORLD,"  Residual norm < 1.e-12\n");
353:           } else {
354:             PetscPrintf(PETSC_COMM_WORLD,"  Residual norm %g\n",(double)norm);
355:           }
356:         }
357:       }
358:     }   /* while (num_rhs--) */

360:     /* - - - - - - - - - - - New Stage - - - - - - - - - - - - -
361:           Check error, print output, free data structures.
362:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

364:     /*
365:        Check error
366:     */
367:     if (trans) {
368:       MatMultTranspose(A,x,u);
369:     } else {
370:       MatMult(A,x,u);
371:     }
372:     VecAXPY(u,-1.0,b);
373:     VecNorm(u,NORM_2,&norm);
374:     /*
375:      Write output (optinally using table for solver details).
376:       - PetscPrintf() handles output for multiprocessor jobs
377:         by printing from only one processor in the communicator.
378:       - KSPView() prints information about the linear solver.
379:     */
380:     if (table) {
381:       char        *matrixname,kspinfo[120];
382:       PetscViewer viewer;

384:       /*
385:        Open a string viewer; then write info to it.
386:       */
387:       PetscViewerStringOpen(PETSC_COMM_WORLD,kspinfo,120,&viewer);
388:       KSPView(ksp,viewer);
389:       PetscStrrchr(file[PetscPreLoadIt],'/',&matrixname);
390:       PetscPrintf(PETSC_COMM_WORLD,"%-8.8s %3D %2.0e %s \n",matrixname,its,norm,kspinfo);

392:       /*
393:         Destroy the viewer
394:       */
395:       PetscViewerDestroy(&viewer);
396:     } else {
397:       PetscPrintf(PETSC_COMM_WORLD,"Number of iterations = %3D\n",its);
398:       if (!PetscIsNanScalar(norm)) {
399:         if (norm < 1.e-12 && !PetscIsNanScalar((PetscScalar)norm)) {
400:           PetscPrintf(PETSC_COMM_WORLD,"  Residual norm < 1.e-12\n");
401:         } else {
402:           PetscPrintf(PETSC_COMM_WORLD,"Residual norm %g\n",(double)norm);
403:         }
404:       }
405:     }
406:     PetscOptionsGetString(NULL,NULL,"-solution",file[3],PETSC_MAX_PATH_LEN,&flg);
407:     if (flg) {
408:       PetscViewer viewer;
409:       Vec         xstar;
410:       PetscReal   norm;

412:       PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[3],FILE_MODE_READ,&viewer);
413:       VecCreate(PETSC_COMM_WORLD,&xstar);
414:       VecLoad(xstar,viewer);
415:       VecAXPY(xstar, -1.0, x);
416:       VecNorm(xstar, NORM_2, &norm);
417:       PetscPrintf(PETSC_COMM_WORLD, "Error norm %g\n", (double)norm);
418:       VecDestroy(&xstar);
419:       PetscViewerDestroy(&viewer);
420:     }
421:     if (outputSoln) {
422:       PetscViewer viewer;

424:       PetscViewerBinaryOpen(PETSC_COMM_WORLD,"solution.petsc",FILE_MODE_WRITE,&viewer);
425:       VecView(x, viewer);
426:       PetscViewerDestroy(&viewer);
427:     }

429:     flg  = PETSC_FALSE;
430:     PetscOptionsGetBool(NULL,NULL, "-ksp_reason", &flg,NULL);
431:     if (flg) {
432:       KSPConvergedReason reason;
433:       KSPGetConvergedReason(ksp,&reason);
434:       PetscPrintf(PETSC_COMM_WORLD,"KSPConvergedReason: %D\n", reason);
435:     }

437:   }   /* while (num_numfac--) */

439:   /*
440:      Free work space.  All PETSc objects should be destroyed when they
441:      are no longer needed.
442:   */
443:   MatDestroy(&A); VecDestroy(&b);
444:   VecDestroy(&u); VecDestroy(&x);
445:   KSPDestroy(&ksp);
446:   PetscPreLoadEnd();
447:   /* -----------------------------------------------------------
448:                       End of linear solver loop
449:      ----------------------------------------------------------- */

451:   PetscFinalize();
452:   return ierr;
453: }


456: /*TEST
457:    
458:    build:
459:       requires: !complex 

461:    testset:
462:       suffix: 1
463:       nsize: 2
464:       args: -f0 ${wPETSC_DIR}/share/petsc/datafiles/matrices/spd-real-int@PETSC_INDEX_SIZE@-float@PETSC_SCALAR_SIZE@

466:    testset:
467:       suffix: 1a
468:       args: -f0 ${wPETSC_DIR}/share/petsc/datafiles/matrices/spd-real-int@PETSC_INDEX_SIZE@-float@PETSC_SCALAR_SIZE@

470:    testset:
471:       nsize: 2
472:       requires: datafilespath double !define(PETSC_USE_64BIT_INDICES)
473:       args: -f0 ${DATAFILESPATH}/matrices/medium
474:       args:  -ksp_type bicg
475:       test:
476:          suffix: 2
477:       test:
478:          suffix: 3
479:          args: -pc_type asm 
480:       
481:    testset:
482:       requires: datafilespath double !define(PETSC_USE_64BIT_INDICES)
483:       args: -f0 ${DATAFILESPATH}/matrices/medium
484:       args: -ksp_type bicg
485:       test:
486:          suffix: 4
487:          args: -pc_type lu 
488:       test:
489:          suffix: 5
490:    
491:    testset:
492:       suffix: 6
493:       requires: datafilespath double !define(PETSC_USE_64BIT_INDICES)
494:       args: -f0 ${DATAFILESPATH}/matrices/fem1
495:       args: -pc_factor_levels 2 -pc_factor_fill 1.73 -ksp_gmres_cgs_refinement_type refine_always
496:    
497:    testset:
498:       suffix: 7
499:       requires: datafilespath double !define(PETSC_USE_64BIT_INDICES)
500:       args: -f0 ${DATAFILESPATH}/matrices/medium 
501:       args: -viewer_binary_skip_info -mat_type seqbaij 
502:       args: -matload_block_size {{2 3 4 5 6 7 8}separate output}
503:       args: -ksp_max_it 100 -ksp_gmres_cgs_refinement_type refine_always 
504:       args: -ksp_rtol 1.0e-15 -ksp_monitor_short
505:       test:
506:          suffix: a
507:       test:
508:          suffix: b
509:          args: -pc_factor_mat_ordering_type nd
510:       test:
511:          suffix: c
512:          args: -pc_factor_levels 1

514:    testset:
515:       suffix: 7_d
516:       requires: datafilespath double !define(PETSC_USE_64BIT_INDICES)
517:       args: -f0 ${DATAFILESPATH}/matrices/medium 
518:       args: -viewer_binary_skip_info -mat_type seqbaij 
519:       args: -matload_block_size {{2 3 4 5 6 7 8}shared output}
520:       args: -ksp_type preonly -pc_type lu
521:    
522:    testset:
523:       suffix: 8
524:       requires: datafilespath double !define(PETSC_USE_64BIT_INDICES)
525:       args: -f0 ${DATAFILESPATH}/matrices/medium  
526:       args: -ksp_diagonal_scale -pc_type eisenstat -ksp_monitor_short -ksp_diagonal_scale_fix -ksp_gmres_cgs_refinement_type refine_always -mat_no_inode
527:    
528:    testset:
529:       suffix: 9
530:       requires: datafilespath double !define(PETSC_USE_64BIT_INDICES)
531:       args: -f0 ${DATAFILESPATH}/matrices/medium 
532:       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
533:       test:
534:          suffix: a
535:          args: -mat_type seqbaij
536:       test:
537:          suffix: b
538:          args: -mat_type seqbaij -trans
539:       test:
540:          suffix: c
541:          nsize: 2
542:          args: -mat_type mpibaij 
543:       test:
544:          suffix: d
545:          nsize: 2
546:          args: -mat_type mpibaij -trans
547:       test:
548:          suffix: e
549:          nsize: 3
550:          args: -mat_type mpibaij 
551:       test:
552:          suffix: f
553:          nsize: 3
554:          args: -mat_type mpibaij -trans
555:    

557:    testset:
558:       suffix: 10
559:       nsize: 2
560:       requires: datafilespath double !define(PETSC_USE_64BIT_INDICES)
561:       args: -ksp_type fgmres -pc_type ksp -f0 ${DATAFILESPATH}/matrices/medium -ksp_fgmres_modifypcksp -ksp_monitor_short
562:    
563:    testset:
564:       TODO: Need to determine goal of this test
565:       suffix: 11
566:       nsize: 2
567:       args: -f0 http://ftp.mcs.anl.gov/pub/petsc/matrices/testmatrix.gz
568:    
569:    testset:
570:       suffix: 12
571:       requires: matlab
572:       args: -pc_type lu -pc_factor_mat_solver_package matlab -f0 ${DATAFILESPATH}/matrices/arco1
573:    
574:    testset:
575:       suffix: 13
576:       requires: lusol
577:       args: -f0 ${DATAFILESPATH}/matrices/arco1
578:       args: -mat_type lusol -pc_type lu
579:    
580:    testset:
581:       nsize: 3
582:       args: -f0 ${DATAFILESPATH}/matrices/medium
583:       requires: datafilespath double !define(PETSC_USE_64BIT_INDICES)
584:       test:
585:          suffix: 14
586:          requires: spai 
587:          args: -pc_type spai 
588:       test:
589:          suffix: 15
590:          requires: hypre
591:          args: -pc_type hypre -pc_hypre_type pilut 
592:       test:
593:          suffix: 16
594:          requires: hypre
595:          args: -pc_type hypre -pc_hypre_type parasails
596:       test:
597:          suffix: 17
598:          requires: hypre
599:          args: -pc_type hypre -pc_hypre_type boomeramg 
600:    
601:    testset:
602:       suffix: 19
603:       requires: datafilespath double !define(PETSC_USE_64BIT_INDICES)
604:       args: -f0 ${DATAFILESPATH}/matrices/poisson1
605:       args: -ksp_type cg -pc_type icc 
606:       args: -pc_factor_levels {{0 2 4}separate output}
607:       test:
608:       test:
609:          args: -mat_type seqsbaij
610:    

612:    testset:
613:       suffix: ILU
614:       requires: datafilespath double !define(PETSC_USE_64BIT_INDICES)
615:       args: -f0 ${DATAFILESPATH}/matrices/small
616:       args: -pc_factor_levels 1
617:       test:
618:       test:
619:          # This is tested against regular ILU (used to be denoted ILUBAIJ)
620:          args: -mat_type baij
621:    
622:    testset:
623:       suffix: aijcusparse
624:       requires: datafilespath double !define(PETSC_USE_64BIT_INDICES) cusparse
625:       args: -f0 ${DATAFILESPATH}/matrices/medium -ksp_monitor_short -ksp_view -mat_view ascii::ascii_info -mat_type aijcusparse -pc_factor_mat_solver_package cusparse -pc_type ilu
626:    
627:    testset:
628:       TODO: No output file. Need to determine if deprecated
629:       suffix: asm_viennacl
630:       nsize: 2
631:       requires: viennacl
632:       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}
633:    
634:    testset:
635:       nsize: 2
636:       requires: datafilespath double !define(PETSC_USE_64BIT_INDICES) hypre
637:       args: -f0 ${DATAFILESPATH}/matrices/poisson2.gz -ksp_monitor_short -ksp_rtol 1.E-9 -pc_type hypre -pc_hypre_type boomeramg 
638:       test:
639:          suffix: boomeramg_euclid
640:          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 
641:          TODO: Need to determine if deprecated
642:       test:
643:          suffix: boomeramg_euclid_bj
644:          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 
645:          TODO: Need to determine if deprecated
646:       test:
647:          suffix: boomeramg_parasails
648:          args: -pc_hypre_boomeramg_smooth_type ParaSails -pc_hypre_boomeramg_smooth_num_levels 2 
649:       test:
650:          suffix: boomeramg_pilut
651:          args: -pc_hypre_boomeramg_smooth_type Pilut -pc_hypre_boomeramg_smooth_num_levels 2 
652:       test:
653:          suffix: boomeramg_schwarz
654:          args: -pc_hypre_boomeramg_smooth_type Schwarz-smoothers 
655:       
656:    testset:
657:       suffix: cg_singlereduction
658:       requires: datafilespath double !define(PETSC_USE_64BIT_INDICES)
659:       args: -f0 ${DATAFILESPATH}/matrices/small
660:       args: -mat_type mpisbaij -ksp_type cg -pc_type eisenstat -ksp_monitor_short -ksp_converged_reason
661:       test:
662:       test:
663:          args: -ksp_cg_single_reduction
664:       
665:    testset:
666:       requires: datafilespath double !define(PETSC_USE_64BIT_INDICES)
667:       args: -f0 ${DATAFILESPATH}/matrices/poisson2.gz
668:       args: -ksp_monitor_short -pc_type icc
669:       test:
670:          suffix: cr
671:          args: -ksp_type cr
672:       test:
673:          suffix: lcd
674:          args: -ksp_type lcd
675:    
676:    testset:
677:       requires: datafilespath double !define(PETSC_USE_64BIT_INDICES)
678:       args: -f0 ${DATAFILESPATH}/matrices/small
679:       args: -ksp_monitor_short -ksp_view -mat_view ascii::ascii_info 
680:       test:
681:          suffix: seqaijcrl
682:          args: -mat_type seqaijcrl
683:       test:
684:          suffix: seqaijperm
685:          args: -mat_type seqaijperm

687:    testset:
688:       nsize: 2
689:       requires: datafilespath double !define(PETSC_USE_64BIT_INDICES)
690:       args: -f0 ${DATAFILESPATH}/matrices/small
691:       args: -ksp_monitor_short -ksp_view 
692:       # Different output files
693:       test:
694:          suffix: mpiaijcrl
695:          args: -mat_type mpiaijcrl
696:       test:
697:          suffix: mpiaijperm
698:          args: -mat_type mpiaijperm
699:    
700:    testset:
701:       nsize: 8
702:       requires: datafilespath double !define(PETSC_USE_64BIT_INDICES)
703:       args: -ksp_monitor_short -ksp_view
704:       test:
705:          suffix: xxt
706:          args: -f0 ${DATAFILESPATH}/matrices/poisson1 -check_symmetry -ksp_type cg -pc_type tfs
707:       test:
708:          suffix: xyt
709:          args: -f0 ${DATAFILESPATH}/matrices/arco1 -ksp_type gmres -pc_type tfs

711:    testset:
712:       # The output file here is the same as mumps
713:       suffix: mumps_cholesky
714:       output_file: output/ex10_mumps.out
715:       requires: datafilespath double !define(PETSC_USE_64BIT_INDICES) mumps
716:       args: -f0 ${DATAFILESPATH}/matrices/small -ksp_type preonly -pc_type cholesky -pc_factor_mat_solver_package mumps -num_numfac 2 -num_rhs 2
717:       nsize: {{1 2}}
718:       test:
719:          args: -mat_type sbaij -mat_ignore_lower_triangular
720:       test:
721:          args: -mat_type aij
722:       test:
723:          args: -mat_type aij -matload_spd
724:       
725:    testset:
726:       # The output file here is the same as mumps
727:       suffix: mumps_lu
728:       output_file: output/ex10_mumps.out
729:       requires: datafilespath double !define(PETSC_USE_64BIT_INDICES) mumps
730:       args: -f0 ${DATAFILESPATH}/matrices/small -ksp_type preonly -pc_type lu -pc_factor_mat_solver_package mumps -num_numfac 2 -num_rhs 2
731:       test:
732:          args: -mat_type seqaij
733:       test:
734:          nsize: 2
735:          args: -mat_type mpiaij
736:       test:
737:          args: -mat_type seqbaij -matload_block_size 2
738:       test:
739:          nsize: 2
740:          args: -mat_type mpibaij -matload_block_size 2
741:       test:
742:          args: -mat_type aij -mat_mumps_icntl_7 5
743:          TODO: Need to determine if deprecated
744:       test:
745:          nsize: 2
746:          args: -mat_type mpiaij -mat_mumps_icntl_28 2 -mat_mumps_icntl_29 2
747:          TODO: Need to determine if deprecated
748:       
749:    testset:
750:       # The output file here is the same as mumps
751:       suffix: mumps_redundant
752:       output_file: output/ex10_mumps_redundant.out
753:       nsize: 8
754:       requires: datafilespath double !define(PETSC_USE_64BIT_INDICES) mumps
755:       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_package mumps -num_numfac 2 -num_rhs 2
756:    
757:    testset:
758:       suffix: pastix_cholesky
759:       requires: datafilespath double !define(PETSC_USE_64BIT_INDICES) pastix
760:       output_file: output/ex10_mumps.out
761:       nsize: {{1 2}}
762:       args: -f0 ${DATAFILESPATH}/matrices/small -ksp_type preonly -pc_factor_mat_solver_package pastix -num_numfac 2 -num_rhs 2 -pc_type cholesky -mat_type sbaij -mat_ignore_lower_triangular
763:       
764:    testset:
765:       suffix: pastix_lu
766:       requires: datafilespath double !define(PETSC_USE_64BIT_INDICES) pastix
767:       args: -f0 ${DATAFILESPATH}/matrices/small -ksp_type preonly -pc_type lu -pc_factor_mat_solver_package pastix -num_numfac 2 -num_rhs 2
768:       output_file: output/ex10_mumps.out
769:       test:
770:          args: -mat_type seqaij
771:       test:
772:          nsize: 2
773:          args: -mat_type mpiaij
774:    
775:    testset:
776:       suffix: pastix_redundant
777:       output_file: output/ex10_mumps_redundant.out
778:       nsize: 8
779:       requires: datafilespath double !define(PETSC_USE_64BIT_INDICES) pastix
780:       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_package pastix -num_numfac 2 -num_rhs 2
781:    
782:      
783:    testset:
784:       suffix: superlu_dist_lu
785:       requires: datafilespath double !define(PETSC_USE_64BIT_INDICES) superlu_dist
786:       output_file: output/ex10_mumps.out
787:       args: -f0 ${DATAFILESPATH}/matrices/small -ksp_type preonly -pc_type lu -pc_factor_mat_solver_package superlu_dist -num_numfac 2 -num_rhs 2
788:       nsize: {{1 2}}
789:       
790:    testset:
791:       suffix: superlu_dist_redundant
792:       nsize: 8
793:       output_file: output/ex10_mumps_redundant.out
794:       requires: datafilespath double !define(PETSC_USE_64BIT_INDICES) superlu
795:       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_package superlu_dist -num_numfac 2 -num_rhs 2
796:    
797:    testset:
798:       suffix: superlu_lu
799:       output_file: output/ex10_mumps.out
800:       requires: datafilespath double !define(PETSC_USE_64BIT_INDICES) superlu
801:       args: -f0 ${DATAFILESPATH}/matrices/small -ksp_type preonly -pc_type lu -pc_factor_mat_solver_package superlu -num_numfac 2 -num_rhs 2
802:    
803:    testset:
804:       suffix: umfpack
805:       requires: datafilespath double !define(PETSC_USE_64BIT_INDICES) suitesparse
806:       args: -f0 ${DATAFILESPATH}/matrices/small -ksp_type preonly -pc_type lu -mat_type seqaij -pc_factor_mat_solver_package umfpack -num_numfac 2 -num_rhs 2
807:    
808:      
809:    testset:
810:       suffix: zeropivot
811:       requires: datafilespath double !define(PETSC_USE_64BIT_INDICES) mumps
812:       args: -f0 ${DATAFILESPATH}/matrices/small -test_zeropivot -ksp_converged_reason -ksp_type fgmres -pc_type ksp
813:       test:
814:          nsize: 3
815:          args: -ksp_pc_type bjacobi
816:       test:
817:          nsize: 2
818:          args: -ksp_ksp_type cg -ksp_pc_type bjacobi -ksp_pc_bjacobi_blocks 1
819:       #test:
820:          #nsize: 3
821:          #args: -ksp_ksp_converged_reason -ksp_pc_type bjacobi -ksp_sub_ksp_converged_reason
822:          #TODO: Need to determine if deprecated
823: TEST*/