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; 56*39989c4fSStefano Zampini PetscBool chol = PETSC_FALSE, view = PETSC_FALSE, matsolvexx = PETSC_FALSE, test_inertia; 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 93*39989c4fSStefano Zampini test_inertia = symm; 94*39989c4fSStefano Zampini PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_inertia", &test_inertia, NULL)); 95*39989c4fSStefano Zampini 965d955bbbSStefano Zampini PetscCall(PetscOptionsGetBool(NULL, NULL, "-cholesky", &chol, NULL)); 975d955bbbSStefano Zampini 9862671d91SStefano Zampini /* test MATNEST support */ 9962671d91SStefano Zampini flg = PETSC_FALSE; 10062671d91SStefano Zampini PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_nest", &flg, NULL)); 10162671d91SStefano Zampini if (flg) { 1025d955bbbSStefano Zampini Mat B; 10362671d91SStefano Zampini 1045d955bbbSStefano Zampini flg = PETSC_FALSE; 1055d955bbbSStefano Zampini PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_nest_bordered", &flg, NULL)); 1065d955bbbSStefano Zampini if (!flg) { 1075d955bbbSStefano Zampini Mat mats[9] = {NULL, NULL, A, NULL, A, NULL, A, NULL, NULL}; 1085d955bbbSStefano Zampini 1095d955bbbSStefano Zampini /* Create a nested matrix representing 11062671d91SStefano Zampini | 0 0 A | 11162671d91SStefano Zampini | 0 A 0 | 11262671d91SStefano Zampini | A 0 0 | 11362671d91SStefano Zampini */ 11462671d91SStefano Zampini PetscCall(MatCreateNest(PETSC_COMM_WORLD, 3, NULL, 3, NULL, mats, &B)); 11553587d93SPierre Jolivet flg = PETSC_TRUE; 1165d955bbbSStefano Zampini } else { 1175d955bbbSStefano Zampini Mat mats[4]; 1185d955bbbSStefano Zampini 1195d955bbbSStefano Zampini /* Create a nested matrix representing 1205d955bbbSStefano Zampini | Id R | 1215d955bbbSStefano Zampini | R^t A | 1225d955bbbSStefano Zampini */ 1235d955bbbSStefano Zampini PetscCall(MatGetSize(A, NULL, &n)); 1245d955bbbSStefano Zampini m = n + 12; 1255d955bbbSStefano Zampini PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL)); 1265d955bbbSStefano Zampini PetscCall(CreateIdentity(m, &mats[0])); 1275d955bbbSStefano Zampini PetscCall(CreateRandom(m, n, &mats[1])); 1285d955bbbSStefano Zampini mats[3] = A; 1295d955bbbSStefano Zampini 1305d955bbbSStefano Zampini /* use CreateTranspose/CreateHermitianTranspose or explicit matrix for debugging purposes */ 1315d955bbbSStefano Zampini flg = PETSC_FALSE; 1325d955bbbSStefano Zampini PetscCall(PetscOptionsGetBool(NULL, NULL, "-expl", &flg, NULL)); 13353587d93SPierre Jolivet #if PetscDefined(USE_COMPLEX) 1345d955bbbSStefano Zampini if (chol) { /* Hermitian transpose not supported by MUMPS Cholesky factor */ 1355d955bbbSStefano Zampini if (!flg) PetscCall(MatCreateTranspose(mats[1], &mats[2])); 1365d955bbbSStefano Zampini else PetscCall(MatTranspose(mats[1], MAT_INITIAL_MATRIX, &mats[2])); 13753587d93SPierre Jolivet flg = PETSC_TRUE; 1385d955bbbSStefano Zampini } else { 13953587d93SPierre Jolivet if (!flg) { 14053587d93SPierre Jolivet Mat B; 14153587d93SPierre Jolivet 14253587d93SPierre Jolivet PetscCall(MatDuplicate(mats[1], MAT_COPY_VALUES, &B)); 14353587d93SPierre Jolivet PetscCall(MatCreateHermitianTranspose(B, &mats[2])); 14453587d93SPierre Jolivet PetscCall(MatDestroy(&B)); 14553587d93SPierre Jolivet if (n == m) { 14653587d93SPierre Jolivet PetscCall(MatScale(mats[2], PetscCMPLX(4.0, -2.0))); 14753587d93SPierre Jolivet PetscCall(MatShift(mats[2], PetscCMPLX(-2.0, 1.0))); // mats[2] = (4 - 2i) B* - (2 - i) I 14853587d93SPierre Jolivet PetscCall(MatCreateHermitianTranspose(mats[2], &B)); 14953587d93SPierre Jolivet PetscCall(MatDestroy(mats + 2)); 15053587d93SPierre Jolivet PetscCall(MatScale(B, 0.5)); 15153587d93SPierre Jolivet PetscCall(MatShift(B, PetscCMPLX(1.0, 0.5))); // B = 0.5 mats[2]* - (1 - 0.5i) I = (2 + i) B - (1 + 0.5i) I + (1 + 0.5i) I = (2 + i) B 15253587d93SPierre Jolivet PetscCall(MatCreateHermitianTranspose(B, &mats[2])); // mats[2] = B* = (2 - i) B* 15353587d93SPierre Jolivet PetscCall(MatDestroy(&B)); 15453587d93SPierre Jolivet PetscCall(MatScale(mats[1], PetscCMPLX(2.0, 1.0))); // mats[1] = (2 + i) B = mats[2]* 15553587d93SPierre Jolivet } else flg = PETSC_TRUE; 15653587d93SPierre Jolivet } else PetscCall(MatHermitianTranspose(mats[1], MAT_INITIAL_MATRIX, &mats[2])); 1575d955bbbSStefano Zampini } 15853587d93SPierre Jolivet #else 15953587d93SPierre Jolivet if (!flg) { 16053587d93SPierre Jolivet Mat B; 16153587d93SPierre Jolivet 16253587d93SPierre Jolivet PetscCall(MatDuplicate(mats[1], MAT_COPY_VALUES, &B)); 16353587d93SPierre Jolivet PetscCall(MatCreateTranspose(B, &mats[2])); 16453587d93SPierre Jolivet PetscCall(MatDestroy(&B)); 16553587d93SPierre Jolivet if (n == m) { 16653587d93SPierre Jolivet PetscCall(MatScale(mats[2], 4.0)); 16753587d93SPierre Jolivet PetscCall(MatShift(mats[2], -2.0)); // mats[2] = 4 B' - 2 I 16853587d93SPierre Jolivet PetscCall(MatCreateTranspose(mats[2], &B)); 16953587d93SPierre Jolivet PetscCall(MatDestroy(mats + 2)); 17053587d93SPierre Jolivet PetscCall(MatScale(B, 0.5)); 17153587d93SPierre Jolivet PetscCall(MatShift(B, 1.0)); // B = 0.5 mats[2]' + I = 0.5 (4 B' - 2 I)' + I = 2 B 17253587d93SPierre Jolivet PetscCall(MatCreateTranspose(B, &mats[2])); // mats[2] = B' = 2 B' 17353587d93SPierre Jolivet PetscCall(MatDestroy(&B)); 17453587d93SPierre Jolivet PetscCall(MatScale(mats[1], 2.0)); // mats[1] = 2 B = mats[2]' 17553587d93SPierre Jolivet } else flg = PETSC_TRUE; 17653587d93SPierre Jolivet } else PetscCall(MatTranspose(mats[1], MAT_INITIAL_MATRIX, &mats[2])); 17753587d93SPierre Jolivet #endif 1785d955bbbSStefano Zampini PetscCall(MatCreateNest(PETSC_COMM_WORLD, 2, NULL, 2, NULL, mats, &B)); 1795d955bbbSStefano Zampini PetscCall(MatDestroy(&mats[0])); 1805d955bbbSStefano Zampini PetscCall(MatDestroy(&mats[1])); 1815d955bbbSStefano Zampini PetscCall(MatDestroy(&mats[2])); 1825d955bbbSStefano Zampini } 18362671d91SStefano Zampini PetscCall(MatDestroy(&A)); 18462671d91SStefano Zampini A = B; 18562671d91SStefano Zampini PetscCall(MatSetOption(A, MAT_SYMMETRIC, symm)); 1865d955bbbSStefano Zampini 1875d955bbbSStefano Zampini /* not all the combinations of MatMat operations are supported by MATNEST. */ 1885d955bbbSStefano Zampini PetscCall(MatComputeOperator(A, MATAIJ, &Ae)); 1895d955bbbSStefano Zampini } else { 1905d955bbbSStefano Zampini PetscCall(PetscObjectReference((PetscObject)A)); 1915d955bbbSStefano Zampini Ae = A; 19253587d93SPierre Jolivet flg = PETSC_TRUE; 19362671d91SStefano Zampini } 19462671d91SStefano Zampini PetscCall(MatGetLocalSize(A, &m, &n)); 19562671d91SStefano 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); 19662671d91SStefano Zampini 1979566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A, NULL, "-A_view")); 1985d955bbbSStefano Zampini PetscCall(MatViewFromOptions(Ae, NULL, "-A_view_expl")); 199c4762a1bSJed Brown 200a5b23f4aSJose E. Roman /* Create dense matrix C and X; C holds true solution with identical columns */ 201c4762a1bSJed Brown nrhs = 2; 2029566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-nrhs", &nrhs, NULL)); 2039566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "ex125: nrhs %" PetscInt_FMT "\n", nrhs)); 2049566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &C)); 2059566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(C, "rhs_")); 2069566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, m, PETSC_DECIDE, PETSC_DECIDE, nrhs)); 2079566063dSJacob Faibussowitsch PetscCall(MatSetType(C, MATDENSE)); 2089566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C)); 2099566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 210c4762a1bSJed Brown 2119566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-view_factor", &view, NULL)); 2129566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_matmatsolve", &testMatMatSolve, NULL)); 21378bc9606SHong Zhang PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_matmatsolvetranspose", &testMatMatSolveTranspose, NULL)); 21478bc9606SHong Zhang PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_matsolvetranspose", &testMatSolveTranspose, NULL)); 215c4762a1bSJed Brown #if defined(PETSC_HAVE_MUMPS) 2169566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_mumps_opts", &test_mumps_opts, NULL)); 217c4762a1bSJed Brown #endif 218c4762a1bSJed Brown 2199566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rand)); 2209566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rand)); 2219566063dSJacob Faibussowitsch PetscCall(MatSetRandom(C, rand)); 2229566063dSJacob Faibussowitsch PetscCall(MatDuplicate(C, MAT_DO_NOT_COPY_VALUES, &X)); 223c4762a1bSJed Brown 224c4762a1bSJed Brown /* Create vectors */ 2259566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &x, &b)); 2269566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &u)); /* save the true solution */ 227c4762a1bSJed Brown 228c4762a1bSJed Brown /* Test Factorization */ 22953587d93SPierre Jolivet if (flg) PetscCall(MatGetOrdering(A, MATORDERINGND, &perm, &iperm)); // TODO FIXME: MatConvert_Nest_AIJ() does not support chained MatCreate[Hermitian]Transpose() 230c4762a1bSJed Brown 23178bc9606SHong Zhang PetscCall(PetscOptionsGetString(NULL, NULL, "-mat_solver_type", pack, sizeof(pack), NULL)); 232c4762a1bSJed Brown #if defined(PETSC_HAVE_SUPERLU) 23378bc9606SHong Zhang PetscCall(PetscStrcmp(MATSOLVERSUPERLU, pack, &match)); 23478bc9606SHong Zhang if (match) { 23528b400f6SJacob Faibussowitsch PetscCheck(!chol, PETSC_COMM_WORLD, PETSC_ERR_SUP, "SuperLU does not provide Cholesky!"); 2369566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " SUPERLU LU:\n")); 2379566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A, MATSOLVERSUPERLU, MAT_FACTOR_LU, &F)); 23878bc9606SHong Zhang matsolvexx = PETSC_FALSE; /* Test MatMatSolve(F,RHS,RHS), RHS is a dense matrix, need further work */ 23978bc9606SHong Zhang ipack = 0; 24078bc9606SHong Zhang goto skipoptions; 24178bc9606SHong Zhang } 242c4762a1bSJed Brown #endif 243c4762a1bSJed Brown #if defined(PETSC_HAVE_SUPERLU_DIST) 24478bc9606SHong Zhang PetscCall(PetscStrcmp(MATSOLVERSUPERLU_DIST, pack, &match)); 24578bc9606SHong Zhang if (match) { 24628b400f6SJacob Faibussowitsch PetscCheck(!chol, PETSC_COMM_WORLD, PETSC_ERR_SUP, "SuperLU does not provide Cholesky!"); 2479566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " SUPERLU_DIST LU:\n")); 2489566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A, MATSOLVERSUPERLU_DIST, MAT_FACTOR_LU, &F)); 249c4762a1bSJed Brown matsolvexx = PETSC_TRUE; 25078bc9606SHong Zhang if (symm) { /* A is symmetric */ 25178bc9606SHong Zhang testMatMatSolveTranspose = PETSC_TRUE; 25278bc9606SHong Zhang testMatSolveTranspose = PETSC_TRUE; 25378bc9606SHong Zhang } else { /* superlu_dist does not support solving A^t x = rhs yet */ 25478bc9606SHong Zhang testMatMatSolveTranspose = PETSC_FALSE; 25578bc9606SHong Zhang testMatSolveTranspose = PETSC_FALSE; 25678bc9606SHong Zhang } 25778bc9606SHong Zhang ipack = 1; 25878bc9606SHong Zhang goto skipoptions; 25978bc9606SHong Zhang } 260c4762a1bSJed Brown #endif 261c4762a1bSJed Brown #if defined(PETSC_HAVE_MUMPS) 26278bc9606SHong Zhang PetscCall(PetscStrcmp(MATSOLVERMUMPS, pack, &match)); 26378bc9606SHong Zhang if (match) { 264c4762a1bSJed Brown if (chol) { 2659566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " MUMPS CHOLESKY:\n")); 2669566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A, MATSOLVERMUMPS, MAT_FACTOR_CHOLESKY, &F)); 267c4762a1bSJed Brown } else { 2689566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " MUMPS LU:\n")); 2699566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A, MATSOLVERMUMPS, MAT_FACTOR_LU, &F)); 270c4762a1bSJed Brown } 271c4762a1bSJed Brown matsolvexx = PETSC_TRUE; 272c4762a1bSJed Brown if (test_mumps_opts) { 273c4762a1bSJed Brown /* test mumps options */ 274c4762a1bSJed Brown PetscInt icntl; 275c4762a1bSJed Brown PetscReal cntl; 276c4762a1bSJed Brown 277c4762a1bSJed Brown icntl = 2; /* sequential matrix ordering */ 2789566063dSJacob Faibussowitsch PetscCall(MatMumpsSetIcntl(F, 7, icntl)); 279c4762a1bSJed Brown 280c4762a1bSJed Brown cntl = 1.e-6; /* threshold for row pivot detection */ 2819566063dSJacob Faibussowitsch PetscCall(MatMumpsSetIcntl(F, 24, 1)); 2829566063dSJacob Faibussowitsch PetscCall(MatMumpsSetCntl(F, 3, cntl)); 283c4762a1bSJed Brown } 28478bc9606SHong Zhang ipack = 2; 28578bc9606SHong Zhang goto skipoptions; 28678bc9606SHong Zhang } 287c4762a1bSJed Brown #endif 288c4762a1bSJed Brown #if defined(PETSC_HAVE_MKL_PARDISO) 28978bc9606SHong Zhang PetscCall(PetscStrcmp(MATSOLVERMKL_PARDISO, pack, &match)); 29078bc9606SHong Zhang if (match) { 291c4762a1bSJed Brown if (chol) { 2929566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " MKL_PARDISO CHOLESKY:\n")); 2939566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A, MATSOLVERMKL_PARDISO, MAT_FACTOR_CHOLESKY, &F)); 294c4762a1bSJed Brown } else { 2959566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " MKL_PARDISO LU:\n")); 2969566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A, MATSOLVERMKL_PARDISO, MAT_FACTOR_LU, &F)); 297c4762a1bSJed Brown } 29878bc9606SHong Zhang ipack = 3; 29978bc9606SHong Zhang goto skipoptions; 30078bc9606SHong Zhang } 301c4762a1bSJed Brown #endif 30238a8e8c1SStefano Zampini #if defined(PETSC_HAVE_CUDA) 30378bc9606SHong Zhang PetscCall(PetscStrcmp(MATSOLVERCUSPARSE, pack, &match)); 30478bc9606SHong Zhang if (match) { 30538a8e8c1SStefano Zampini if (chol) { 3069566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " CUSPARSE CHOLESKY:\n")); 3079566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A, MATSOLVERCUSPARSE, MAT_FACTOR_CHOLESKY, &F)); 30838a8e8c1SStefano Zampini } else { 3099566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " CUSPARSE LU:\n")); 3109566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A, MATSOLVERCUSPARSE, MAT_FACTOR_LU, &F)); 31138a8e8c1SStefano Zampini } 31278bc9606SHong Zhang testMatSolveTranspose = PETSC_FALSE; 31378bc9606SHong Zhang testMatMatSolveTranspose = PETSC_FALSE; 31478bc9606SHong Zhang ipack = 4; 31578bc9606SHong Zhang goto skipoptions; 31678bc9606SHong Zhang } 31738a8e8c1SStefano Zampini #endif 31878bc9606SHong Zhang /* PETSc */ 31978bc9606SHong Zhang match = PETSC_TRUE; 32078bc9606SHong Zhang if (match) { 321c4762a1bSJed Brown if (chol) { 3229566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " PETSC CHOLESKY:\n")); 3239566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A, MATSOLVERPETSC, MAT_FACTOR_CHOLESKY, &F)); 324c4762a1bSJed Brown } else { 3259566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " PETSC LU:\n")); 3269566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A, MATSOLVERPETSC, MAT_FACTOR_LU, &F)); 327c4762a1bSJed Brown } 328c4762a1bSJed Brown matsolvexx = PETSC_TRUE; 32978bc9606SHong Zhang ipack = 5; 33078bc9606SHong Zhang goto skipoptions; 331c4762a1bSJed Brown } 332c4762a1bSJed Brown 33378bc9606SHong Zhang skipoptions: 3349566063dSJacob Faibussowitsch PetscCall(MatFactorInfoInitialize(&info)); 335c4762a1bSJed Brown info.fill = 5.0; 336c4762a1bSJed Brown info.shifttype = (PetscReal)MAT_SHIFT_NONE; 337c4762a1bSJed Brown if (chol) { 3389566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorSymbolic(F, A, perm, &info)); 339c4762a1bSJed Brown } else { 3409566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(F, A, perm, iperm, &info)); 341c4762a1bSJed Brown } 342c4762a1bSJed Brown 343c4762a1bSJed Brown for (nfact = 0; nfact < 2; nfact++) { 344c4762a1bSJed Brown if (chol) { 3459566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " %" PetscInt_FMT "-the CHOLESKY numfactorization \n", nfact)); 3469566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorNumeric(F, A, &info)); 347c4762a1bSJed Brown } else { 3489566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " %" PetscInt_FMT "-the LU numfactorization \n", nfact)); 3499566063dSJacob Faibussowitsch PetscCall(MatLUFactorNumeric(F, A, &info)); 350c4762a1bSJed Brown } 351c4762a1bSJed Brown if (view) { 3529566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_INFO)); 3539566063dSJacob Faibussowitsch PetscCall(MatView(F, PETSC_VIEWER_STDOUT_WORLD)); 3549566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD)); 355c4762a1bSJed Brown view = PETSC_FALSE; 356c4762a1bSJed Brown } 357c4762a1bSJed Brown 358c4762a1bSJed Brown #if defined(PETSC_HAVE_SUPERLU_DIST) 359c4762a1bSJed Brown if (ipack == 1) { /* Test MatSuperluDistGetDiagU() 360c4762a1bSJed Brown -- input: matrix factor F; output: main diagonal of matrix U on all processes */ 361c4762a1bSJed Brown PetscInt M; 362c4762a1bSJed Brown PetscScalar *diag; 363c4762a1bSJed Brown #if !defined(PETSC_USE_COMPLEX) 364c4762a1bSJed Brown PetscInt nneg, nzero, npos; 365c4762a1bSJed Brown #endif 366c4762a1bSJed Brown 3679566063dSJacob Faibussowitsch PetscCall(MatGetSize(F, &M, NULL)); 3689566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(M, &diag)); 3699566063dSJacob Faibussowitsch PetscCall(MatSuperluDistGetDiagU(F, diag)); 3709566063dSJacob Faibussowitsch PetscCall(PetscFree(diag)); 371c4762a1bSJed Brown 372c4762a1bSJed Brown #if !defined(PETSC_USE_COMPLEX) 373c4762a1bSJed Brown /* Test MatGetInertia() */ 374*39989c4fSStefano Zampini if (test_inertia) { /* A is symmetric */ 3759566063dSJacob Faibussowitsch PetscCall(MatGetInertia(F, &nneg, &nzero, &npos)); 3769566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, " MatInertia: nneg: %" PetscInt_FMT ", nzero: %" PetscInt_FMT ", npos: %" PetscInt_FMT "\n", nneg, nzero, npos)); 37778bc9606SHong Zhang } 378c4762a1bSJed Brown #endif 379c4762a1bSJed Brown } 380c4762a1bSJed Brown #endif 381c4762a1bSJed Brown 382d47f36abSHong Zhang #if defined(PETSC_HAVE_MUMPS) 383d47f36abSHong Zhang /* mumps interface allows repeated call of MatCholeskyFactorSymbolic(), while the succession calls do nothing */ 384d47f36abSHong Zhang if (ipack == 2) { 385d47f36abSHong Zhang if (chol) { 3869566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorSymbolic(F, A, perm, &info)); 3879566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorNumeric(F, A, &info)); 388d47f36abSHong Zhang } else { 3899566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(F, A, perm, iperm, &info)); 3909566063dSJacob Faibussowitsch PetscCall(MatLUFactorNumeric(F, A, &info)); 391d47f36abSHong Zhang } 392d47f36abSHong Zhang } 393d47f36abSHong Zhang #endif 394d47f36abSHong Zhang 395b18964edSHong Zhang /* Test MatMatSolve(), A X = B, where B can be dense or sparse */ 396c4762a1bSJed Brown if (testMatMatSolve) { 397c4762a1bSJed Brown if (!nfact) { 3985d955bbbSStefano Zampini PetscCall(MatMatMult(Ae, C, MAT_INITIAL_MATRIX, 2.0, &RHS)); 399c4762a1bSJed Brown } else { 4005d955bbbSStefano Zampini PetscCall(MatMatMult(Ae, C, MAT_REUSE_MATRIX, 2.0, &RHS)); 401c4762a1bSJed Brown } 402c4762a1bSJed Brown for (nsolve = 0; nsolve < 2; nsolve++) { 4039566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " %" PetscInt_FMT "-the MatMatSolve \n", nsolve)); 4049566063dSJacob Faibussowitsch PetscCall(MatMatSolve(F, RHS, X)); 405c4762a1bSJed Brown 406c4762a1bSJed Brown /* Check the error */ 4079566063dSJacob Faibussowitsch PetscCall(MatAXPY(X, -1.0, C, SAME_NONZERO_PATTERN)); 4089566063dSJacob Faibussowitsch PetscCall(MatNorm(X, NORM_FROBENIUS, &norm)); 409b18964edSHong 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)); 410c4762a1bSJed Brown } 411b18964edSHong Zhang 412c4762a1bSJed Brown if (matsolvexx) { 413c4762a1bSJed Brown /* Test MatMatSolve(F,RHS,RHS), RHS is a dense matrix */ 4149566063dSJacob Faibussowitsch PetscCall(MatCopy(RHS, X, SAME_NONZERO_PATTERN)); 4159566063dSJacob Faibussowitsch PetscCall(MatMatSolve(F, X, X)); 416c4762a1bSJed Brown /* Check the error */ 4179566063dSJacob Faibussowitsch PetscCall(MatAXPY(X, -1.0, C, SAME_NONZERO_PATTERN)); 4189566063dSJacob Faibussowitsch PetscCall(MatNorm(X, NORM_FROBENIUS, &norm)); 419b18964edSHong Zhang if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatMatSolve(F,RHS,RHS): Norm of error %g\n", (double)norm)); 420c4762a1bSJed Brown } 421c4762a1bSJed Brown 422c4762a1bSJed Brown if (ipack == 2 && size == 1) { 423c4762a1bSJed Brown Mat spRHS, spRHST, RHST; 424c4762a1bSJed Brown 4259566063dSJacob Faibussowitsch PetscCall(MatTranspose(RHS, MAT_INITIAL_MATRIX, &RHST)); 4269566063dSJacob Faibussowitsch PetscCall(MatConvert(RHST, MATAIJ, MAT_INITIAL_MATRIX, &spRHST)); 4279566063dSJacob Faibussowitsch PetscCall(MatCreateTranspose(spRHST, &spRHS)); 428c4762a1bSJed Brown for (nsolve = 0; nsolve < 2; nsolve++) { 4299566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " %" PetscInt_FMT "-the sparse MatMatSolve \n", nsolve)); 4309566063dSJacob Faibussowitsch PetscCall(MatMatSolve(F, spRHS, X)); 431c4762a1bSJed Brown 432c4762a1bSJed Brown /* Check the error */ 4339566063dSJacob Faibussowitsch PetscCall(MatAXPY(X, -1.0, C, SAME_NONZERO_PATTERN)); 4349566063dSJacob Faibussowitsch PetscCall(MatNorm(X, NORM_FROBENIUS, &norm)); 435b18964edSHong 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)); 436c4762a1bSJed Brown } 437b18964edSHong Zhang PetscCall(MatDestroy(&spRHST)); 438b18964edSHong Zhang PetscCall(MatDestroy(&spRHS)); 439b18964edSHong Zhang PetscCall(MatDestroy(&RHST)); 440b18964edSHong Zhang } 441b18964edSHong Zhang } 442b18964edSHong Zhang 443b18964edSHong Zhang /* Test testMatMatSolveTranspose(), A^T X = B, where B can be dense or sparse */ 444b18964edSHong Zhang if (testMatMatSolveTranspose) { 445b18964edSHong Zhang if (!nfact) { 4465d955bbbSStefano Zampini PetscCall(MatTransposeMatMult(Ae, C, MAT_INITIAL_MATRIX, 2.0, &RHS1)); 447b18964edSHong Zhang } else { 4485d955bbbSStefano Zampini PetscCall(MatTransposeMatMult(Ae, C, MAT_REUSE_MATRIX, 2.0, &RHS1)); 449b18964edSHong Zhang } 450b18964edSHong Zhang 451b18964edSHong Zhang for (nsolve = 0; nsolve < 2; nsolve++) { 45278bc9606SHong Zhang PetscCall(PetscPrintf(PETSC_COMM_WORLD, " %" PetscInt_FMT "-the MatMatSolveTranspose\n", nsolve)); 453b18964edSHong Zhang PetscCall(MatMatSolveTranspose(F, RHS1, X)); 454b18964edSHong Zhang 455b18964edSHong Zhang /* Check the error */ 456b18964edSHong Zhang PetscCall(MatAXPY(X, -1.0, C, SAME_NONZERO_PATTERN)); 457b18964edSHong Zhang PetscCall(MatNorm(X, NORM_FROBENIUS, &norm)); 458b18964edSHong 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)); 459b18964edSHong Zhang } 460b18964edSHong Zhang 461b18964edSHong Zhang if (ipack == 2 && size == 1) { 462b18964edSHong Zhang Mat spRHS, spRHST, RHST; 463b18964edSHong Zhang 464b18964edSHong Zhang PetscCall(MatTranspose(RHS1, MAT_INITIAL_MATRIX, &RHST)); 465b18964edSHong Zhang PetscCall(MatConvert(RHST, MATAIJ, MAT_INITIAL_MATRIX, &spRHST)); 466b18964edSHong Zhang PetscCall(MatCreateTranspose(spRHST, &spRHS)); 467b18964edSHong Zhang for (nsolve = 0; nsolve < 2; nsolve++) { 468b18964edSHong Zhang PetscCall(MatMatSolveTranspose(F, spRHS, X)); 469b18964edSHong Zhang 470b18964edSHong Zhang /* Check the error */ 471b18964edSHong Zhang PetscCall(MatAXPY(X, -1.0, C, SAME_NONZERO_PATTERN)); 472b18964edSHong Zhang PetscCall(MatNorm(X, NORM_FROBENIUS, &norm)); 473b18964edSHong 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)); 474c4762a1bSJed Brown } 4759566063dSJacob Faibussowitsch PetscCall(MatDestroy(&spRHST)); 4769566063dSJacob Faibussowitsch PetscCall(MatDestroy(&spRHS)); 4779566063dSJacob Faibussowitsch PetscCall(MatDestroy(&RHST)); 478c4762a1bSJed Brown } 479c4762a1bSJed Brown } 480c4762a1bSJed Brown 481c4762a1bSJed Brown /* Test MatSolve() */ 482c4762a1bSJed Brown if (testMatSolve) { 483c4762a1bSJed Brown for (nsolve = 0; nsolve < 2; nsolve++) { 4849566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rand)); 4859566063dSJacob Faibussowitsch PetscCall(VecCopy(x, u)); 4865d955bbbSStefano Zampini PetscCall(MatMult(Ae, x, b)); 487c4762a1bSJed Brown 4889566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " %" PetscInt_FMT "-the MatSolve \n", nsolve)); 4899566063dSJacob Faibussowitsch PetscCall(MatSolve(F, b, x)); 490c4762a1bSJed Brown 491c4762a1bSJed Brown /* Check the error */ 4929566063dSJacob Faibussowitsch PetscCall(VecAXPY(u, -1.0, x)); /* u <- (-1.0)x + u */ 4939566063dSJacob Faibussowitsch PetscCall(VecNorm(u, NORM_2, &norm)); 494c4762a1bSJed Brown if (norm > tol) { 495c4762a1bSJed Brown PetscReal resi; 4965d955bbbSStefano Zampini PetscCall(MatMult(Ae, x, u)); /* u = A*x */ 4979566063dSJacob Faibussowitsch PetscCall(VecAXPY(u, -1.0, b)); /* u <- (-1.0)b + u */ 4989566063dSJacob Faibussowitsch PetscCall(VecNorm(u, NORM_2, &resi)); 4999566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatSolve: Norm of error %g, resi %g, numfact %" PetscInt_FMT "\n", (double)norm, (double)resi, nfact)); 500c4762a1bSJed Brown } 501c4762a1bSJed Brown } 502c4762a1bSJed Brown } 50378bc9606SHong Zhang 50478bc9606SHong Zhang /* Test MatSolveTranspose() */ 50578bc9606SHong Zhang if (testMatSolveTranspose) { 50678bc9606SHong Zhang for (nsolve = 0; nsolve < 2; nsolve++) { 50778bc9606SHong Zhang PetscCall(VecSetRandom(x, rand)); 50878bc9606SHong Zhang PetscCall(VecCopy(x, u)); 5095d955bbbSStefano Zampini PetscCall(MatMultTranspose(Ae, x, b)); 51078bc9606SHong Zhang 51178bc9606SHong Zhang PetscCall(PetscPrintf(PETSC_COMM_WORLD, " %" PetscInt_FMT "-the MatSolveTranspose\n", nsolve)); 51278bc9606SHong Zhang PetscCall(MatSolveTranspose(F, b, x)); 51378bc9606SHong Zhang 51478bc9606SHong Zhang /* Check the error */ 51578bc9606SHong Zhang PetscCall(VecAXPY(u, -1.0, x)); /* u <- (-1.0)x + u */ 51678bc9606SHong Zhang PetscCall(VecNorm(u, NORM_2, &norm)); 51778bc9606SHong Zhang if (norm > tol) { 51878bc9606SHong Zhang PetscReal resi; 519070736c7SStefano Zampini PetscCall(MatMultTranspose(Ae, x, u)); /* u = A*x */ 52078bc9606SHong Zhang PetscCall(VecAXPY(u, -1.0, b)); /* u <- (-1.0)b + u */ 52178bc9606SHong Zhang PetscCall(VecNorm(u, NORM_2, &resi)); 52278bc9606SHong Zhang PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatSolveTranspose: Norm of error %g, resi %g, numfact %" PetscInt_FMT "\n", (double)norm, (double)resi, nfact)); 52378bc9606SHong Zhang } 52478bc9606SHong Zhang } 52578bc9606SHong Zhang } 526c4762a1bSJed Brown } 527c4762a1bSJed Brown 528c4762a1bSJed Brown /* Free data structures */ 5295d955bbbSStefano Zampini PetscCall(MatDestroy(&Ae)); 5309566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 5319566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 5329566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F)); 5339566063dSJacob Faibussowitsch PetscCall(MatDestroy(&X)); 5349566063dSJacob Faibussowitsch PetscCall(MatDestroy(&RHS)); 535b18964edSHong Zhang PetscCall(MatDestroy(&RHS1)); 536c4762a1bSJed Brown 5379566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand)); 5389566063dSJacob Faibussowitsch PetscCall(ISDestroy(&perm)); 5399566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iperm)); 5409566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 5419566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 5429566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 5439566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 544b122ec5aSJacob Faibussowitsch return 0; 545c4762a1bSJed Brown } 546c4762a1bSJed Brown 547c4762a1bSJed Brown /*TEST 548c4762a1bSJed Brown 549c4762a1bSJed Brown test: 550dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 55178bc9606SHong Zhang args: -f ${DATAFILESPATH}/matrices/medium -mat_solver_type petsc 552c4762a1bSJed Brown output_file: output/ex125.out 553c4762a1bSJed Brown 554c4762a1bSJed Brown test: 5559a14fc28SStefano Zampini suffix: 2 55678bc9606SHong Zhang args: -mat_solver_type petsc 5579a14fc28SStefano Zampini output_file: output/ex125.out 5589a14fc28SStefano Zampini 5599a14fc28SStefano Zampini test: 560c4762a1bSJed Brown suffix: mkl_pardiso 561dfd57a17SPierre Jolivet requires: mkl_pardiso datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 56278bc9606SHong Zhang args: -f ${DATAFILESPATH}/matrices/small -mat_solver_type mkl_pardiso 563c4762a1bSJed Brown 564c4762a1bSJed Brown test: 5659a14fc28SStefano Zampini suffix: mkl_pardiso_2 5669a14fc28SStefano Zampini requires: mkl_pardiso 56778bc9606SHong Zhang args: -mat_solver_type mkl_pardiso 5689a14fc28SStefano Zampini output_file: output/ex125_mkl_pardiso.out 5699a14fc28SStefano Zampini 5709a14fc28SStefano Zampini test: 571c4762a1bSJed Brown suffix: mumps 572dfd57a17SPierre Jolivet requires: mumps datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 57378bc9606SHong Zhang args: -f ${DATAFILESPATH}/matrices/small -mat_solver_type mumps 574c4762a1bSJed Brown output_file: output/ex125_mumps_seq.out 575c4762a1bSJed Brown 576c4762a1bSJed Brown test: 57762671d91SStefano Zampini suffix: mumps_nest 57862671d91SStefano Zampini requires: mumps datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 5795d955bbbSStefano Zampini args: -f ${DATAFILESPATH}/matrices/small -mat_solver_type mumps -test_nest -test_nest_bordered {{0 1}} 58062671d91SStefano Zampini output_file: output/ex125_mumps_seq.out 58162671d91SStefano Zampini 58262671d91SStefano Zampini test: 583c4762a1bSJed Brown suffix: mumps_2 584c4762a1bSJed Brown nsize: 3 585dfd57a17SPierre Jolivet requires: mumps datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 58678bc9606SHong Zhang args: -f ${DATAFILESPATH}/matrices/small -mat_solver_type mumps 587c4762a1bSJed Brown output_file: output/ex125_mumps_par.out 588c4762a1bSJed Brown 589c4762a1bSJed Brown test: 59062671d91SStefano Zampini suffix: mumps_2_nest 59162671d91SStefano Zampini nsize: 3 59262671d91SStefano Zampini requires: mumps datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 5935d955bbbSStefano Zampini args: -f ${DATAFILESPATH}/matrices/small -mat_solver_type mumps -test_nest -test_nest_bordered {{0 1}} 59462671d91SStefano Zampini output_file: output/ex125_mumps_par.out 59562671d91SStefano Zampini 59662671d91SStefano Zampini test: 5979a14fc28SStefano Zampini suffix: mumps_3 5989a14fc28SStefano Zampini requires: mumps 59978bc9606SHong Zhang args: -mat_solver_type mumps 6009a14fc28SStefano Zampini output_file: output/ex125_mumps_seq.out 6019a14fc28SStefano Zampini 6029a14fc28SStefano Zampini test: 60362671d91SStefano Zampini suffix: mumps_3_nest 60462671d91SStefano Zampini requires: mumps 6055d955bbbSStefano Zampini args: -mat_solver_type mumps -test_nest -test_nest_bordered {{0 1}} 60662671d91SStefano Zampini output_file: output/ex125_mumps_seq.out 60762671d91SStefano Zampini 60862671d91SStefano Zampini test: 6099a14fc28SStefano Zampini suffix: mumps_4 6109a14fc28SStefano Zampini nsize: 3 6119a14fc28SStefano Zampini requires: mumps 61278bc9606SHong Zhang args: -mat_solver_type mumps 6139a14fc28SStefano Zampini output_file: output/ex125_mumps_par.out 6149a14fc28SStefano Zampini 6159a14fc28SStefano Zampini test: 61662671d91SStefano Zampini suffix: mumps_4_nest 61762671d91SStefano Zampini nsize: 3 61862671d91SStefano Zampini requires: mumps 6195d955bbbSStefano Zampini args: -mat_solver_type mumps -test_nest -test_nest_bordered {{0 1}} 62062671d91SStefano Zampini output_file: output/ex125_mumps_par.out 62162671d91SStefano Zampini 62262671d91SStefano Zampini test: 623d47f36abSHong Zhang suffix: mumps_5 624d47f36abSHong Zhang nsize: 3 625d47f36abSHong Zhang requires: mumps 62678bc9606SHong Zhang args: -mat_solver_type mumps -cholesky 627d47f36abSHong Zhang output_file: output/ex125_mumps_par_cholesky.out 628d47f36abSHong Zhang 629d47f36abSHong Zhang test: 63062671d91SStefano Zampini suffix: mumps_5_nest 63162671d91SStefano Zampini nsize: 3 63262671d91SStefano Zampini requires: mumps 6335d955bbbSStefano Zampini args: -mat_solver_type mumps -cholesky -test_nest -test_nest_bordered {{0 1}} 63462671d91SStefano Zampini output_file: output/ex125_mumps_par_cholesky.out 63562671d91SStefano Zampini 63662671d91SStefano Zampini test: 63753587d93SPierre Jolivet suffix: mumps_6 63853587d93SPierre Jolivet nsize: 2 63953587d93SPierre Jolivet requires: mumps 64053587d93SPierre Jolivet args: -mat_solver_type mumps -test_nest -test_nest_bordered -m 13 -n 13 64153587d93SPierre Jolivet output_file: output/ex125_mumps_par.out 64253587d93SPierre Jolivet 64353587d93SPierre Jolivet test: 64478bc9606SHong Zhang suffix: superlu 64578bc9606SHong Zhang requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) superlu 64678bc9606SHong Zhang args: -f ${DATAFILESPATH}/matrices/medium -mat_solver_type superlu 64778bc9606SHong Zhang output_file: output/ex125_superlu.out 64878bc9606SHong Zhang 64978bc9606SHong Zhang test: 650c4762a1bSJed Brown suffix: superlu_dist 6519a14fc28SStefano Zampini nsize: {{1 3}} 652d4783600SBarry Smith requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) superlu_dist 65378bc9606SHong Zhang args: -f ${DATAFILESPATH}/matrices/small -mat_solver_type superlu_dist -mat_superlu_dist_rowperm NOROWPERM 65478bc9606SHong Zhang output_file: output/ex125_superlu_dist.out 655c4762a1bSJed Brown 656c4762a1bSJed Brown test: 657c4762a1bSJed Brown suffix: superlu_dist_2 6589a14fc28SStefano Zampini nsize: {{1 3}} 6599a14fc28SStefano Zampini requires: superlu_dist !complex 66078bc9606SHong Zhang args: -n 36 -mat_solver_type superlu_dist -mat_superlu_dist_rowperm NOROWPERM 661c4762a1bSJed Brown output_file: output/ex125_superlu_dist.out 662c4762a1bSJed Brown 663c4762a1bSJed Brown test: 66478bc9606SHong Zhang suffix: superlu_dist_3 66578bc9606SHong Zhang nsize: {{1 3}} 66678bc9606SHong Zhang requires: superlu_dist !complex 66778bc9606SHong Zhang requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) superlu_dist 66878bc9606SHong Zhang args: -f ${DATAFILESPATH}/matrices/medium -mat_solver_type superlu_dist -mat_superlu_dist_rowperm NOROWPERM 66978bc9606SHong Zhang output_file: output/ex125_superlu_dist_nonsymmetric.out 67078bc9606SHong Zhang 67178bc9606SHong Zhang test: 672c4762a1bSJed Brown suffix: superlu_dist_complex 673c4762a1bSJed Brown nsize: 3 674d4783600SBarry Smith requires: datafilespath double superlu_dist complex !defined(PETSC_USE_64BIT_INDICES) 67578bc9606SHong Zhang args: -f ${DATAFILESPATH}/matrices/farzad_B_rhs -mat_solver_type superlu_dist 676c4762a1bSJed Brown output_file: output/ex125_superlu_dist_complex.out 677c4762a1bSJed Brown 67838a8e8c1SStefano Zampini test: 6799a14fc28SStefano Zampini suffix: superlu_dist_complex_2 6809a14fc28SStefano Zampini nsize: 3 6819a14fc28SStefano Zampini requires: superlu_dist complex 68278bc9606SHong Zhang args: -mat_solver_type superlu_dist 68378bc9606SHong Zhang output_file: output/ex125_superlu_dist_complex_2.out 6849a14fc28SStefano Zampini 6859a14fc28SStefano Zampini test: 68638a8e8c1SStefano Zampini suffix: cusparse 687dfd57a17SPierre Jolivet requires: cuda datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 6880338c944SBarry Smith #TODO: fix the bug with cholesky 68978bc9606SHong Zhang #args: -mat_type aijcusparse -f ${DATAFILESPATH}/matrices/small -mat_solver_type cusparse -cholesky {{0 1}separate output} 69078bc9606SHong Zhang args: -mat_type aijcusparse -f ${DATAFILESPATH}/matrices/small -mat_solver_type cusparse -cholesky {{0}separate output} 69138a8e8c1SStefano Zampini 6929a14fc28SStefano Zampini test: 6939a14fc28SStefano Zampini suffix: cusparse_2 6949a14fc28SStefano Zampini requires: cuda 69578bc9606SHong Zhang args: -mat_type aijcusparse -mat_solver_type cusparse -cholesky {{0 1}separate output} 6969a14fc28SStefano Zampini 697*39989c4fSStefano Zampini testset: 698*39989c4fSStefano Zampini nsize: {{1 2}separate output} 699*39989c4fSStefano Zampini requires: double !defined(PETSC_USE_64BIT_INDICES) datafilespath !complex 700*39989c4fSStefano Zampini args: -f ${DATAFILESPATH}/matrices/mixed_poisson 701*39989c4fSStefano Zampini test: 702*39989c4fSStefano Zampini requires: superlu_dist TODO # superlu_dist is broken 703*39989c4fSStefano Zampini suffix: saddle_point_superlu_dist 704*39989c4fSStefano Zampini args: -mat_solver_type superlu_dist -mat_superlu_dist_rowperm {{norowperm largediag_mc64}} -test_inertia 0 705*39989c4fSStefano Zampini test: 706*39989c4fSStefano Zampini requires: mumps 707*39989c4fSStefano Zampini suffix: saddle_point_mumps_lu 708*39989c4fSStefano Zampini args: -mat_solver_type mumps -mat_mumps_icntl_14 100 709*39989c4fSStefano Zampini test: 710*39989c4fSStefano Zampini requires: mumps 711*39989c4fSStefano Zampini suffix: saddle_point_mumps_cholesky 712*39989c4fSStefano Zampini args: -cholesky -mat_solver_type mumps 713*39989c4fSStefano Zampini 714c4762a1bSJed Brown TEST*/ 715