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