1c4762a1bSJed Brown static char help[] = "Tests MatSolve() and MatMatSolve() with MUMPS or MKL_PARDISO sequential solvers in Schur complement mode.\n\ 2c4762a1bSJed Brown Example: mpiexec -n 1 ./ex192 -f <matrix binary file> -nrhs 4 -symmetric_solve -hermitian_solve -schur_ratio 0.3\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 7d71ae5a4SJacob Faibussowitsch { 8c4762a1bSJed Brown Mat A, RHS, C, F, X, S; 9c4762a1bSJed Brown Vec u, x, b; 10c4762a1bSJed Brown Vec xschur, bschur, uschur; 11c4762a1bSJed Brown IS is_schur; 12c4762a1bSJed Brown PetscMPIInt size; 13c4762a1bSJed Brown PetscInt isolver = 0, size_schur, m, n, nfact, nsolve, nrhs; 14c4762a1bSJed Brown PetscReal norm, tol = PETSC_SQRT_MACHINE_EPSILON; 15c4762a1bSJed Brown PetscRandom rand; 16c4762a1bSJed Brown PetscBool data_provided, herm, symm, use_lu, cuda = PETSC_FALSE; 1779606b51SStefano Zampini PetscBool isdata_provided; 18c4762a1bSJed Brown PetscReal sratio = 5.1 / 12.; 19c4762a1bSJed Brown PetscViewer fd; /* viewer */ 20c4762a1bSJed Brown char solver[256]; 2179606b51SStefano Zampini char file[PETSC_MAX_PATH_LEN]; /* input Mat file name */ 2279606b51SStefano Zampini char isfile[PETSC_MAX_PATH_LEN]; /* input IS file name */ 23c4762a1bSJed Brown 24327415f7SBarry Smith PetscFunctionBeginUser; 25c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help)); 269566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 27be096a46SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor test"); 28c4762a1bSJed Brown /* Determine which type of solver we want to test for */ 29c4762a1bSJed Brown herm = PETSC_FALSE; 30c4762a1bSJed Brown symm = PETSC_FALSE; 319566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-symmetric_solve", &symm, NULL)); 329566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-hermitian_solve", &herm, NULL)); 33c4762a1bSJed Brown if (herm) symm = PETSC_TRUE; 349566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-cuda_solve", &cuda, NULL)); 359566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-tol", &tol, NULL)); 36c4762a1bSJed Brown 37c4762a1bSJed Brown /* Determine file from which we read the matrix A */ 389566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-f", file, sizeof(file), &data_provided)); 39c4762a1bSJed Brown if (!data_provided) { /* get matrices from PETSc distribution */ 409566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(file, "${PETSC_DIR}/share/petsc/datafiles/matrices/", sizeof(file))); 41c4762a1bSJed Brown if (symm) { 42c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX) 439566063dSJacob Faibussowitsch PetscCall(PetscStrlcat(file, "hpd-complex-", sizeof(file))); 44c4762a1bSJed Brown #else 459566063dSJacob Faibussowitsch PetscCall(PetscStrlcat(file, "spd-real-", sizeof(file))); 46c4762a1bSJed Brown #endif 47c4762a1bSJed Brown } else { 48c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX) 499566063dSJacob Faibussowitsch PetscCall(PetscStrlcat(file, "nh-complex-", sizeof(file))); 50c4762a1bSJed Brown #else 519566063dSJacob Faibussowitsch PetscCall(PetscStrlcat(file, "ns-real-", sizeof(file))); 52c4762a1bSJed Brown #endif 53c4762a1bSJed Brown } 54c4762a1bSJed Brown #if defined(PETSC_USE_64BIT_INDICES) 559566063dSJacob Faibussowitsch PetscCall(PetscStrlcat(file, "int64-", sizeof(file))); 56c4762a1bSJed Brown #else 579566063dSJacob Faibussowitsch PetscCall(PetscStrlcat(file, "int32-", sizeof(file))); 58c4762a1bSJed Brown #endif 59c4762a1bSJed Brown #if defined(PETSC_USE_REAL_SINGLE) 609566063dSJacob Faibussowitsch PetscCall(PetscStrlcat(file, "float32", sizeof(file))); 61c4762a1bSJed Brown #else 629566063dSJacob Faibussowitsch PetscCall(PetscStrlcat(file, "float64", sizeof(file))); 63c4762a1bSJed Brown #endif 64c4762a1bSJed Brown } 6579606b51SStefano Zampini 66c4762a1bSJed Brown /* Load matrix A */ 679566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &fd)); 689566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 699566063dSJacob Faibussowitsch PetscCall(MatLoad(A, fd)); 709566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &m, &n)); 7108401ef6SPierre Jolivet PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "This example is not intended for rectangular matrices (%" PetscInt_FMT ", %" PetscInt_FMT ")", m, n); 72c4762a1bSJed Brown 7379606b51SStefano Zampini PetscCall(PetscOptionsGetString(NULL, NULL, "-fis", isfile, sizeof(isfile), &isdata_provided)); 7479606b51SStefano Zampini if (isdata_provided) { 7579606b51SStefano Zampini PetscBool samefile; 7679606b51SStefano Zampini 7779606b51SStefano Zampini PetscCall(PetscStrcmp(isfile, file, &samefile)); 7879606b51SStefano Zampini if (!samefile) { 7979606b51SStefano Zampini PetscCall(PetscViewerDestroy(&fd)); 8079606b51SStefano Zampini PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, isfile, FILE_MODE_READ, &fd)); 8179606b51SStefano Zampini } 8279606b51SStefano Zampini PetscCall(ISCreate(PETSC_COMM_SELF, &is_schur)); 8379606b51SStefano Zampini PetscCall(ISLoad(is_schur, fd)); 8479606b51SStefano Zampini } else { 8579606b51SStefano Zampini PetscCall(PetscOptionsGetReal(NULL, NULL, "-schur_ratio", &sratio, NULL)); 8679606b51SStefano Zampini PetscCheck(sratio >= 0. && sratio <= 1., PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Invalid ratio for schur degrees of freedom %g", (double)sratio); 8779606b51SStefano Zampini size_schur = (PetscInt)(sratio * m); 8879606b51SStefano Zampini PetscCall(ISCreateStride(PETSC_COMM_SELF, size_schur, m - size_schur, 1, &is_schur)); 8979606b51SStefano Zampini } 9079606b51SStefano Zampini PetscCall(ISGetSize(is_schur, &size_schur)); 9179606b51SStefano Zampini PetscCall(PetscViewerDestroy(&fd)); 9279606b51SStefano Zampini 93a5b23f4aSJose E. Roman /* Create dense matrix C and X; C holds true solution with identical columns */ 94c4762a1bSJed Brown nrhs = 2; 959566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-nrhs", &nrhs, NULL)); 969566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &C)); 979566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, m, PETSC_DECIDE, PETSC_DECIDE, nrhs)); 989566063dSJacob Faibussowitsch PetscCall(MatSetType(C, MATDENSE)); 999566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C)); 1009566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 101c4762a1bSJed Brown 1029566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rand)); 1039566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rand)); 1049566063dSJacob Faibussowitsch PetscCall(MatSetRandom(C, rand)); 1059566063dSJacob Faibussowitsch PetscCall(MatDuplicate(C, MAT_DO_NOT_COPY_VALUES, &X)); 106c4762a1bSJed Brown 107c4762a1bSJed Brown /* Create vectors */ 1089566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &x)); 1099566063dSJacob Faibussowitsch PetscCall(VecSetSizes(x, n, PETSC_DECIDE)); 1109566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(x)); 1119566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &b)); 1129566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &u)); /* save the true solution */ 113c4762a1bSJed Brown 1149566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-solver", &isolver, NULL)); 115c4762a1bSJed Brown switch (isolver) { 116c4762a1bSJed Brown #if defined(PETSC_HAVE_MUMPS) 117d71ae5a4SJacob Faibussowitsch case 0: 118c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(solver, MATSOLVERMUMPS, sizeof(solver))); 119d71ae5a4SJacob Faibussowitsch break; 120c4762a1bSJed Brown #endif 121c4762a1bSJed Brown #if defined(PETSC_HAVE_MKL_PARDISO) 122d71ae5a4SJacob Faibussowitsch case 1: 123c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(solver, MATSOLVERMKL_PARDISO, sizeof(solver))); 124d71ae5a4SJacob Faibussowitsch break; 125c4762a1bSJed Brown #endif 126d71ae5a4SJacob Faibussowitsch default: 127c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(solver, MATSOLVERPETSC, sizeof(solver))); 128d71ae5a4SJacob Faibussowitsch break; 129c4762a1bSJed Brown } 130c4762a1bSJed Brown 131c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX) 132*b0c98d1dSPierre Jolivet if (isolver == 0 && symm && !data_provided) { /* MUMPS (5.0.0) does not have support for Hermitian matrices, so make them symmetric */ 133c4762a1bSJed Brown PetscScalar im = PetscSqrtScalar((PetscScalar)-1.); 134c4762a1bSJed Brown PetscScalar val = -1.0; 135c4762a1bSJed Brown val = val + im; 1369566063dSJacob Faibussowitsch PetscCall(MatSetValue(A, 1, 0, val, INSERT_VALUES)); 1379566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1389566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 139c4762a1bSJed Brown } 140c4762a1bSJed Brown #endif 141c4762a1bSJed Brown 1429566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Solving with %s: nrhs %" PetscInt_FMT ", sym %d, herm %d, size schur %" PetscInt_FMT ", size mat %" PetscInt_FMT "\n", solver, nrhs, symm, herm, size_schur, m)); 143c4762a1bSJed Brown 144c4762a1bSJed Brown /* Test LU/Cholesky Factorization */ 145c4762a1bSJed Brown use_lu = PETSC_FALSE; 146c4762a1bSJed Brown if (!symm) use_lu = PETSC_TRUE; 147c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX) 148c4762a1bSJed Brown if (isolver == 1) use_lu = PETSC_TRUE; 149c4762a1bSJed Brown #endif 150c4762a1bSJed Brown if (cuda && symm && !herm) use_lu = PETSC_TRUE; 151c4762a1bSJed Brown 152c4762a1bSJed Brown if (herm && !use_lu) { /* test also conversion routines inside the solver packages */ 1539566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_SYMMETRIC, PETSC_TRUE)); 1549566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATSEQSBAIJ, MAT_INPLACE_MATRIX, &A)); 155c4762a1bSJed Brown } 156c4762a1bSJed Brown 157c4762a1bSJed Brown if (use_lu) { 1589566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A, solver, MAT_FACTOR_LU, &F)); 159c4762a1bSJed Brown } else { 160c4762a1bSJed Brown if (herm) { 1619566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_SPD, PETSC_TRUE)); 162c4762a1bSJed Brown } else { 1639566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_SYMMETRIC, PETSC_TRUE)); 1649566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_SPD, PETSC_FALSE)); 165c4762a1bSJed Brown } 1669566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A, solver, MAT_FACTOR_CHOLESKY, &F)); 167c4762a1bSJed Brown } 168c4762a1bSJed Brown 16979606b51SStefano Zampini /* Set Schur complement indices */ 17079606b51SStefano Zampini PetscCall(MatFactorSetSchurIS(F, is_schur)); 1719566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_schur)); 17279606b51SStefano Zampini 173c4762a1bSJed Brown if (use_lu) { 1749566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(F, A, NULL, NULL, NULL)); 175c4762a1bSJed Brown } else { 1769566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorSymbolic(F, A, NULL, NULL)); 177c4762a1bSJed Brown } 178c4762a1bSJed Brown 179c4762a1bSJed Brown for (nfact = 0; nfact < 3; nfact++) { 180c4762a1bSJed Brown Mat AD; 181c4762a1bSJed Brown 18279606b51SStefano Zampini if (nfact == 1) { 1839566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rand)); 18448a46eb9SPierre Jolivet if (symm && herm) PetscCall(VecAbs(x)); 1859566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(A, x, ADD_VALUES)); 186c4762a1bSJed Brown } 187c4762a1bSJed Brown if (use_lu) { 1889566063dSJacob Faibussowitsch PetscCall(MatLUFactorNumeric(F, A, NULL)); 189c4762a1bSJed Brown } else { 1909566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorNumeric(F, A, NULL)); 191c4762a1bSJed Brown } 19279606b51SStefano Zampini 193c4762a1bSJed Brown if (cuda) { 1949566063dSJacob Faibussowitsch PetscCall(MatFactorGetSchurComplement(F, &S, NULL)); 1959566063dSJacob Faibussowitsch PetscCall(MatSetType(S, MATSEQDENSECUDA)); 1969566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(S, &xschur, &bschur)); 1979566063dSJacob Faibussowitsch PetscCall(MatFactorRestoreSchurComplement(F, &S, MAT_FACTOR_SCHUR_UNFACTORED)); 198c4762a1bSJed Brown } 1999566063dSJacob Faibussowitsch PetscCall(MatFactorCreateSchurComplement(F, &S, NULL)); 20048a46eb9SPierre Jolivet if (!cuda) PetscCall(MatCreateVecs(S, &xschur, &bschur)); 2019566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xschur, &uschur)); 20248a46eb9SPierre Jolivet if (nfact == 1 && (!cuda || (herm && symm))) PetscCall(MatFactorInvertSchurComplement(F)); 203c4762a1bSJed Brown for (nsolve = 0; nsolve < 2; nsolve++) { 2049566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rand)); 2059566063dSJacob Faibussowitsch PetscCall(VecCopy(x, u)); 206c4762a1bSJed Brown 207c4762a1bSJed Brown if (nsolve) { 2089566063dSJacob Faibussowitsch PetscCall(MatMult(A, x, b)); 2099566063dSJacob Faibussowitsch PetscCall(MatSolve(F, b, x)); 210c4762a1bSJed Brown } else { 2119566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, x, b)); 2129566063dSJacob Faibussowitsch PetscCall(MatSolveTranspose(F, b, x)); 213c4762a1bSJed Brown } 214c4762a1bSJed Brown /* Check the error */ 2159566063dSJacob Faibussowitsch PetscCall(VecAXPY(u, -1.0, x)); /* u <- (-1.0)x + u */ 2169566063dSJacob Faibussowitsch PetscCall(VecNorm(u, NORM_2, &norm)); 217c4762a1bSJed Brown if (norm > tol) { 218c4762a1bSJed Brown PetscReal resi; 219c4762a1bSJed Brown if (nsolve) { 2209566063dSJacob Faibussowitsch PetscCall(MatMult(A, x, u)); /* u = A*x */ 221c4762a1bSJed Brown } else { 2229566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, x, u)); /* u = A*x */ 223c4762a1bSJed Brown } 2249566063dSJacob Faibussowitsch PetscCall(VecAXPY(u, -1.0, b)); /* u <- (-1.0)b + u */ 2259566063dSJacob Faibussowitsch PetscCall(VecNorm(u, NORM_2, &resi)); 226c4762a1bSJed Brown if (nsolve) { 2279566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "(f %" PetscInt_FMT ", s %" PetscInt_FMT ") MatSolve error: Norm of error %g, residual %g\n", nfact, nsolve, (double)norm, (double)resi)); 228c4762a1bSJed Brown } else { 2299566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "(f %" PetscInt_FMT ", s %" PetscInt_FMT ") MatSolveTranspose error: Norm of error %g, residual %f\n", nfact, nsolve, (double)norm, (double)resi)); 230c4762a1bSJed Brown } 231c4762a1bSJed Brown } 2329566063dSJacob Faibussowitsch PetscCall(VecSetRandom(xschur, rand)); 2339566063dSJacob Faibussowitsch PetscCall(VecCopy(xschur, uschur)); 234c4762a1bSJed Brown if (nsolve) { 2359566063dSJacob Faibussowitsch PetscCall(MatMult(S, xschur, bschur)); 2369566063dSJacob Faibussowitsch PetscCall(MatFactorSolveSchurComplement(F, bschur, xschur)); 237c4762a1bSJed Brown } else { 2389566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(S, xschur, bschur)); 2399566063dSJacob Faibussowitsch PetscCall(MatFactorSolveSchurComplementTranspose(F, bschur, xschur)); 240c4762a1bSJed Brown } 241c4762a1bSJed Brown /* Check the error */ 2429566063dSJacob Faibussowitsch PetscCall(VecAXPY(uschur, -1.0, xschur)); /* u <- (-1.0)x + u */ 2439566063dSJacob Faibussowitsch PetscCall(VecNorm(uschur, NORM_2, &norm)); 244c4762a1bSJed Brown if (norm > tol) { 245c4762a1bSJed Brown PetscReal resi; 246c4762a1bSJed Brown if (nsolve) { 2479566063dSJacob Faibussowitsch PetscCall(MatMult(S, xschur, uschur)); /* u = A*x */ 248c4762a1bSJed Brown } else { 2499566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(S, xschur, uschur)); /* u = A*x */ 250c4762a1bSJed Brown } 2519566063dSJacob Faibussowitsch PetscCall(VecAXPY(uschur, -1.0, bschur)); /* u <- (-1.0)b + u */ 2529566063dSJacob Faibussowitsch PetscCall(VecNorm(uschur, NORM_2, &resi)); 253c4762a1bSJed Brown if (nsolve) { 2549566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "(f %" PetscInt_FMT ", s %" PetscInt_FMT ") MatFactorSolveSchurComplement error: Norm of error %g, residual %g\n", nfact, nsolve, (double)norm, (double)resi)); 255c4762a1bSJed Brown } else { 2569566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "(f %" PetscInt_FMT ", s %" PetscInt_FMT ") MatFactorSolveSchurComplementTranspose error: Norm of error %g, residual %f\n", nfact, nsolve, (double)norm, (double)resi)); 257c4762a1bSJed Brown } 258c4762a1bSJed Brown } 259c4762a1bSJed Brown } 2609566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATSEQAIJ, MAT_INITIAL_MATRIX, &AD)); 261c4762a1bSJed Brown if (!nfact) { 2629566063dSJacob Faibussowitsch PetscCall(MatMatMult(AD, C, MAT_INITIAL_MATRIX, 2.0, &RHS)); 263c4762a1bSJed Brown } else { 2649566063dSJacob Faibussowitsch PetscCall(MatMatMult(AD, C, MAT_REUSE_MATRIX, 2.0, &RHS)); 265c4762a1bSJed Brown } 2669566063dSJacob Faibussowitsch PetscCall(MatDestroy(&AD)); 267c4762a1bSJed Brown for (nsolve = 0; nsolve < 2; nsolve++) { 2689566063dSJacob Faibussowitsch PetscCall(MatMatSolve(F, RHS, X)); 269c4762a1bSJed Brown 270c4762a1bSJed Brown /* Check the error */ 2719566063dSJacob Faibussowitsch PetscCall(MatAXPY(X, -1.0, C, SAME_NONZERO_PATTERN)); 2729566063dSJacob Faibussowitsch PetscCall(MatNorm(X, NORM_FROBENIUS, &norm)); 27348a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "(f %" PetscInt_FMT ", s %" PetscInt_FMT ") MatMatSolve: Norm of error %g\n", nfact, nsolve, (double)norm)); 2743e5b40d0SPierre Jolivet #if PetscDefined(HAVE_MUMPS) 2753e5b40d0SPierre Jolivet PetscCall(MatMumpsSetIcntl(F, 26, 1)); 2763e5b40d0SPierre Jolivet PetscCall(MatMatSolve(F, RHS, X)); 2773e5b40d0SPierre Jolivet PetscCall(MatMumpsSetIcntl(F, 26, 2)); 2783e5b40d0SPierre Jolivet PetscCall(MatMatSolve(F, RHS, X)); 2793e5b40d0SPierre Jolivet PetscCall(MatMumpsSetIcntl(F, 26, -1)); 2803e5b40d0SPierre Jolivet 2813e5b40d0SPierre Jolivet /* Check the error */ 2823e5b40d0SPierre Jolivet PetscCall(MatAXPY(X, -1.0, C, SAME_NONZERO_PATTERN)); 2833e5b40d0SPierre Jolivet PetscCall(MatNorm(X, NORM_FROBENIUS, &norm)); 2843e5b40d0SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "(f %" PetscInt_FMT ", s %" PetscInt_FMT ") MatMatSolve: Norm of error %g\n", nfact, nsolve, (double)norm)); 2853e5b40d0SPierre Jolivet #endif 286c4762a1bSJed Brown } 287c4762a1bSJed Brown if (isolver == 0) { 288c4762a1bSJed Brown Mat spRHS, spRHST, RHST; 289c4762a1bSJed Brown 2909566063dSJacob Faibussowitsch PetscCall(MatTranspose(RHS, MAT_INITIAL_MATRIX, &RHST)); 2919566063dSJacob Faibussowitsch PetscCall(MatConvert(RHST, MATSEQAIJ, MAT_INITIAL_MATRIX, &spRHST)); 2929566063dSJacob Faibussowitsch PetscCall(MatCreateTranspose(spRHST, &spRHS)); 293c4762a1bSJed Brown for (nsolve = 0; nsolve < 2; nsolve++) { 2949566063dSJacob Faibussowitsch PetscCall(MatMatSolve(F, spRHS, X)); 295c4762a1bSJed Brown 296c4762a1bSJed Brown /* Check the error */ 2979566063dSJacob Faibussowitsch PetscCall(MatAXPY(X, -1.0, C, SAME_NONZERO_PATTERN)); 2989566063dSJacob Faibussowitsch PetscCall(MatNorm(X, NORM_FROBENIUS, &norm)); 29948a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "(f %" PetscInt_FMT ", s %" PetscInt_FMT ") sparse MatMatSolve: Norm of error %g\n", nfact, nsolve, (double)norm)); 300c4762a1bSJed Brown } 3019566063dSJacob Faibussowitsch PetscCall(MatDestroy(&spRHST)); 3029566063dSJacob Faibussowitsch PetscCall(MatDestroy(&spRHS)); 3039566063dSJacob Faibussowitsch PetscCall(MatDestroy(&RHST)); 304c4762a1bSJed Brown } 3059566063dSJacob Faibussowitsch PetscCall(MatDestroy(&S)); 3069566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xschur)); 3079566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bschur)); 3089566063dSJacob Faibussowitsch PetscCall(VecDestroy(&uschur)); 309c4762a1bSJed Brown } 310c4762a1bSJed Brown /* Free data structures */ 3119566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 3129566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 3139566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F)); 3149566063dSJacob Faibussowitsch PetscCall(MatDestroy(&X)); 3159566063dSJacob Faibussowitsch PetscCall(MatDestroy(&RHS)); 3169566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand)); 3179566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 3189566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 3199566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 3209566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 321b122ec5aSJacob Faibussowitsch return 0; 322c4762a1bSJed Brown } 323c4762a1bSJed Brown 324c4762a1bSJed Brown /*TEST 325c4762a1bSJed Brown 326c4762a1bSJed Brown testset: 32738cf7977SStefano Zampini requires: mkl_pardiso double !complex 328c4762a1bSJed Brown args: -solver 1 329c4762a1bSJed Brown 330c4762a1bSJed Brown test: 331c4762a1bSJed Brown suffix: mkl_pardiso 332c4762a1bSJed Brown test: 333c4762a1bSJed Brown requires: cuda 334c4762a1bSJed Brown suffix: mkl_pardiso_cuda 335c4762a1bSJed Brown args: -cuda_solve 336c4762a1bSJed Brown output_file: output/ex192_mkl_pardiso.out 337c4762a1bSJed Brown test: 338c4762a1bSJed Brown suffix: mkl_pardiso_1 339c4762a1bSJed Brown args: -symmetric_solve 340c4762a1bSJed Brown output_file: output/ex192_mkl_pardiso_1.out 341c4762a1bSJed Brown test: 342c4762a1bSJed Brown requires: cuda 343c4762a1bSJed Brown suffix: mkl_pardiso_cuda_1 344c4762a1bSJed Brown args: -symmetric_solve -cuda_solve 345c4762a1bSJed Brown output_file: output/ex192_mkl_pardiso_1.out 346c4762a1bSJed Brown test: 347c4762a1bSJed Brown suffix: mkl_pardiso_3 348c4762a1bSJed Brown args: -symmetric_solve -hermitian_solve 349c4762a1bSJed Brown output_file: output/ex192_mkl_pardiso_3.out 350c4762a1bSJed Brown test: 351dfd57a17SPierre Jolivet requires: cuda defined(PETSC_HAVE_CUSOLVERDNDPOTRI) 352c4762a1bSJed Brown suffix: mkl_pardiso_cuda_3 353c4762a1bSJed Brown args: -symmetric_solve -hermitian_solve -cuda_solve 354c4762a1bSJed Brown output_file: output/ex192_mkl_pardiso_3.out 355c4762a1bSJed Brown 356c4762a1bSJed Brown testset: 357c4762a1bSJed Brown requires: mumps double !complex 358c4762a1bSJed Brown args: -solver 0 359c4762a1bSJed Brown 360c4762a1bSJed Brown test: 361c4762a1bSJed Brown suffix: mumps 362c4762a1bSJed Brown test: 363c4762a1bSJed Brown requires: cuda 364c4762a1bSJed Brown suffix: mumps_cuda 365c4762a1bSJed Brown args: -cuda_solve 366c4762a1bSJed Brown output_file: output/ex192_mumps.out 367c4762a1bSJed Brown test: 368c4762a1bSJed Brown suffix: mumps_2 369c4762a1bSJed Brown args: -symmetric_solve 370c4762a1bSJed Brown output_file: output/ex192_mumps_2.out 371c4762a1bSJed Brown test: 372c4762a1bSJed Brown requires: cuda 373c4762a1bSJed Brown suffix: mumps_cuda_2 374c4762a1bSJed Brown args: -symmetric_solve -cuda_solve 375c4762a1bSJed Brown output_file: output/ex192_mumps_2.out 376c4762a1bSJed Brown test: 377c4762a1bSJed Brown suffix: mumps_3 378c4762a1bSJed Brown args: -symmetric_solve -hermitian_solve 379c4762a1bSJed Brown output_file: output/ex192_mumps_3.out 380c4762a1bSJed Brown test: 381dfd57a17SPierre Jolivet requires: cuda defined(PETSC_HAVE_CUSOLVERDNDPOTRI) 382c4762a1bSJed Brown suffix: mumps_cuda_3 383c4762a1bSJed Brown args: -symmetric_solve -hermitian_solve -cuda_solve 384c4762a1bSJed Brown output_file: output/ex192_mumps_3.out 385c4762a1bSJed Brown 386cf053153SJunchao Zhang testset: 387cf053153SJunchao Zhang requires: mumps double !complex defined(PETSC_HAVE_MUMPS_MIXED_PRECISION) 388cf053153SJunchao Zhang args: -solver 0 -pc_precision single -tol 3.4e-4 389cf053153SJunchao Zhang 390cf053153SJunchao Zhang test: 391cf053153SJunchao Zhang suffix: mumps_s 392cf053153SJunchao Zhang output_file: output/ex192_mumps.out 393cf053153SJunchao Zhang 394cf053153SJunchao Zhang test: 395cf053153SJunchao Zhang requires: cuda 396cf053153SJunchao Zhang suffix: mumps_cuda_s 397cf053153SJunchao Zhang args: -cuda_solve 398cf053153SJunchao Zhang output_file: output/ex192_mumps.out 399cf053153SJunchao Zhang test: 400cf053153SJunchao Zhang suffix: mumps_2_s 401cf053153SJunchao Zhang args: -symmetric_solve 402cf053153SJunchao Zhang output_file: output/ex192_mumps_2.out 403cf053153SJunchao Zhang test: 404cf053153SJunchao Zhang requires: cuda 405cf053153SJunchao Zhang suffix: mumps_cuda_2_s 406cf053153SJunchao Zhang args: -symmetric_solve -cuda_solve 407cf053153SJunchao Zhang output_file: output/ex192_mumps_2.out 408cf053153SJunchao Zhang test: 409cf053153SJunchao Zhang suffix: mumps_3_s 410cf053153SJunchao Zhang args: -symmetric_solve -hermitian_solve 411cf053153SJunchao Zhang output_file: output/ex192_mumps_3.out 412cf053153SJunchao Zhang test: 413cf053153SJunchao Zhang requires: cuda defined(PETSC_HAVE_CUSOLVERDNDPOTRI) 414cf053153SJunchao Zhang suffix: mumps_cuda_3_s 415cf053153SJunchao Zhang args: -symmetric_solve -hermitian_solve -cuda_solve 416cf053153SJunchao Zhang output_file: output/ex192_mumps_3.out 417cf053153SJunchao Zhang 418c4762a1bSJed Brown TEST*/ 419