xref: /petsc/src/mat/tests/ex125.c (revision 82cad2ff0f8a4693401f7b4c96d2f3af9a0b620a)
1 static char help[] = "Tests MatSolve() and MatMatSolve() (interface to superlu_dist, mumps and mkl_pardiso).\n\
2 Example: mpiexec -n <np> ./ex125 -f <matrix binary file> -nrhs 4 -mat_solver_type <>\n\n";
3 
4 /*
5 -mat_solver_type:
6   superlu
7   superlu_dist
8   mumps
9   mkl_pardiso
10   cusparse
11   petsc
12 */
13 
14 #include <petscmat.h>
15 
16 int main(int argc, char **args)
17 {
18   Mat           A, RHS = NULL, RHS1 = NULL, C, F, X;
19   Vec           u, x, b;
20   PetscMPIInt   size;
21   PetscInt      m, n, nfact, nsolve, nrhs, ipack = 5;
22   PetscReal     norm, tol = 1.e-10;
23   IS            perm = NULL, iperm = NULL;
24   MatFactorInfo info;
25   PetscRandom   rand;
26   PetscBool     flg, symm, testMatSolve = PETSC_TRUE, testMatMatSolve = PETSC_TRUE, testMatMatSolveTranspose = PETSC_TRUE, testMatSolveTranspose = PETSC_TRUE, match = PETSC_FALSE;
27   PetscBool     chol = PETSC_FALSE, view = PETSC_FALSE, matsolvexx = PETSC_FALSE;
28 #if defined(PETSC_HAVE_MUMPS)
29   PetscBool test_mumps_opts = PETSC_FALSE;
30 #endif
31   PetscViewer fd;                       /* viewer */
32   char        file[PETSC_MAX_PATH_LEN]; /* input file name */
33   char        pack[PETSC_MAX_PATH_LEN];
34 
35   PetscFunctionBeginUser;
36   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
37   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
38 
39   /* Determine file from which we read the matrix A */
40   PetscCall(PetscOptionsGetString(NULL, NULL, "-f", file, sizeof(file), &flg));
41   if (flg) { /* Load matrix A */
42     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &fd));
43     PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
44     PetscCall(MatSetFromOptions(A));
45     PetscCall(MatLoad(A, fd));
46     PetscCall(PetscViewerDestroy(&fd));
47   } else {
48     n = 13;
49     PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL));
50     PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
51     PetscCall(MatSetType(A, MATAIJ));
52     PetscCall(MatSetFromOptions(A));
53     PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, n, n));
54     PetscCall(MatSetUp(A));
55     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
56     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
57     PetscCall(MatShift(A, 1.0));
58   }
59 
60   /* if A is symmetric, set its flag -- required by MatGetInertia() */
61   PetscCall(MatIsSymmetric(A, 0.0, &symm));
62   PetscCall(MatSetOption(A, MAT_SYMMETRIC, symm));
63 
64   /* test MATNEST support */
65   flg = PETSC_FALSE;
66   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_nest", &flg, NULL));
67   if (flg) {
68     /* Create a nested matrix representing
69 
70        | 0 0 A |
71        | 0 A 0 |
72        | A 0 0 |
73 
74     */
75     Mat B, mats[9] = {NULL, NULL, A, NULL, A, NULL, A, NULL, NULL};
76 
77     PetscCall(MatCreateNest(PETSC_COMM_WORLD, 3, NULL, 3, NULL, mats, &B));
78     PetscCall(MatDestroy(&A));
79     A = B;
80     PetscCall(MatSetOption(A, MAT_SYMMETRIC, symm));
81   }
82   PetscCall(MatGetLocalSize(A, &m, &n));
83   PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "This example is not intended for rectangular matrices (%" PetscInt_FMT ", %" PetscInt_FMT ")", m, n);
84 
85   PetscCall(MatViewFromOptions(A, NULL, "-A_view"));
86 
87   /* Create dense matrix C and X; C holds true solution with identical columns */
88   nrhs = 2;
89   PetscCall(PetscOptionsGetInt(NULL, NULL, "-nrhs", &nrhs, NULL));
90   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "ex125: nrhs %" PetscInt_FMT "\n", nrhs));
91   PetscCall(MatCreate(PETSC_COMM_WORLD, &C));
92   PetscCall(MatSetOptionsPrefix(C, "rhs_"));
93   PetscCall(MatSetSizes(C, m, PETSC_DECIDE, PETSC_DECIDE, nrhs));
94   PetscCall(MatSetType(C, MATDENSE));
95   PetscCall(MatSetFromOptions(C));
96   PetscCall(MatSetUp(C));
97 
98   PetscCall(PetscOptionsGetBool(NULL, NULL, "-view_factor", &view, NULL));
99   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_matmatsolve", &testMatMatSolve, NULL));
100   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_matmatsolvetranspose", &testMatMatSolveTranspose, NULL));
101   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_matsolvetranspose", &testMatSolveTranspose, NULL));
102   PetscCall(PetscOptionsGetBool(NULL, NULL, "-cholesky", &chol, NULL));
103 #if defined(PETSC_HAVE_MUMPS)
104   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_mumps_opts", &test_mumps_opts, NULL));
105 #endif
106 
107   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rand));
108   PetscCall(PetscRandomSetFromOptions(rand));
109   PetscCall(MatSetRandom(C, rand));
110   PetscCall(MatDuplicate(C, MAT_DO_NOT_COPY_VALUES, &X));
111 
112   /* Create vectors */
113   PetscCall(MatCreateVecs(A, &x, &b));
114   PetscCall(VecDuplicate(x, &u)); /* save the true solution */
115 
116   /* Test Factorization */
117   PetscCall(MatGetOrdering(A, MATORDERINGND, &perm, &iperm));
118 
119   PetscCall(PetscOptionsGetString(NULL, NULL, "-mat_solver_type", pack, sizeof(pack), NULL));
120 #if defined(PETSC_HAVE_SUPERLU)
121   PetscCall(PetscStrcmp(MATSOLVERSUPERLU, pack, &match));
122   if (match) {
123     PetscCheck(!chol, PETSC_COMM_WORLD, PETSC_ERR_SUP, "SuperLU does not provide Cholesky!");
124     PetscCall(PetscPrintf(PETSC_COMM_WORLD, " SUPERLU LU:\n"));
125     PetscCall(MatGetFactor(A, MATSOLVERSUPERLU, MAT_FACTOR_LU, &F));
126     matsolvexx = PETSC_FALSE; /* Test MatMatSolve(F,RHS,RHS), RHS is a dense matrix, need further work */
127     ipack      = 0;
128     goto skipoptions;
129   }
130 #endif
131 #if defined(PETSC_HAVE_SUPERLU_DIST)
132   PetscCall(PetscStrcmp(MATSOLVERSUPERLU_DIST, pack, &match));
133   if (match) {
134     PetscCheck(!chol, PETSC_COMM_WORLD, PETSC_ERR_SUP, "SuperLU does not provide Cholesky!");
135     PetscCall(PetscPrintf(PETSC_COMM_WORLD, " SUPERLU_DIST LU:\n"));
136     PetscCall(MatGetFactor(A, MATSOLVERSUPERLU_DIST, MAT_FACTOR_LU, &F));
137     matsolvexx = PETSC_TRUE;
138     if (symm) { /* A is symmetric */
139       testMatMatSolveTranspose = PETSC_TRUE;
140       testMatSolveTranspose    = PETSC_TRUE;
141     } else { /* superlu_dist does not support solving A^t x = rhs yet */
142       testMatMatSolveTranspose = PETSC_FALSE;
143       testMatSolveTranspose    = PETSC_FALSE;
144     }
145     ipack = 1;
146     goto skipoptions;
147   }
148 #endif
149 #if defined(PETSC_HAVE_MUMPS)
150   PetscCall(PetscStrcmp(MATSOLVERMUMPS, pack, &match));
151   if (match) {
152     if (chol) {
153       PetscCall(PetscPrintf(PETSC_COMM_WORLD, " MUMPS CHOLESKY:\n"));
154       PetscCall(MatGetFactor(A, MATSOLVERMUMPS, MAT_FACTOR_CHOLESKY, &F));
155     } else {
156       PetscCall(PetscPrintf(PETSC_COMM_WORLD, " MUMPS LU:\n"));
157       PetscCall(MatGetFactor(A, MATSOLVERMUMPS, MAT_FACTOR_LU, &F));
158     }
159     matsolvexx = PETSC_TRUE;
160     if (test_mumps_opts) {
161       /* test mumps options */
162       PetscInt  icntl;
163       PetscReal cntl;
164 
165       icntl = 2; /* sequential matrix ordering */
166       PetscCall(MatMumpsSetIcntl(F, 7, icntl));
167 
168       cntl = 1.e-6; /* threshold for row pivot detection */
169       PetscCall(MatMumpsSetIcntl(F, 24, 1));
170       PetscCall(MatMumpsSetCntl(F, 3, cntl));
171     }
172     ipack = 2;
173     goto skipoptions;
174   }
175 #endif
176 #if defined(PETSC_HAVE_MKL_PARDISO)
177   PetscCall(PetscStrcmp(MATSOLVERMKL_PARDISO, pack, &match));
178   if (match) {
179     if (chol) {
180       PetscCall(PetscPrintf(PETSC_COMM_WORLD, " MKL_PARDISO CHOLESKY:\n"));
181       PetscCall(MatGetFactor(A, MATSOLVERMKL_PARDISO, MAT_FACTOR_CHOLESKY, &F));
182     } else {
183       PetscCall(PetscPrintf(PETSC_COMM_WORLD, " MKL_PARDISO LU:\n"));
184       PetscCall(MatGetFactor(A, MATSOLVERMKL_PARDISO, MAT_FACTOR_LU, &F));
185     }
186     ipack = 3;
187     goto skipoptions;
188   }
189 #endif
190 #if defined(PETSC_HAVE_CUDA)
191   PetscCall(PetscStrcmp(MATSOLVERCUSPARSE, pack, &match));
192   if (match) {
193     if (chol) {
194       PetscCall(PetscPrintf(PETSC_COMM_WORLD, " CUSPARSE CHOLESKY:\n"));
195       PetscCall(MatGetFactor(A, MATSOLVERCUSPARSE, MAT_FACTOR_CHOLESKY, &F));
196     } else {
197       PetscCall(PetscPrintf(PETSC_COMM_WORLD, " CUSPARSE LU:\n"));
198       PetscCall(MatGetFactor(A, MATSOLVERCUSPARSE, MAT_FACTOR_LU, &F));
199     }
200     testMatSolveTranspose    = PETSC_FALSE;
201     testMatMatSolveTranspose = PETSC_FALSE;
202     ipack                    = 4;
203     goto skipoptions;
204   }
205 #endif
206   /* PETSc */
207   match = PETSC_TRUE;
208   if (match) {
209     if (chol) {
210       PetscCall(PetscPrintf(PETSC_COMM_WORLD, " PETSC CHOLESKY:\n"));
211       PetscCall(MatGetFactor(A, MATSOLVERPETSC, MAT_FACTOR_CHOLESKY, &F));
212     } else {
213       PetscCall(PetscPrintf(PETSC_COMM_WORLD, " PETSC LU:\n"));
214       PetscCall(MatGetFactor(A, MATSOLVERPETSC, MAT_FACTOR_LU, &F));
215     }
216     matsolvexx = PETSC_TRUE;
217     ipack      = 5;
218     goto skipoptions;
219   }
220 
221 skipoptions:
222   PetscCall(MatFactorInfoInitialize(&info));
223   info.fill      = 5.0;
224   info.shifttype = (PetscReal)MAT_SHIFT_NONE;
225   if (chol) {
226     PetscCall(MatCholeskyFactorSymbolic(F, A, perm, &info));
227   } else {
228     PetscCall(MatLUFactorSymbolic(F, A, perm, iperm, &info));
229   }
230 
231   for (nfact = 0; nfact < 2; nfact++) {
232     if (chol) {
233       PetscCall(PetscPrintf(PETSC_COMM_WORLD, " %" PetscInt_FMT "-the CHOLESKY numfactorization \n", nfact));
234       PetscCall(MatCholeskyFactorNumeric(F, A, &info));
235     } else {
236       PetscCall(PetscPrintf(PETSC_COMM_WORLD, " %" PetscInt_FMT "-the LU numfactorization \n", nfact));
237       PetscCall(MatLUFactorNumeric(F, A, &info));
238     }
239     if (view) {
240       PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_INFO));
241       PetscCall(MatView(F, PETSC_VIEWER_STDOUT_WORLD));
242       PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
243       view = PETSC_FALSE;
244     }
245 
246 #if defined(PETSC_HAVE_SUPERLU_DIST)
247     if (ipack == 1) { /* Test MatSuperluDistGetDiagU()
248        -- input: matrix factor F; output: main diagonal of matrix U on all processes */
249       PetscInt     M;
250       PetscScalar *diag;
251   #if !defined(PETSC_USE_COMPLEX)
252       PetscInt nneg, nzero, npos;
253   #endif
254 
255       PetscCall(MatGetSize(F, &M, NULL));
256       PetscCall(PetscMalloc1(M, &diag));
257       PetscCall(MatSuperluDistGetDiagU(F, diag));
258       PetscCall(PetscFree(diag));
259 
260   #if !defined(PETSC_USE_COMPLEX)
261       /* Test MatGetInertia() */
262       if (symm) { /* A is symmetric */
263         PetscCall(MatGetInertia(F, &nneg, &nzero, &npos));
264         PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, " MatInertia: nneg: %" PetscInt_FMT ", nzero: %" PetscInt_FMT ", npos: %" PetscInt_FMT "\n", nneg, nzero, npos));
265       }
266   #endif
267     }
268 #endif
269 
270 #if defined(PETSC_HAVE_MUMPS)
271     /* mumps interface allows repeated call of MatCholeskyFactorSymbolic(), while the succession calls do nothing */
272     if (ipack == 2) {
273       if (chol) {
274         PetscCall(MatCholeskyFactorSymbolic(F, A, perm, &info));
275         PetscCall(MatCholeskyFactorNumeric(F, A, &info));
276       } else {
277         PetscCall(MatLUFactorSymbolic(F, A, perm, iperm, &info));
278         PetscCall(MatLUFactorNumeric(F, A, &info));
279       }
280     }
281 #endif
282 
283     /* Test MatMatSolve(), A X = B, where B can be dense or sparse */
284     if (testMatMatSolve) {
285       if (!nfact) {
286         PetscCall(MatMatMult(A, C, MAT_INITIAL_MATRIX, 2.0, &RHS));
287       } else {
288         PetscCall(MatMatMult(A, C, MAT_REUSE_MATRIX, 2.0, &RHS));
289       }
290       for (nsolve = 0; nsolve < 2; nsolve++) {
291         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "   %" PetscInt_FMT "-the MatMatSolve \n", nsolve));
292         PetscCall(MatMatSolve(F, RHS, X));
293 
294         /* Check the error */
295         PetscCall(MatAXPY(X, -1.0, C, SAME_NONZERO_PATTERN));
296         PetscCall(MatNorm(X, NORM_FROBENIUS, &norm));
297         if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT "-the MatMatSolve: Norm of error %g, nsolve %" PetscInt_FMT "\n", nsolve, (double)norm, nsolve));
298       }
299 
300       if (matsolvexx) {
301         /* Test MatMatSolve(F,RHS,RHS), RHS is a dense matrix */
302         PetscCall(MatCopy(RHS, X, SAME_NONZERO_PATTERN));
303         PetscCall(MatMatSolve(F, X, X));
304         /* Check the error */
305         PetscCall(MatAXPY(X, -1.0, C, SAME_NONZERO_PATTERN));
306         PetscCall(MatNorm(X, NORM_FROBENIUS, &norm));
307         if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatMatSolve(F,RHS,RHS): Norm of error %g\n", (double)norm));
308       }
309 
310       if (ipack == 2 && size == 1) {
311         Mat spRHS, spRHST, RHST;
312 
313         PetscCall(MatTranspose(RHS, MAT_INITIAL_MATRIX, &RHST));
314         PetscCall(MatConvert(RHST, MATAIJ, MAT_INITIAL_MATRIX, &spRHST));
315         PetscCall(MatCreateTranspose(spRHST, &spRHS));
316         for (nsolve = 0; nsolve < 2; nsolve++) {
317           PetscCall(PetscPrintf(PETSC_COMM_WORLD, "   %" PetscInt_FMT "-the sparse MatMatSolve \n", nsolve));
318           PetscCall(MatMatSolve(F, spRHS, X));
319 
320           /* Check the error */
321           PetscCall(MatAXPY(X, -1.0, C, SAME_NONZERO_PATTERN));
322           PetscCall(MatNorm(X, NORM_FROBENIUS, &norm));
323           if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT "-the sparse MatMatSolve: Norm of error %g, nsolve %" PetscInt_FMT "\n", nsolve, (double)norm, nsolve));
324         }
325         PetscCall(MatDestroy(&spRHST));
326         PetscCall(MatDestroy(&spRHS));
327         PetscCall(MatDestroy(&RHST));
328       }
329     }
330 
331     /* Test testMatMatSolveTranspose(), A^T X = B, where B can be dense or sparse */
332     if (testMatMatSolveTranspose) {
333       if (!nfact) {
334         PetscCall(MatTransposeMatMult(A, C, MAT_INITIAL_MATRIX, 2.0, &RHS1));
335       } else {
336         PetscCall(MatTransposeMatMult(A, C, MAT_REUSE_MATRIX, 2.0, &RHS1));
337       }
338 
339       for (nsolve = 0; nsolve < 2; nsolve++) {
340         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "   %" PetscInt_FMT "-the MatMatSolveTranspose\n", nsolve));
341         PetscCall(MatMatSolveTranspose(F, RHS1, X));
342 
343         /* Check the error */
344         PetscCall(MatAXPY(X, -1.0, C, SAME_NONZERO_PATTERN));
345         PetscCall(MatNorm(X, NORM_FROBENIUS, &norm));
346         if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT "-the MatMatSolveTranspose: Norm of error %g, nsolve %" PetscInt_FMT "\n", nsolve, (double)norm, nsolve));
347       }
348 
349       if (ipack == 2 && size == 1) {
350         Mat spRHS, spRHST, RHST;
351 
352         PetscCall(MatTranspose(RHS1, MAT_INITIAL_MATRIX, &RHST));
353         PetscCall(MatConvert(RHST, MATAIJ, MAT_INITIAL_MATRIX, &spRHST));
354         PetscCall(MatCreateTranspose(spRHST, &spRHS));
355         for (nsolve = 0; nsolve < 2; nsolve++) {
356           PetscCall(MatMatSolveTranspose(F, spRHS, X));
357 
358           /* Check the error */
359           PetscCall(MatAXPY(X, -1.0, C, SAME_NONZERO_PATTERN));
360           PetscCall(MatNorm(X, NORM_FROBENIUS, &norm));
361           if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT "-the sparse MatMatSolveTranspose: Norm of error %g, nsolve %" PetscInt_FMT "\n", nsolve, (double)norm, nsolve));
362         }
363         PetscCall(MatDestroy(&spRHST));
364         PetscCall(MatDestroy(&spRHS));
365         PetscCall(MatDestroy(&RHST));
366       }
367     }
368 
369     /* Test MatSolve() */
370     if (testMatSolve) {
371       for (nsolve = 0; nsolve < 2; nsolve++) {
372         PetscCall(VecSetRandom(x, rand));
373         PetscCall(VecCopy(x, u));
374         PetscCall(MatMult(A, x, b));
375 
376         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "   %" PetscInt_FMT "-the MatSolve \n", nsolve));
377         PetscCall(MatSolve(F, b, x));
378 
379         /* Check the error */
380         PetscCall(VecAXPY(u, -1.0, x)); /* u <- (-1.0)x + u */
381         PetscCall(VecNorm(u, NORM_2, &norm));
382         if (norm > tol) {
383           PetscReal resi;
384           PetscCall(MatMult(A, x, u));    /* u = A*x */
385           PetscCall(VecAXPY(u, -1.0, b)); /* u <- (-1.0)b + u */
386           PetscCall(VecNorm(u, NORM_2, &resi));
387           PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatSolve: Norm of error %g, resi %g, numfact %" PetscInt_FMT "\n", (double)norm, (double)resi, nfact));
388         }
389       }
390     }
391 
392     /* Test MatSolveTranspose() */
393     if (testMatSolveTranspose) {
394       for (nsolve = 0; nsolve < 2; nsolve++) {
395         PetscCall(VecSetRandom(x, rand));
396         PetscCall(VecCopy(x, u));
397         PetscCall(MatMultTranspose(A, x, b));
398 
399         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "   %" PetscInt_FMT "-the MatSolveTranspose\n", nsolve));
400         PetscCall(MatSolveTranspose(F, b, x));
401 
402         /* Check the error */
403         PetscCall(VecAXPY(u, -1.0, x)); /* u <- (-1.0)x + u */
404         PetscCall(VecNorm(u, NORM_2, &norm));
405         if (norm > tol) {
406           PetscReal resi;
407           PetscCall(VecAXPY(u, -1.0, b)); /* u <- (-1.0)b + u */
408           PetscCall(VecNorm(u, NORM_2, &resi));
409           PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatSolveTranspose: Norm of error %g, resi %g, numfact %" PetscInt_FMT "\n", (double)norm, (double)resi, nfact));
410         }
411       }
412     }
413   }
414 
415   /* Free data structures */
416   PetscCall(MatDestroy(&A));
417   PetscCall(MatDestroy(&C));
418   PetscCall(MatDestroy(&F));
419   PetscCall(MatDestroy(&X));
420   PetscCall(MatDestroy(&RHS));
421   PetscCall(MatDestroy(&RHS1));
422 
423   PetscCall(PetscRandomDestroy(&rand));
424   PetscCall(ISDestroy(&perm));
425   PetscCall(ISDestroy(&iperm));
426   PetscCall(VecDestroy(&x));
427   PetscCall(VecDestroy(&b));
428   PetscCall(VecDestroy(&u));
429   PetscCall(PetscFinalize());
430   return 0;
431 }
432 
433 /*TEST
434 
435    test:
436       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
437       args: -f ${DATAFILESPATH}/matrices/medium -mat_solver_type petsc
438       output_file: output/ex125.out
439 
440    test:
441       suffix: 2
442       args: -mat_solver_type petsc
443       output_file: output/ex125.out
444 
445    test:
446       suffix: mkl_pardiso
447       requires: mkl_pardiso datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
448       args: -f ${DATAFILESPATH}/matrices/small -mat_solver_type mkl_pardiso
449 
450    test:
451       suffix: mkl_pardiso_2
452       requires: mkl_pardiso
453       args: -mat_solver_type mkl_pardiso
454       output_file: output/ex125_mkl_pardiso.out
455 
456    test:
457       suffix: mumps
458       requires: mumps datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
459       args: -f ${DATAFILESPATH}/matrices/small -mat_solver_type mumps
460       output_file: output/ex125_mumps_seq.out
461 
462    test:
463       suffix: mumps_nest
464       requires: mumps datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
465       args: -f ${DATAFILESPATH}/matrices/small -mat_solver_type mumps -test_nest
466       output_file: output/ex125_mumps_seq.out
467 
468    test:
469       suffix: mumps_2
470       nsize: 3
471       requires: mumps datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
472       args: -f ${DATAFILESPATH}/matrices/small -mat_solver_type mumps
473       output_file: output/ex125_mumps_par.out
474 
475    test:
476       suffix: mumps_2_nest
477       nsize: 3
478       requires: mumps datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
479       args: -f ${DATAFILESPATH}/matrices/small -mat_solver_type mumps -test_nest
480       output_file: output/ex125_mumps_par.out
481 
482    test:
483       suffix: mumps_3
484       requires: mumps
485       args: -mat_solver_type mumps
486       output_file: output/ex125_mumps_seq.out
487 
488    test:
489       suffix: mumps_3_nest
490       requires: mumps
491       args: -mat_solver_type mumps -test_nest
492       output_file: output/ex125_mumps_seq.out
493 
494    test:
495       suffix: mumps_4
496       nsize: 3
497       requires: mumps
498       args: -mat_solver_type mumps
499       output_file: output/ex125_mumps_par.out
500 
501    test:
502       suffix: mumps_4_nest
503       nsize: 3
504       requires: mumps
505       args: -mat_solver_type mumps -test_nest
506       output_file: output/ex125_mumps_par.out
507 
508    test:
509       suffix: mumps_5
510       nsize: 3
511       requires: mumps
512       args: -mat_solver_type mumps -cholesky
513       output_file: output/ex125_mumps_par_cholesky.out
514 
515    test:
516       suffix: mumps_5_nest
517       nsize: 3
518       requires: mumps
519       args: -mat_solver_type mumps -cholesky -test_nest
520       output_file: output/ex125_mumps_par_cholesky.out
521 
522    test:
523       suffix: superlu
524       requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) superlu
525       args: -f ${DATAFILESPATH}/matrices/medium -mat_solver_type superlu
526       output_file: output/ex125_superlu.out
527 
528    test:
529       suffix: superlu_dist
530       nsize: {{1 3}}
531       requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) superlu_dist
532       args: -f ${DATAFILESPATH}/matrices/small -mat_solver_type superlu_dist -mat_superlu_dist_rowperm NOROWPERM
533       output_file: output/ex125_superlu_dist.out
534 
535    test:
536       suffix: superlu_dist_2
537       nsize: {{1 3}}
538       requires: superlu_dist !complex
539       args: -n 36 -mat_solver_type superlu_dist -mat_superlu_dist_rowperm NOROWPERM
540       output_file: output/ex125_superlu_dist.out
541 
542    test:
543       suffix: superlu_dist_3
544       nsize: {{1 3}}
545       requires: superlu_dist !complex
546       requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) superlu_dist
547       args: -f ${DATAFILESPATH}/matrices/medium -mat_solver_type superlu_dist -mat_superlu_dist_rowperm NOROWPERM
548       output_file: output/ex125_superlu_dist_nonsymmetric.out
549 
550    test:
551       suffix: superlu_dist_complex
552       nsize: 3
553       requires: datafilespath double superlu_dist complex !defined(PETSC_USE_64BIT_INDICES)
554       args: -f ${DATAFILESPATH}/matrices/farzad_B_rhs -mat_solver_type superlu_dist
555       output_file: output/ex125_superlu_dist_complex.out
556 
557    test:
558       suffix: superlu_dist_complex_2
559       nsize: 3
560       requires: superlu_dist complex
561       args: -mat_solver_type superlu_dist
562       output_file: output/ex125_superlu_dist_complex_2.out
563 
564    test:
565       suffix: cusparse
566       requires: cuda datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
567       #todo: fix the bug with cholesky
568       #args: -mat_type aijcusparse -f ${DATAFILESPATH}/matrices/small -mat_solver_type cusparse -cholesky {{0 1}separate output}
569       args: -mat_type aijcusparse -f ${DATAFILESPATH}/matrices/small -mat_solver_type cusparse -cholesky {{0}separate output}
570 
571    test:
572       suffix: cusparse_2
573       requires: cuda
574       args: -mat_type aijcusparse -mat_solver_type cusparse -cholesky {{0 1}separate output}
575 
576 TEST*/
577