xref: /petsc/src/mat/tests/ex125.c (revision 0338c94495b54376573e2dd674640a2b674ea157)
1c4762a1bSJed Brown static char help[] = "Tests MatSolve() and MatMatSolve() (interface to superlu_dist, mumps and mkl_pardiso).\n\
278bc9606SHong Zhang Example: mpiexec -n <np> ./ex125 -f <matrix binary file> -nrhs 4 -mat_solver_type <>\n\n";
378bc9606SHong Zhang 
478bc9606SHong Zhang /*
578bc9606SHong Zhang -mat_solver_type:
678bc9606SHong Zhang   superlu
778bc9606SHong Zhang   superlu_dist
878bc9606SHong Zhang   mumps
978bc9606SHong Zhang   mkl_pardiso
1078bc9606SHong Zhang   cusparse
1178bc9606SHong Zhang   petsc
1278bc9606SHong Zhang */
13c4762a1bSJed Brown 
14c4762a1bSJed Brown #include <petscmat.h>
15c4762a1bSJed Brown 
165d955bbbSStefano Zampini PetscErrorCode CreateRandom(PetscInt n, PetscInt m, Mat *A)
175d955bbbSStefano Zampini {
185d955bbbSStefano Zampini   PetscFunctionBeginUser;
195d955bbbSStefano Zampini   PetscCall(MatCreate(PETSC_COMM_WORLD, A));
205d955bbbSStefano Zampini   PetscCall(MatSetType(*A, MATAIJ));
215d955bbbSStefano Zampini   PetscCall(MatSetFromOptions(*A));
225d955bbbSStefano Zampini   PetscCall(MatSetSizes(*A, PETSC_DECIDE, PETSC_DECIDE, n, m));
235d955bbbSStefano Zampini   PetscCall(MatSeqAIJSetPreallocation(*A, 5, NULL));
245d955bbbSStefano Zampini   PetscCall(MatMPIAIJSetPreallocation(*A, 5, NULL, 5, NULL));
255d955bbbSStefano Zampini   PetscCall(MatSetRandom(*A, NULL));
265d955bbbSStefano Zampini   PetscCall(MatAssemblyBegin(*A, MAT_FINAL_ASSEMBLY));
275d955bbbSStefano Zampini   PetscCall(MatAssemblyEnd(*A, MAT_FINAL_ASSEMBLY));
285d955bbbSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
295d955bbbSStefano Zampini }
305d955bbbSStefano Zampini 
315d955bbbSStefano Zampini PetscErrorCode CreateIdentity(PetscInt n, Mat *A)
325d955bbbSStefano Zampini {
335d955bbbSStefano Zampini   PetscFunctionBeginUser;
345d955bbbSStefano Zampini   PetscCall(MatCreate(PETSC_COMM_WORLD, A));
355d955bbbSStefano Zampini   PetscCall(MatSetType(*A, MATAIJ));
365d955bbbSStefano Zampini   PetscCall(MatSetFromOptions(*A));
375d955bbbSStefano Zampini   PetscCall(MatSetSizes(*A, PETSC_DECIDE, PETSC_DECIDE, n, n));
385d955bbbSStefano Zampini   PetscCall(MatSetUp(*A));
395d955bbbSStefano Zampini   PetscCall(MatAssemblyBegin(*A, MAT_FINAL_ASSEMBLY));
405d955bbbSStefano Zampini   PetscCall(MatAssemblyEnd(*A, MAT_FINAL_ASSEMBLY));
415d955bbbSStefano Zampini   PetscCall(MatShift(*A, 1.0));
425d955bbbSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
435d955bbbSStefano Zampini }
445d955bbbSStefano Zampini 
45d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
46d71ae5a4SJacob Faibussowitsch {
475d955bbbSStefano Zampini   Mat           A, Ae, RHS = NULL, RHS1 = NULL, C, F, X;
48c4762a1bSJed Brown   Vec           u, x, b;
49c4762a1bSJed Brown   PetscMPIInt   size;
5078bc9606SHong Zhang   PetscInt      m, n, nfact, nsolve, nrhs, ipack = 5;
51070736c7SStefano Zampini   PetscReal     norm, tol = 10 * PETSC_SQRT_MACHINE_EPSILON;
5262671d91SStefano Zampini   IS            perm = NULL, iperm = NULL;
53c4762a1bSJed Brown   MatFactorInfo info;
54c4762a1bSJed Brown   PetscRandom   rand;
5578bc9606SHong Zhang   PetscBool     flg, symm, testMatSolve = PETSC_TRUE, testMatMatSolve = PETSC_TRUE, testMatMatSolveTranspose = PETSC_TRUE, testMatSolveTranspose = PETSC_TRUE, match = PETSC_FALSE;
56c4762a1bSJed Brown   PetscBool     chol = PETSC_FALSE, view = PETSC_FALSE, matsolvexx = PETSC_FALSE;
57c4762a1bSJed Brown #if defined(PETSC_HAVE_MUMPS)
58c4762a1bSJed Brown   PetscBool test_mumps_opts = PETSC_FALSE;
59c4762a1bSJed Brown #endif
60c4762a1bSJed Brown   PetscViewer fd;                       /* viewer */
61c4762a1bSJed Brown   char        file[PETSC_MAX_PATH_LEN]; /* input file name */
6278bc9606SHong Zhang   char        pack[PETSC_MAX_PATH_LEN];
63c4762a1bSJed Brown 
64327415f7SBarry Smith   PetscFunctionBeginUser;
65c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &args, NULL, help));
669566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   /* Determine file from which we read the matrix A */
699566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, NULL, "-f", file, sizeof(file), &flg));
709a14fc28SStefano Zampini   if (flg) { /* Load matrix A */
719566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &fd));
729566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
739566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(A));
749566063dSJacob Faibussowitsch     PetscCall(MatLoad(A, fd));
759566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&fd));
769a14fc28SStefano Zampini   } else {
779a14fc28SStefano Zampini     n = 13;
789566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL));
799566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
809566063dSJacob Faibussowitsch     PetscCall(MatSetType(A, MATAIJ));
819566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(A));
829566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, n, n));
839566063dSJacob Faibussowitsch     PetscCall(MatSetUp(A));
849566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
859566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
869566063dSJacob Faibussowitsch     PetscCall(MatShift(A, 1.0));
879a14fc28SStefano Zampini   }
88c4762a1bSJed Brown 
89c4762a1bSJed Brown   /* if A is symmetric, set its flag -- required by MatGetInertia() */
9078bc9606SHong Zhang   PetscCall(MatIsSymmetric(A, 0.0, &symm));
9178bc9606SHong Zhang   PetscCall(MatSetOption(A, MAT_SYMMETRIC, symm));
92c4762a1bSJed Brown 
935d955bbbSStefano Zampini   PetscCall(PetscOptionsGetBool(NULL, NULL, "-cholesky", &chol, NULL));
945d955bbbSStefano Zampini 
9562671d91SStefano Zampini   /* test MATNEST support */
9662671d91SStefano Zampini   flg = PETSC_FALSE;
9762671d91SStefano Zampini   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_nest", &flg, NULL));
9862671d91SStefano Zampini   if (flg) {
995d955bbbSStefano Zampini     Mat B;
10062671d91SStefano Zampini 
1015d955bbbSStefano Zampini     flg = PETSC_FALSE;
1025d955bbbSStefano Zampini     PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_nest_bordered", &flg, NULL));
1035d955bbbSStefano Zampini     if (!flg) {
1045d955bbbSStefano Zampini       Mat mats[9] = {NULL, NULL, A, NULL, A, NULL, A, NULL, NULL};
1055d955bbbSStefano Zampini 
1065d955bbbSStefano Zampini       /* Create a nested matrix representing
10762671d91SStefano Zampini               | 0 0 A |
10862671d91SStefano Zampini               | 0 A 0 |
10962671d91SStefano Zampini               | A 0 0 |
11062671d91SStefano Zampini       */
11162671d91SStefano Zampini       PetscCall(MatCreateNest(PETSC_COMM_WORLD, 3, NULL, 3, NULL, mats, &B));
1125d955bbbSStefano Zampini     } else {
1135d955bbbSStefano Zampini       Mat mats[4];
1145d955bbbSStefano Zampini 
1155d955bbbSStefano Zampini       /* Create a nested matrix representing
1165d955bbbSStefano Zampini               | Id  R |
1175d955bbbSStefano Zampini               | R^t A |
1185d955bbbSStefano Zampini       */
1195d955bbbSStefano Zampini       PetscCall(MatGetSize(A, NULL, &n));
1205d955bbbSStefano Zampini       m = n + 12;
1215d955bbbSStefano Zampini       PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL));
1225d955bbbSStefano Zampini       PetscCall(CreateIdentity(m, &mats[0]));
1235d955bbbSStefano Zampini       PetscCall(CreateRandom(m, n, &mats[1]));
1245d955bbbSStefano Zampini       mats[3] = A;
1255d955bbbSStefano Zampini 
1265d955bbbSStefano Zampini       /* use CreateTranspose/CreateHermitianTranspose or explicit matrix for debugging purposes */
1275d955bbbSStefano Zampini       flg = PETSC_FALSE;
1285d955bbbSStefano Zampini       PetscCall(PetscOptionsGetBool(NULL, NULL, "-expl", &flg, NULL));
1295d955bbbSStefano Zampini       if (PetscDefined(USE_COMPLEX)) {
1305d955bbbSStefano Zampini         if (chol) { /* Hermitian transpose not supported by MUMPS Cholesky factor */
1315d955bbbSStefano Zampini           if (!flg) PetscCall(MatCreateTranspose(mats[1], &mats[2]));
1325d955bbbSStefano Zampini           else PetscCall(MatTranspose(mats[1], MAT_INITIAL_MATRIX, &mats[2]));
1335d955bbbSStefano Zampini         } else {
1345d955bbbSStefano Zampini           if (!flg) PetscCall(MatCreateHermitianTranspose(mats[1], &mats[2]));
1355d955bbbSStefano Zampini           else PetscCall(MatHermitianTranspose(mats[1], MAT_INITIAL_MATRIX, &mats[2]));
1365d955bbbSStefano Zampini         }
1375d955bbbSStefano Zampini       } else {
1385d955bbbSStefano Zampini         if (!flg) PetscCall(MatCreateTranspose(mats[1], &mats[2]));
1395d955bbbSStefano Zampini         else PetscCall(MatTranspose(mats[1], MAT_INITIAL_MATRIX, &mats[2]));
1405d955bbbSStefano Zampini       }
1415d955bbbSStefano Zampini       PetscCall(MatCreateNest(PETSC_COMM_WORLD, 2, NULL, 2, NULL, mats, &B));
1425d955bbbSStefano Zampini       PetscCall(MatDestroy(&mats[0]));
1435d955bbbSStefano Zampini       PetscCall(MatDestroy(&mats[1]));
1445d955bbbSStefano Zampini       PetscCall(MatDestroy(&mats[2]));
1455d955bbbSStefano Zampini     }
14662671d91SStefano Zampini     PetscCall(MatDestroy(&A));
14762671d91SStefano Zampini     A = B;
14862671d91SStefano Zampini     PetscCall(MatSetOption(A, MAT_SYMMETRIC, symm));
1495d955bbbSStefano Zampini 
1505d955bbbSStefano Zampini     /* not all the combinations of MatMat operations are supported by MATNEST. */
1515d955bbbSStefano Zampini     PetscCall(MatComputeOperator(A, MATAIJ, &Ae));
1525d955bbbSStefano Zampini   } else {
1535d955bbbSStefano Zampini     PetscCall(PetscObjectReference((PetscObject)A));
1545d955bbbSStefano Zampini     Ae = A;
15562671d91SStefano Zampini   }
15662671d91SStefano Zampini   PetscCall(MatGetLocalSize(A, &m, &n));
15762671d91SStefano Zampini   PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "This example is not intended for rectangular matrices (%" PetscInt_FMT ", %" PetscInt_FMT ")", m, n);
15862671d91SStefano Zampini 
1599566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(A, NULL, "-A_view"));
1605d955bbbSStefano Zampini   PetscCall(MatViewFromOptions(Ae, NULL, "-A_view_expl"));
161c4762a1bSJed Brown 
162a5b23f4aSJose E. Roman   /* Create dense matrix C and X; C holds true solution with identical columns */
163c4762a1bSJed Brown   nrhs = 2;
1649566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-nrhs", &nrhs, NULL));
1659566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "ex125: nrhs %" PetscInt_FMT "\n", nrhs));
1669566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &C));
1679566063dSJacob Faibussowitsch   PetscCall(MatSetOptionsPrefix(C, "rhs_"));
1689566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(C, m, PETSC_DECIDE, PETSC_DECIDE, nrhs));
1699566063dSJacob Faibussowitsch   PetscCall(MatSetType(C, MATDENSE));
1709566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(C));
1719566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
172c4762a1bSJed Brown 
1739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-view_factor", &view, NULL));
1749566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_matmatsolve", &testMatMatSolve, NULL));
17578bc9606SHong Zhang   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_matmatsolvetranspose", &testMatMatSolveTranspose, NULL));
17678bc9606SHong Zhang   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_matsolvetranspose", &testMatSolveTranspose, NULL));
177c4762a1bSJed Brown #if defined(PETSC_HAVE_MUMPS)
1789566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_mumps_opts", &test_mumps_opts, NULL));
179c4762a1bSJed Brown #endif
180c4762a1bSJed Brown 
1819566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rand));
1829566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rand));
1839566063dSJacob Faibussowitsch   PetscCall(MatSetRandom(C, rand));
1849566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(C, MAT_DO_NOT_COPY_VALUES, &X));
185c4762a1bSJed Brown 
186c4762a1bSJed Brown   /* Create vectors */
1879566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &x, &b));
1889566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &u)); /* save the true solution */
189c4762a1bSJed Brown 
190c4762a1bSJed Brown   /* Test Factorization */
1919566063dSJacob Faibussowitsch   PetscCall(MatGetOrdering(A, MATORDERINGND, &perm, &iperm));
192c4762a1bSJed Brown 
19378bc9606SHong Zhang   PetscCall(PetscOptionsGetString(NULL, NULL, "-mat_solver_type", pack, sizeof(pack), NULL));
194c4762a1bSJed Brown #if defined(PETSC_HAVE_SUPERLU)
19578bc9606SHong Zhang   PetscCall(PetscStrcmp(MATSOLVERSUPERLU, pack, &match));
19678bc9606SHong Zhang   if (match) {
19728b400f6SJacob Faibussowitsch     PetscCheck(!chol, PETSC_COMM_WORLD, PETSC_ERR_SUP, "SuperLU does not provide Cholesky!");
1989566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, " SUPERLU LU:\n"));
1999566063dSJacob Faibussowitsch     PetscCall(MatGetFactor(A, MATSOLVERSUPERLU, MAT_FACTOR_LU, &F));
20078bc9606SHong Zhang     matsolvexx = PETSC_FALSE; /* Test MatMatSolve(F,RHS,RHS), RHS is a dense matrix, need further work */
20178bc9606SHong Zhang     ipack      = 0;
20278bc9606SHong Zhang     goto skipoptions;
20378bc9606SHong Zhang   }
204c4762a1bSJed Brown #endif
205c4762a1bSJed Brown #if defined(PETSC_HAVE_SUPERLU_DIST)
20678bc9606SHong Zhang   PetscCall(PetscStrcmp(MATSOLVERSUPERLU_DIST, pack, &match));
20778bc9606SHong Zhang   if (match) {
20828b400f6SJacob Faibussowitsch     PetscCheck(!chol, PETSC_COMM_WORLD, PETSC_ERR_SUP, "SuperLU does not provide Cholesky!");
2099566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, " SUPERLU_DIST LU:\n"));
2109566063dSJacob Faibussowitsch     PetscCall(MatGetFactor(A, MATSOLVERSUPERLU_DIST, MAT_FACTOR_LU, &F));
211c4762a1bSJed Brown     matsolvexx = PETSC_TRUE;
21278bc9606SHong Zhang     if (symm) { /* A is symmetric */
21378bc9606SHong Zhang       testMatMatSolveTranspose = PETSC_TRUE;
21478bc9606SHong Zhang       testMatSolveTranspose    = PETSC_TRUE;
21578bc9606SHong Zhang     } else { /* superlu_dist does not support solving A^t x = rhs yet */
21678bc9606SHong Zhang       testMatMatSolveTranspose = PETSC_FALSE;
21778bc9606SHong Zhang       testMatSolveTranspose    = PETSC_FALSE;
21878bc9606SHong Zhang     }
21978bc9606SHong Zhang     ipack = 1;
22078bc9606SHong Zhang     goto skipoptions;
22178bc9606SHong Zhang   }
222c4762a1bSJed Brown #endif
223c4762a1bSJed Brown #if defined(PETSC_HAVE_MUMPS)
22478bc9606SHong Zhang   PetscCall(PetscStrcmp(MATSOLVERMUMPS, pack, &match));
22578bc9606SHong Zhang   if (match) {
226c4762a1bSJed Brown     if (chol) {
2279566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, " MUMPS CHOLESKY:\n"));
2289566063dSJacob Faibussowitsch       PetscCall(MatGetFactor(A, MATSOLVERMUMPS, MAT_FACTOR_CHOLESKY, &F));
229c4762a1bSJed Brown     } else {
2309566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, " MUMPS LU:\n"));
2319566063dSJacob Faibussowitsch       PetscCall(MatGetFactor(A, MATSOLVERMUMPS, MAT_FACTOR_LU, &F));
232c4762a1bSJed Brown     }
233c4762a1bSJed Brown     matsolvexx = PETSC_TRUE;
234c4762a1bSJed Brown     if (test_mumps_opts) {
235c4762a1bSJed Brown       /* test mumps options */
236c4762a1bSJed Brown       PetscInt  icntl;
237c4762a1bSJed Brown       PetscReal cntl;
238c4762a1bSJed Brown 
239c4762a1bSJed Brown       icntl = 2; /* sequential matrix ordering */
2409566063dSJacob Faibussowitsch       PetscCall(MatMumpsSetIcntl(F, 7, icntl));
241c4762a1bSJed Brown 
242c4762a1bSJed Brown       cntl = 1.e-6; /* threshold for row pivot detection */
2439566063dSJacob Faibussowitsch       PetscCall(MatMumpsSetIcntl(F, 24, 1));
2449566063dSJacob Faibussowitsch       PetscCall(MatMumpsSetCntl(F, 3, cntl));
245c4762a1bSJed Brown     }
24678bc9606SHong Zhang     ipack = 2;
24778bc9606SHong Zhang     goto skipoptions;
24878bc9606SHong Zhang   }
249c4762a1bSJed Brown #endif
250c4762a1bSJed Brown #if defined(PETSC_HAVE_MKL_PARDISO)
25178bc9606SHong Zhang   PetscCall(PetscStrcmp(MATSOLVERMKL_PARDISO, pack, &match));
25278bc9606SHong Zhang   if (match) {
253c4762a1bSJed Brown     if (chol) {
2549566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, " MKL_PARDISO CHOLESKY:\n"));
2559566063dSJacob Faibussowitsch       PetscCall(MatGetFactor(A, MATSOLVERMKL_PARDISO, MAT_FACTOR_CHOLESKY, &F));
256c4762a1bSJed Brown     } else {
2579566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, " MKL_PARDISO LU:\n"));
2589566063dSJacob Faibussowitsch       PetscCall(MatGetFactor(A, MATSOLVERMKL_PARDISO, MAT_FACTOR_LU, &F));
259c4762a1bSJed Brown     }
26078bc9606SHong Zhang     ipack = 3;
26178bc9606SHong Zhang     goto skipoptions;
26278bc9606SHong Zhang   }
263c4762a1bSJed Brown #endif
26438a8e8c1SStefano Zampini #if defined(PETSC_HAVE_CUDA)
26578bc9606SHong Zhang   PetscCall(PetscStrcmp(MATSOLVERCUSPARSE, pack, &match));
26678bc9606SHong Zhang   if (match) {
26738a8e8c1SStefano Zampini     if (chol) {
2689566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, " CUSPARSE CHOLESKY:\n"));
2699566063dSJacob Faibussowitsch       PetscCall(MatGetFactor(A, MATSOLVERCUSPARSE, MAT_FACTOR_CHOLESKY, &F));
27038a8e8c1SStefano Zampini     } else {
2719566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, " CUSPARSE LU:\n"));
2729566063dSJacob Faibussowitsch       PetscCall(MatGetFactor(A, MATSOLVERCUSPARSE, MAT_FACTOR_LU, &F));
27338a8e8c1SStefano Zampini     }
27478bc9606SHong Zhang     testMatSolveTranspose    = PETSC_FALSE;
27578bc9606SHong Zhang     testMatMatSolveTranspose = PETSC_FALSE;
27678bc9606SHong Zhang     ipack                    = 4;
27778bc9606SHong Zhang     goto skipoptions;
27878bc9606SHong Zhang   }
27938a8e8c1SStefano Zampini #endif
28078bc9606SHong Zhang   /* PETSc */
28178bc9606SHong Zhang   match = PETSC_TRUE;
28278bc9606SHong Zhang   if (match) {
283c4762a1bSJed Brown     if (chol) {
2849566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, " PETSC CHOLESKY:\n"));
2859566063dSJacob Faibussowitsch       PetscCall(MatGetFactor(A, MATSOLVERPETSC, MAT_FACTOR_CHOLESKY, &F));
286c4762a1bSJed Brown     } else {
2879566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, " PETSC LU:\n"));
2889566063dSJacob Faibussowitsch       PetscCall(MatGetFactor(A, MATSOLVERPETSC, MAT_FACTOR_LU, &F));
289c4762a1bSJed Brown     }
290c4762a1bSJed Brown     matsolvexx = PETSC_TRUE;
29178bc9606SHong Zhang     ipack      = 5;
29278bc9606SHong Zhang     goto skipoptions;
293c4762a1bSJed Brown   }
294c4762a1bSJed Brown 
29578bc9606SHong Zhang skipoptions:
2969566063dSJacob Faibussowitsch   PetscCall(MatFactorInfoInitialize(&info));
297c4762a1bSJed Brown   info.fill      = 5.0;
298c4762a1bSJed Brown   info.shifttype = (PetscReal)MAT_SHIFT_NONE;
299c4762a1bSJed Brown   if (chol) {
3009566063dSJacob Faibussowitsch     PetscCall(MatCholeskyFactorSymbolic(F, A, perm, &info));
301c4762a1bSJed Brown   } else {
3029566063dSJacob Faibussowitsch     PetscCall(MatLUFactorSymbolic(F, A, perm, iperm, &info));
303c4762a1bSJed Brown   }
304c4762a1bSJed Brown 
305c4762a1bSJed Brown   for (nfact = 0; nfact < 2; nfact++) {
306c4762a1bSJed Brown     if (chol) {
3079566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, " %" PetscInt_FMT "-the CHOLESKY numfactorization \n", nfact));
3089566063dSJacob Faibussowitsch       PetscCall(MatCholeskyFactorNumeric(F, A, &info));
309c4762a1bSJed Brown     } else {
3109566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, " %" PetscInt_FMT "-the LU numfactorization \n", nfact));
3119566063dSJacob Faibussowitsch       PetscCall(MatLUFactorNumeric(F, A, &info));
312c4762a1bSJed Brown     }
313c4762a1bSJed Brown     if (view) {
3149566063dSJacob Faibussowitsch       PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_INFO));
3159566063dSJacob Faibussowitsch       PetscCall(MatView(F, PETSC_VIEWER_STDOUT_WORLD));
3169566063dSJacob Faibussowitsch       PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
317c4762a1bSJed Brown       view = PETSC_FALSE;
318c4762a1bSJed Brown     }
319c4762a1bSJed Brown 
320c4762a1bSJed Brown #if defined(PETSC_HAVE_SUPERLU_DIST)
321c4762a1bSJed Brown     if (ipack == 1) { /* Test MatSuperluDistGetDiagU()
322c4762a1bSJed Brown        -- input: matrix factor F; output: main diagonal of matrix U on all processes */
323c4762a1bSJed Brown       PetscInt     M;
324c4762a1bSJed Brown       PetscScalar *diag;
325c4762a1bSJed Brown   #if !defined(PETSC_USE_COMPLEX)
326c4762a1bSJed Brown       PetscInt nneg, nzero, npos;
327c4762a1bSJed Brown   #endif
328c4762a1bSJed Brown 
3299566063dSJacob Faibussowitsch       PetscCall(MatGetSize(F, &M, NULL));
3309566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(M, &diag));
3319566063dSJacob Faibussowitsch       PetscCall(MatSuperluDistGetDiagU(F, diag));
3329566063dSJacob Faibussowitsch       PetscCall(PetscFree(diag));
333c4762a1bSJed Brown 
334c4762a1bSJed Brown   #if !defined(PETSC_USE_COMPLEX)
335c4762a1bSJed Brown       /* Test MatGetInertia() */
33678bc9606SHong Zhang       if (symm) { /* A is symmetric */
3379566063dSJacob Faibussowitsch         PetscCall(MatGetInertia(F, &nneg, &nzero, &npos));
3389566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, " MatInertia: nneg: %" PetscInt_FMT ", nzero: %" PetscInt_FMT ", npos: %" PetscInt_FMT "\n", nneg, nzero, npos));
33978bc9606SHong Zhang       }
340c4762a1bSJed Brown   #endif
341c4762a1bSJed Brown     }
342c4762a1bSJed Brown #endif
343c4762a1bSJed Brown 
344d47f36abSHong Zhang #if defined(PETSC_HAVE_MUMPS)
345d47f36abSHong Zhang     /* mumps interface allows repeated call of MatCholeskyFactorSymbolic(), while the succession calls do nothing */
346d47f36abSHong Zhang     if (ipack == 2) {
347d47f36abSHong Zhang       if (chol) {
3489566063dSJacob Faibussowitsch         PetscCall(MatCholeskyFactorSymbolic(F, A, perm, &info));
3499566063dSJacob Faibussowitsch         PetscCall(MatCholeskyFactorNumeric(F, A, &info));
350d47f36abSHong Zhang       } else {
3519566063dSJacob Faibussowitsch         PetscCall(MatLUFactorSymbolic(F, A, perm, iperm, &info));
3529566063dSJacob Faibussowitsch         PetscCall(MatLUFactorNumeric(F, A, &info));
353d47f36abSHong Zhang       }
354d47f36abSHong Zhang     }
355d47f36abSHong Zhang #endif
356d47f36abSHong Zhang 
357b18964edSHong Zhang     /* Test MatMatSolve(), A X = B, where B can be dense or sparse */
358c4762a1bSJed Brown     if (testMatMatSolve) {
359c4762a1bSJed Brown       if (!nfact) {
3605d955bbbSStefano Zampini         PetscCall(MatMatMult(Ae, C, MAT_INITIAL_MATRIX, 2.0, &RHS));
361c4762a1bSJed Brown       } else {
3625d955bbbSStefano Zampini         PetscCall(MatMatMult(Ae, C, MAT_REUSE_MATRIX, 2.0, &RHS));
363c4762a1bSJed Brown       }
364c4762a1bSJed Brown       for (nsolve = 0; nsolve < 2; nsolve++) {
3659566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "   %" PetscInt_FMT "-the MatMatSolve \n", nsolve));
3669566063dSJacob Faibussowitsch         PetscCall(MatMatSolve(F, RHS, X));
367c4762a1bSJed Brown 
368c4762a1bSJed Brown         /* Check the error */
3699566063dSJacob Faibussowitsch         PetscCall(MatAXPY(X, -1.0, C, SAME_NONZERO_PATTERN));
3709566063dSJacob Faibussowitsch         PetscCall(MatNorm(X, NORM_FROBENIUS, &norm));
371b18964edSHong Zhang         if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT "-the MatMatSolve: Norm of error %g, nsolve %" PetscInt_FMT "\n", nsolve, (double)norm, nsolve));
372c4762a1bSJed Brown       }
373b18964edSHong Zhang 
374c4762a1bSJed Brown       if (matsolvexx) {
375c4762a1bSJed Brown         /* Test MatMatSolve(F,RHS,RHS), RHS is a dense matrix */
3769566063dSJacob Faibussowitsch         PetscCall(MatCopy(RHS, X, SAME_NONZERO_PATTERN));
3779566063dSJacob Faibussowitsch         PetscCall(MatMatSolve(F, X, X));
378c4762a1bSJed Brown         /* Check the error */
3799566063dSJacob Faibussowitsch         PetscCall(MatAXPY(X, -1.0, C, SAME_NONZERO_PATTERN));
3809566063dSJacob Faibussowitsch         PetscCall(MatNorm(X, NORM_FROBENIUS, &norm));
381b18964edSHong Zhang         if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatMatSolve(F,RHS,RHS): Norm of error %g\n", (double)norm));
382c4762a1bSJed Brown       }
383c4762a1bSJed Brown 
384c4762a1bSJed Brown       if (ipack == 2 && size == 1) {
385c4762a1bSJed Brown         Mat spRHS, spRHST, RHST;
386c4762a1bSJed Brown 
3879566063dSJacob Faibussowitsch         PetscCall(MatTranspose(RHS, MAT_INITIAL_MATRIX, &RHST));
3889566063dSJacob Faibussowitsch         PetscCall(MatConvert(RHST, MATAIJ, MAT_INITIAL_MATRIX, &spRHST));
3899566063dSJacob Faibussowitsch         PetscCall(MatCreateTranspose(spRHST, &spRHS));
390c4762a1bSJed Brown         for (nsolve = 0; nsolve < 2; nsolve++) {
3919566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_WORLD, "   %" PetscInt_FMT "-the sparse MatMatSolve \n", nsolve));
3929566063dSJacob Faibussowitsch           PetscCall(MatMatSolve(F, spRHS, X));
393c4762a1bSJed Brown 
394c4762a1bSJed Brown           /* Check the error */
3959566063dSJacob Faibussowitsch           PetscCall(MatAXPY(X, -1.0, C, SAME_NONZERO_PATTERN));
3969566063dSJacob Faibussowitsch           PetscCall(MatNorm(X, NORM_FROBENIUS, &norm));
397b18964edSHong Zhang           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));
398c4762a1bSJed Brown         }
399b18964edSHong Zhang         PetscCall(MatDestroy(&spRHST));
400b18964edSHong Zhang         PetscCall(MatDestroy(&spRHS));
401b18964edSHong Zhang         PetscCall(MatDestroy(&RHST));
402b18964edSHong Zhang       }
403b18964edSHong Zhang     }
404b18964edSHong Zhang 
405b18964edSHong Zhang     /* Test testMatMatSolveTranspose(), A^T X = B, where B can be dense or sparse */
406b18964edSHong Zhang     if (testMatMatSolveTranspose) {
407b18964edSHong Zhang       if (!nfact) {
4085d955bbbSStefano Zampini         PetscCall(MatTransposeMatMult(Ae, C, MAT_INITIAL_MATRIX, 2.0, &RHS1));
409b18964edSHong Zhang       } else {
4105d955bbbSStefano Zampini         PetscCall(MatTransposeMatMult(Ae, C, MAT_REUSE_MATRIX, 2.0, &RHS1));
411b18964edSHong Zhang       }
412b18964edSHong Zhang 
413b18964edSHong Zhang       for (nsolve = 0; nsolve < 2; nsolve++) {
41478bc9606SHong Zhang         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "   %" PetscInt_FMT "-the MatMatSolveTranspose\n", nsolve));
415b18964edSHong Zhang         PetscCall(MatMatSolveTranspose(F, RHS1, X));
416b18964edSHong Zhang 
417b18964edSHong Zhang         /* Check the error */
418b18964edSHong Zhang         PetscCall(MatAXPY(X, -1.0, C, SAME_NONZERO_PATTERN));
419b18964edSHong Zhang         PetscCall(MatNorm(X, NORM_FROBENIUS, &norm));
420b18964edSHong Zhang         if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT "-the MatMatSolveTranspose: Norm of error %g, nsolve %" PetscInt_FMT "\n", nsolve, (double)norm, nsolve));
421b18964edSHong Zhang       }
422b18964edSHong Zhang 
423b18964edSHong Zhang       if (ipack == 2 && size == 1) {
424b18964edSHong Zhang         Mat spRHS, spRHST, RHST;
425b18964edSHong Zhang 
426b18964edSHong Zhang         PetscCall(MatTranspose(RHS1, MAT_INITIAL_MATRIX, &RHST));
427b18964edSHong Zhang         PetscCall(MatConvert(RHST, MATAIJ, MAT_INITIAL_MATRIX, &spRHST));
428b18964edSHong Zhang         PetscCall(MatCreateTranspose(spRHST, &spRHS));
429b18964edSHong Zhang         for (nsolve = 0; nsolve < 2; nsolve++) {
430b18964edSHong Zhang           PetscCall(MatMatSolveTranspose(F, spRHS, X));
431b18964edSHong Zhang 
432b18964edSHong Zhang           /* Check the error */
433b18964edSHong Zhang           PetscCall(MatAXPY(X, -1.0, C, SAME_NONZERO_PATTERN));
434b18964edSHong Zhang           PetscCall(MatNorm(X, NORM_FROBENIUS, &norm));
435b18964edSHong Zhang           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));
436c4762a1bSJed Brown         }
4379566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&spRHST));
4389566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&spRHS));
4399566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&RHST));
440c4762a1bSJed Brown       }
441c4762a1bSJed Brown     }
442c4762a1bSJed Brown 
443c4762a1bSJed Brown     /* Test MatSolve() */
444c4762a1bSJed Brown     if (testMatSolve) {
445c4762a1bSJed Brown       for (nsolve = 0; nsolve < 2; nsolve++) {
4469566063dSJacob Faibussowitsch         PetscCall(VecSetRandom(x, rand));
4479566063dSJacob Faibussowitsch         PetscCall(VecCopy(x, u));
4485d955bbbSStefano Zampini         PetscCall(MatMult(Ae, x, b));
449c4762a1bSJed Brown 
4509566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "   %" PetscInt_FMT "-the MatSolve \n", nsolve));
4519566063dSJacob Faibussowitsch         PetscCall(MatSolve(F, b, x));
452c4762a1bSJed Brown 
453c4762a1bSJed Brown         /* Check the error */
4549566063dSJacob Faibussowitsch         PetscCall(VecAXPY(u, -1.0, x)); /* u <- (-1.0)x + u */
4559566063dSJacob Faibussowitsch         PetscCall(VecNorm(u, NORM_2, &norm));
456c4762a1bSJed Brown         if (norm > tol) {
457c4762a1bSJed Brown           PetscReal resi;
4585d955bbbSStefano Zampini           PetscCall(MatMult(Ae, x, u));   /* u = A*x */
4599566063dSJacob Faibussowitsch           PetscCall(VecAXPY(u, -1.0, b)); /* u <- (-1.0)b + u */
4609566063dSJacob Faibussowitsch           PetscCall(VecNorm(u, NORM_2, &resi));
4619566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatSolve: Norm of error %g, resi %g, numfact %" PetscInt_FMT "\n", (double)norm, (double)resi, nfact));
462c4762a1bSJed Brown         }
463c4762a1bSJed Brown       }
464c4762a1bSJed Brown     }
46578bc9606SHong Zhang 
46678bc9606SHong Zhang     /* Test MatSolveTranspose() */
46778bc9606SHong Zhang     if (testMatSolveTranspose) {
46878bc9606SHong Zhang       for (nsolve = 0; nsolve < 2; nsolve++) {
46978bc9606SHong Zhang         PetscCall(VecSetRandom(x, rand));
47078bc9606SHong Zhang         PetscCall(VecCopy(x, u));
4715d955bbbSStefano Zampini         PetscCall(MatMultTranspose(Ae, x, b));
47278bc9606SHong Zhang 
47378bc9606SHong Zhang         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "   %" PetscInt_FMT "-the MatSolveTranspose\n", nsolve));
47478bc9606SHong Zhang         PetscCall(MatSolveTranspose(F, b, x));
47578bc9606SHong Zhang 
47678bc9606SHong Zhang         /* Check the error */
47778bc9606SHong Zhang         PetscCall(VecAXPY(u, -1.0, x)); /* u <- (-1.0)x + u */
47878bc9606SHong Zhang         PetscCall(VecNorm(u, NORM_2, &norm));
47978bc9606SHong Zhang         if (norm > tol) {
48078bc9606SHong Zhang           PetscReal resi;
481070736c7SStefano Zampini           PetscCall(MatMultTranspose(Ae, x, u)); /* u = A*x */
48278bc9606SHong Zhang           PetscCall(VecAXPY(u, -1.0, b));        /* u <- (-1.0)b + u */
48378bc9606SHong Zhang           PetscCall(VecNorm(u, NORM_2, &resi));
48478bc9606SHong Zhang           PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatSolveTranspose: Norm of error %g, resi %g, numfact %" PetscInt_FMT "\n", (double)norm, (double)resi, nfact));
48578bc9606SHong Zhang         }
48678bc9606SHong Zhang       }
48778bc9606SHong Zhang     }
488c4762a1bSJed Brown   }
489c4762a1bSJed Brown 
490c4762a1bSJed Brown   /* Free data structures */
4915d955bbbSStefano Zampini   PetscCall(MatDestroy(&Ae));
4929566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
4939566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
4949566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&F));
4959566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&X));
4969566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&RHS));
497b18964edSHong Zhang   PetscCall(MatDestroy(&RHS1));
498c4762a1bSJed Brown 
4999566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rand));
5009566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&perm));
5019566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iperm));
5029566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
5039566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b));
5049566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
5059566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
506b122ec5aSJacob Faibussowitsch   return 0;
507c4762a1bSJed Brown }
508c4762a1bSJed Brown 
509c4762a1bSJed Brown /*TEST
510c4762a1bSJed Brown 
511c4762a1bSJed Brown    test:
512dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
51378bc9606SHong Zhang       args: -f ${DATAFILESPATH}/matrices/medium -mat_solver_type petsc
514c4762a1bSJed Brown       output_file: output/ex125.out
515c4762a1bSJed Brown 
516c4762a1bSJed Brown    test:
5179a14fc28SStefano Zampini       suffix: 2
51878bc9606SHong Zhang       args: -mat_solver_type petsc
5199a14fc28SStefano Zampini       output_file: output/ex125.out
5209a14fc28SStefano Zampini 
5219a14fc28SStefano Zampini    test:
522c4762a1bSJed Brown       suffix: mkl_pardiso
523dfd57a17SPierre Jolivet       requires: mkl_pardiso datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
52478bc9606SHong Zhang       args: -f ${DATAFILESPATH}/matrices/small -mat_solver_type mkl_pardiso
525c4762a1bSJed Brown 
526c4762a1bSJed Brown    test:
5279a14fc28SStefano Zampini       suffix: mkl_pardiso_2
5289a14fc28SStefano Zampini       requires: mkl_pardiso
52978bc9606SHong Zhang       args: -mat_solver_type mkl_pardiso
5309a14fc28SStefano Zampini       output_file: output/ex125_mkl_pardiso.out
5319a14fc28SStefano Zampini 
5329a14fc28SStefano Zampini    test:
533c4762a1bSJed Brown       suffix: mumps
534dfd57a17SPierre Jolivet       requires: mumps datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
53578bc9606SHong Zhang       args: -f ${DATAFILESPATH}/matrices/small -mat_solver_type mumps
536c4762a1bSJed Brown       output_file: output/ex125_mumps_seq.out
537c4762a1bSJed Brown 
538c4762a1bSJed Brown    test:
53962671d91SStefano Zampini       suffix: mumps_nest
54062671d91SStefano Zampini       requires: mumps datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
5415d955bbbSStefano Zampini       args: -f ${DATAFILESPATH}/matrices/small -mat_solver_type mumps -test_nest -test_nest_bordered {{0 1}}
54262671d91SStefano Zampini       output_file: output/ex125_mumps_seq.out
54362671d91SStefano Zampini 
54462671d91SStefano Zampini    test:
545c4762a1bSJed Brown       suffix: mumps_2
546c4762a1bSJed Brown       nsize: 3
547dfd57a17SPierre Jolivet       requires: mumps datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
54878bc9606SHong Zhang       args: -f ${DATAFILESPATH}/matrices/small -mat_solver_type mumps
549c4762a1bSJed Brown       output_file: output/ex125_mumps_par.out
550c4762a1bSJed Brown 
551c4762a1bSJed Brown    test:
55262671d91SStefano Zampini       suffix: mumps_2_nest
55362671d91SStefano Zampini       nsize: 3
55462671d91SStefano Zampini       requires: mumps datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
5555d955bbbSStefano Zampini       args: -f ${DATAFILESPATH}/matrices/small -mat_solver_type mumps -test_nest -test_nest_bordered {{0 1}}
55662671d91SStefano Zampini       output_file: output/ex125_mumps_par.out
55762671d91SStefano Zampini 
55862671d91SStefano Zampini    test:
5599a14fc28SStefano Zampini       suffix: mumps_3
5609a14fc28SStefano Zampini       requires: mumps
56178bc9606SHong Zhang       args: -mat_solver_type mumps
5629a14fc28SStefano Zampini       output_file: output/ex125_mumps_seq.out
5639a14fc28SStefano Zampini 
5649a14fc28SStefano Zampini    test:
56562671d91SStefano Zampini       suffix: mumps_3_nest
56662671d91SStefano Zampini       requires: mumps
5675d955bbbSStefano Zampini       args: -mat_solver_type mumps -test_nest -test_nest_bordered {{0 1}}
56862671d91SStefano Zampini       output_file: output/ex125_mumps_seq.out
56962671d91SStefano Zampini 
57062671d91SStefano Zampini    test:
5719a14fc28SStefano Zampini       suffix: mumps_4
5729a14fc28SStefano Zampini       nsize: 3
5739a14fc28SStefano Zampini       requires: mumps
57478bc9606SHong Zhang       args: -mat_solver_type mumps
5759a14fc28SStefano Zampini       output_file: output/ex125_mumps_par.out
5769a14fc28SStefano Zampini 
5779a14fc28SStefano Zampini    test:
57862671d91SStefano Zampini       suffix: mumps_4_nest
57962671d91SStefano Zampini       nsize: 3
58062671d91SStefano Zampini       requires: mumps
5815d955bbbSStefano Zampini       args: -mat_solver_type mumps -test_nest -test_nest_bordered {{0 1}}
58262671d91SStefano Zampini       output_file: output/ex125_mumps_par.out
58362671d91SStefano Zampini 
58462671d91SStefano Zampini    test:
585d47f36abSHong Zhang       suffix: mumps_5
586d47f36abSHong Zhang       nsize: 3
587d47f36abSHong Zhang       requires: mumps
58878bc9606SHong Zhang       args: -mat_solver_type mumps -cholesky
589d47f36abSHong Zhang       output_file: output/ex125_mumps_par_cholesky.out
590d47f36abSHong Zhang 
591d47f36abSHong Zhang    test:
59262671d91SStefano Zampini       suffix: mumps_5_nest
59362671d91SStefano Zampini       nsize: 3
59462671d91SStefano Zampini       requires: mumps
5955d955bbbSStefano Zampini       args: -mat_solver_type mumps -cholesky -test_nest -test_nest_bordered {{0 1}}
59662671d91SStefano Zampini       output_file: output/ex125_mumps_par_cholesky.out
59762671d91SStefano Zampini 
59862671d91SStefano Zampini    test:
59978bc9606SHong Zhang       suffix: superlu
60078bc9606SHong Zhang       requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) superlu
60178bc9606SHong Zhang       args: -f ${DATAFILESPATH}/matrices/medium -mat_solver_type superlu
60278bc9606SHong Zhang       output_file: output/ex125_superlu.out
60378bc9606SHong Zhang 
60478bc9606SHong Zhang    test:
605c4762a1bSJed Brown       suffix: superlu_dist
6069a14fc28SStefano Zampini       nsize: {{1 3}}
607d4783600SBarry Smith       requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) superlu_dist
60878bc9606SHong Zhang       args: -f ${DATAFILESPATH}/matrices/small -mat_solver_type superlu_dist -mat_superlu_dist_rowperm NOROWPERM
60978bc9606SHong Zhang       output_file: output/ex125_superlu_dist.out
610c4762a1bSJed Brown 
611c4762a1bSJed Brown    test:
612c4762a1bSJed Brown       suffix: superlu_dist_2
6139a14fc28SStefano Zampini       nsize: {{1 3}}
6149a14fc28SStefano Zampini       requires: superlu_dist !complex
61578bc9606SHong Zhang       args: -n 36 -mat_solver_type superlu_dist -mat_superlu_dist_rowperm NOROWPERM
616c4762a1bSJed Brown       output_file: output/ex125_superlu_dist.out
617c4762a1bSJed Brown 
618c4762a1bSJed Brown    test:
61978bc9606SHong Zhang       suffix: superlu_dist_3
62078bc9606SHong Zhang       nsize: {{1 3}}
62178bc9606SHong Zhang       requires: superlu_dist !complex
62278bc9606SHong Zhang       requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) superlu_dist
62378bc9606SHong Zhang       args: -f ${DATAFILESPATH}/matrices/medium -mat_solver_type superlu_dist -mat_superlu_dist_rowperm NOROWPERM
62478bc9606SHong Zhang       output_file: output/ex125_superlu_dist_nonsymmetric.out
62578bc9606SHong Zhang 
62678bc9606SHong Zhang    test:
627c4762a1bSJed Brown       suffix: superlu_dist_complex
628c4762a1bSJed Brown       nsize: 3
629d4783600SBarry Smith       requires: datafilespath double superlu_dist complex !defined(PETSC_USE_64BIT_INDICES)
63078bc9606SHong Zhang       args: -f ${DATAFILESPATH}/matrices/farzad_B_rhs -mat_solver_type superlu_dist
631c4762a1bSJed Brown       output_file: output/ex125_superlu_dist_complex.out
632c4762a1bSJed Brown 
63338a8e8c1SStefano Zampini    test:
6349a14fc28SStefano Zampini       suffix: superlu_dist_complex_2
6359a14fc28SStefano Zampini       nsize: 3
6369a14fc28SStefano Zampini       requires: superlu_dist complex
63778bc9606SHong Zhang       args: -mat_solver_type superlu_dist
63878bc9606SHong Zhang       output_file: output/ex125_superlu_dist_complex_2.out
6399a14fc28SStefano Zampini 
6409a14fc28SStefano Zampini    test:
64138a8e8c1SStefano Zampini       suffix: cusparse
642dfd57a17SPierre Jolivet       requires: cuda datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
643*0338c944SBarry Smith       #TODO: fix the bug with cholesky
64478bc9606SHong Zhang       #args: -mat_type aijcusparse -f ${DATAFILESPATH}/matrices/small -mat_solver_type cusparse -cholesky {{0 1}separate output}
64578bc9606SHong Zhang       args: -mat_type aijcusparse -f ${DATAFILESPATH}/matrices/small -mat_solver_type cusparse -cholesky {{0}separate output}
64638a8e8c1SStefano Zampini 
6479a14fc28SStefano Zampini    test:
6489a14fc28SStefano Zampini       suffix: cusparse_2
6499a14fc28SStefano Zampini       requires: cuda
65078bc9606SHong Zhang       args: -mat_type aijcusparse -mat_solver_type cusparse -cholesky {{0 1}separate output}
6519a14fc28SStefano Zampini 
652c4762a1bSJed Brown TEST*/
653