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