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; 17c4762a1bSJed Brown PetscReal sratio = 5.1 / 12.; 18c4762a1bSJed Brown PetscViewer fd; /* viewer */ 19c4762a1bSJed Brown char solver[256]; 20c4762a1bSJed Brown char file[PETSC_MAX_PATH_LEN]; /* input file name */ 21c4762a1bSJed Brown 22327415f7SBarry Smith PetscFunctionBeginUser; 239566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 249566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 25be096a46SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor test"); 26c4762a1bSJed Brown /* Determine which type of solver we want to test for */ 27c4762a1bSJed Brown herm = PETSC_FALSE; 28c4762a1bSJed Brown symm = PETSC_FALSE; 299566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-symmetric_solve", &symm, NULL)); 309566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-hermitian_solve", &herm, NULL)); 31c4762a1bSJed Brown if (herm) symm = PETSC_TRUE; 329566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-cuda_solve", &cuda, NULL)); 339566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-tol", &tol, NULL)); 34c4762a1bSJed Brown 35c4762a1bSJed Brown /* Determine file from which we read the matrix A */ 369566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-f", file, sizeof(file), &data_provided)); 37c4762a1bSJed Brown if (!data_provided) { /* get matrices from PETSc distribution */ 389566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(file, "${PETSC_DIR}/share/petsc/datafiles/matrices/", sizeof(file))); 39c4762a1bSJed Brown if (symm) { 40c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX) 419566063dSJacob Faibussowitsch PetscCall(PetscStrlcat(file, "hpd-complex-", sizeof(file))); 42c4762a1bSJed Brown #else 439566063dSJacob Faibussowitsch PetscCall(PetscStrlcat(file, "spd-real-", sizeof(file))); 44c4762a1bSJed Brown #endif 45c4762a1bSJed Brown } else { 46c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX) 479566063dSJacob Faibussowitsch PetscCall(PetscStrlcat(file, "nh-complex-", sizeof(file))); 48c4762a1bSJed Brown #else 499566063dSJacob Faibussowitsch PetscCall(PetscStrlcat(file, "ns-real-", sizeof(file))); 50c4762a1bSJed Brown #endif 51c4762a1bSJed Brown } 52c4762a1bSJed Brown #if defined(PETSC_USE_64BIT_INDICES) 539566063dSJacob Faibussowitsch PetscCall(PetscStrlcat(file, "int64-", sizeof(file))); 54c4762a1bSJed Brown #else 559566063dSJacob Faibussowitsch PetscCall(PetscStrlcat(file, "int32-", sizeof(file))); 56c4762a1bSJed Brown #endif 57c4762a1bSJed Brown #if defined(PETSC_USE_REAL_SINGLE) 589566063dSJacob Faibussowitsch PetscCall(PetscStrlcat(file, "float32", sizeof(file))); 59c4762a1bSJed Brown #else 609566063dSJacob Faibussowitsch PetscCall(PetscStrlcat(file, "float64", sizeof(file))); 61c4762a1bSJed Brown #endif 62c4762a1bSJed Brown } 63c4762a1bSJed Brown /* Load matrix A */ 649566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &fd)); 659566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 669566063dSJacob Faibussowitsch PetscCall(MatLoad(A, fd)); 679566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&fd)); 689566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &m, &n)); 6908401ef6SPierre 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); 70c4762a1bSJed Brown 71a5b23f4aSJose E. Roman /* Create dense matrix C and X; C holds true solution with identical columns */ 72c4762a1bSJed Brown nrhs = 2; 739566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-nrhs", &nrhs, NULL)); 749566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &C)); 759566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, m, PETSC_DECIDE, PETSC_DECIDE, nrhs)); 769566063dSJacob Faibussowitsch PetscCall(MatSetType(C, MATDENSE)); 779566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C)); 789566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 79c4762a1bSJed Brown 809566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rand)); 819566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rand)); 829566063dSJacob Faibussowitsch PetscCall(MatSetRandom(C, rand)); 839566063dSJacob Faibussowitsch PetscCall(MatDuplicate(C, MAT_DO_NOT_COPY_VALUES, &X)); 84c4762a1bSJed Brown 85c4762a1bSJed Brown /* Create vectors */ 869566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &x)); 879566063dSJacob Faibussowitsch PetscCall(VecSetSizes(x, n, PETSC_DECIDE)); 889566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(x)); 899566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &b)); 909566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &u)); /* save the true solution */ 91c4762a1bSJed Brown 929566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-solver", &isolver, NULL)); 93c4762a1bSJed Brown switch (isolver) { 94c4762a1bSJed Brown #if defined(PETSC_HAVE_MUMPS) 95d71ae5a4SJacob Faibussowitsch case 0: 96c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(solver, MATSOLVERMUMPS, sizeof(solver))); 97d71ae5a4SJacob Faibussowitsch break; 98c4762a1bSJed Brown #endif 99c4762a1bSJed Brown #if defined(PETSC_HAVE_MKL_PARDISO) 100d71ae5a4SJacob Faibussowitsch case 1: 101c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(solver, MATSOLVERMKL_PARDISO, sizeof(solver))); 102d71ae5a4SJacob Faibussowitsch break; 103c4762a1bSJed Brown #endif 104d71ae5a4SJacob Faibussowitsch default: 105c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(solver, MATSOLVERPETSC, sizeof(solver))); 106d71ae5a4SJacob Faibussowitsch break; 107c4762a1bSJed Brown } 108c4762a1bSJed Brown 109c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX) 110c4762a1bSJed Brown if (isolver == 0 && symm && !data_provided) { /* MUMPS (5.0.0) does not have support for hermitian matrices, so make them symmetric */ 111c4762a1bSJed Brown PetscScalar im = PetscSqrtScalar((PetscScalar)-1.); 112c4762a1bSJed Brown PetscScalar val = -1.0; 113c4762a1bSJed Brown val = val + im; 1149566063dSJacob Faibussowitsch PetscCall(MatSetValue(A, 1, 0, val, INSERT_VALUES)); 1159566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1169566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 117c4762a1bSJed Brown } 118c4762a1bSJed Brown #endif 119c4762a1bSJed Brown 1209566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-schur_ratio", &sratio, NULL)); 121bc30f867SBarry Smith PetscCheck(sratio >= 0. && sratio <= 1., PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Invalid ratio for schur degrees of freedom %g", (double)sratio); 122c4762a1bSJed Brown size_schur = (PetscInt)(sratio * m); 123c4762a1bSJed Brown 1249566063dSJacob 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)); 125c4762a1bSJed Brown 126c4762a1bSJed Brown /* Test LU/Cholesky Factorization */ 127c4762a1bSJed Brown use_lu = PETSC_FALSE; 128c4762a1bSJed Brown if (!symm) use_lu = PETSC_TRUE; 129c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX) 130c4762a1bSJed Brown if (isolver == 1) use_lu = PETSC_TRUE; 131c4762a1bSJed Brown #endif 132c4762a1bSJed Brown if (cuda && symm && !herm) use_lu = PETSC_TRUE; 133c4762a1bSJed Brown 134c4762a1bSJed Brown if (herm && !use_lu) { /* test also conversion routines inside the solver packages */ 1359566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_SYMMETRIC, PETSC_TRUE)); 1369566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATSEQSBAIJ, MAT_INPLACE_MATRIX, &A)); 137c4762a1bSJed Brown } 138c4762a1bSJed Brown 139c4762a1bSJed Brown if (use_lu) { 1409566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A, solver, MAT_FACTOR_LU, &F)); 141c4762a1bSJed Brown } else { 142c4762a1bSJed Brown if (herm) { 1439566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_SPD, PETSC_TRUE)); 144c4762a1bSJed Brown } else { 1459566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_SYMMETRIC, PETSC_TRUE)); 1469566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_SPD, PETSC_FALSE)); 147c4762a1bSJed Brown } 1489566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A, solver, MAT_FACTOR_CHOLESKY, &F)); 149c4762a1bSJed Brown } 1509566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, size_schur, m - size_schur, 1, &is_schur)); 1519566063dSJacob Faibussowitsch PetscCall(MatFactorSetSchurIS(F, is_schur)); 152c4762a1bSJed Brown 1539566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_schur)); 154c4762a1bSJed Brown if (use_lu) { 1559566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(F, A, NULL, NULL, NULL)); 156c4762a1bSJed Brown } else { 1579566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorSymbolic(F, A, NULL, NULL)); 158c4762a1bSJed Brown } 159c4762a1bSJed Brown 160c4762a1bSJed Brown for (nfact = 0; nfact < 3; nfact++) { 161c4762a1bSJed Brown Mat AD; 162c4762a1bSJed Brown 163c4762a1bSJed Brown if (!nfact) { 1649566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rand)); 16548a46eb9SPierre Jolivet if (symm && herm) PetscCall(VecAbs(x)); 1669566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(A, x, ADD_VALUES)); 167c4762a1bSJed Brown } 168c4762a1bSJed Brown if (use_lu) { 1699566063dSJacob Faibussowitsch PetscCall(MatLUFactorNumeric(F, A, NULL)); 170c4762a1bSJed Brown } else { 1719566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorNumeric(F, A, NULL)); 172c4762a1bSJed Brown } 173c4762a1bSJed Brown if (cuda) { 1749566063dSJacob Faibussowitsch PetscCall(MatFactorGetSchurComplement(F, &S, NULL)); 1759566063dSJacob Faibussowitsch PetscCall(MatSetType(S, MATSEQDENSECUDA)); 1769566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(S, &xschur, &bschur)); 1779566063dSJacob Faibussowitsch PetscCall(MatFactorRestoreSchurComplement(F, &S, MAT_FACTOR_SCHUR_UNFACTORED)); 178c4762a1bSJed Brown } 1799566063dSJacob Faibussowitsch PetscCall(MatFactorCreateSchurComplement(F, &S, NULL)); 18048a46eb9SPierre Jolivet if (!cuda) PetscCall(MatCreateVecs(S, &xschur, &bschur)); 1819566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xschur, &uschur)); 18248a46eb9SPierre Jolivet if (nfact == 1 && (!cuda || (herm && symm))) PetscCall(MatFactorInvertSchurComplement(F)); 183c4762a1bSJed Brown for (nsolve = 0; nsolve < 2; nsolve++) { 1849566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rand)); 1859566063dSJacob Faibussowitsch PetscCall(VecCopy(x, u)); 186c4762a1bSJed Brown 187c4762a1bSJed Brown if (nsolve) { 1889566063dSJacob Faibussowitsch PetscCall(MatMult(A, x, b)); 1899566063dSJacob Faibussowitsch PetscCall(MatSolve(F, b, x)); 190c4762a1bSJed Brown } else { 1919566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, x, b)); 1929566063dSJacob Faibussowitsch PetscCall(MatSolveTranspose(F, b, x)); 193c4762a1bSJed Brown } 194c4762a1bSJed Brown /* Check the error */ 1959566063dSJacob Faibussowitsch PetscCall(VecAXPY(u, -1.0, x)); /* u <- (-1.0)x + u */ 1969566063dSJacob Faibussowitsch PetscCall(VecNorm(u, NORM_2, &norm)); 197c4762a1bSJed Brown if (norm > tol) { 198c4762a1bSJed Brown PetscReal resi; 199c4762a1bSJed Brown if (nsolve) { 2009566063dSJacob Faibussowitsch PetscCall(MatMult(A, x, u)); /* u = A*x */ 201c4762a1bSJed Brown } else { 2029566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, x, u)); /* u = A*x */ 203c4762a1bSJed Brown } 2049566063dSJacob Faibussowitsch PetscCall(VecAXPY(u, -1.0, b)); /* u <- (-1.0)b + u */ 2059566063dSJacob Faibussowitsch PetscCall(VecNorm(u, NORM_2, &resi)); 206c4762a1bSJed Brown if (nsolve) { 2079566063dSJacob 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)); 208c4762a1bSJed Brown } else { 2099566063dSJacob 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)); 210c4762a1bSJed Brown } 211c4762a1bSJed Brown } 2129566063dSJacob Faibussowitsch PetscCall(VecSetRandom(xschur, rand)); 2139566063dSJacob Faibussowitsch PetscCall(VecCopy(xschur, uschur)); 214c4762a1bSJed Brown if (nsolve) { 2159566063dSJacob Faibussowitsch PetscCall(MatMult(S, xschur, bschur)); 2169566063dSJacob Faibussowitsch PetscCall(MatFactorSolveSchurComplement(F, bschur, xschur)); 217c4762a1bSJed Brown } else { 2189566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(S, xschur, bschur)); 2199566063dSJacob Faibussowitsch PetscCall(MatFactorSolveSchurComplementTranspose(F, bschur, xschur)); 220c4762a1bSJed Brown } 221c4762a1bSJed Brown /* Check the error */ 2229566063dSJacob Faibussowitsch PetscCall(VecAXPY(uschur, -1.0, xschur)); /* u <- (-1.0)x + u */ 2239566063dSJacob Faibussowitsch PetscCall(VecNorm(uschur, NORM_2, &norm)); 224c4762a1bSJed Brown if (norm > tol) { 225c4762a1bSJed Brown PetscReal resi; 226c4762a1bSJed Brown if (nsolve) { 2279566063dSJacob Faibussowitsch PetscCall(MatMult(S, xschur, uschur)); /* u = A*x */ 228c4762a1bSJed Brown } else { 2299566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(S, xschur, uschur)); /* u = A*x */ 230c4762a1bSJed Brown } 2319566063dSJacob Faibussowitsch PetscCall(VecAXPY(uschur, -1.0, bschur)); /* u <- (-1.0)b + u */ 2329566063dSJacob Faibussowitsch PetscCall(VecNorm(uschur, NORM_2, &resi)); 233c4762a1bSJed Brown if (nsolve) { 2349566063dSJacob 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)); 235c4762a1bSJed Brown } else { 2369566063dSJacob 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)); 237c4762a1bSJed Brown } 238c4762a1bSJed Brown } 239c4762a1bSJed Brown } 2409566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATSEQAIJ, MAT_INITIAL_MATRIX, &AD)); 241c4762a1bSJed Brown if (!nfact) { 2429566063dSJacob Faibussowitsch PetscCall(MatMatMult(AD, C, MAT_INITIAL_MATRIX, 2.0, &RHS)); 243c4762a1bSJed Brown } else { 2449566063dSJacob Faibussowitsch PetscCall(MatMatMult(AD, C, MAT_REUSE_MATRIX, 2.0, &RHS)); 245c4762a1bSJed Brown } 2469566063dSJacob Faibussowitsch PetscCall(MatDestroy(&AD)); 247c4762a1bSJed Brown for (nsolve = 0; nsolve < 2; nsolve++) { 2489566063dSJacob Faibussowitsch PetscCall(MatMatSolve(F, RHS, X)); 249c4762a1bSJed Brown 250c4762a1bSJed Brown /* Check the error */ 2519566063dSJacob Faibussowitsch PetscCall(MatAXPY(X, -1.0, C, SAME_NONZERO_PATTERN)); 2529566063dSJacob Faibussowitsch PetscCall(MatNorm(X, NORM_FROBENIUS, &norm)); 25348a46eb9SPierre 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)); 254*3e5b40d0SPierre Jolivet #if PetscDefined(HAVE_MUMPS) 255*3e5b40d0SPierre Jolivet PetscCall(MatMumpsSetIcntl(F, 26, 1)); 256*3e5b40d0SPierre Jolivet PetscCall(MatMatSolve(F, RHS, X)); 257*3e5b40d0SPierre Jolivet PetscCall(MatMumpsSetIcntl(F, 26, 2)); 258*3e5b40d0SPierre Jolivet PetscCall(MatMatSolve(F, RHS, X)); 259*3e5b40d0SPierre Jolivet PetscCall(MatMumpsSetIcntl(F, 26, -1)); 260*3e5b40d0SPierre Jolivet 261*3e5b40d0SPierre Jolivet /* Check the error */ 262*3e5b40d0SPierre Jolivet PetscCall(MatAXPY(X, -1.0, C, SAME_NONZERO_PATTERN)); 263*3e5b40d0SPierre Jolivet PetscCall(MatNorm(X, NORM_FROBENIUS, &norm)); 264*3e5b40d0SPierre 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)); 265*3e5b40d0SPierre Jolivet #endif 266c4762a1bSJed Brown } 267c4762a1bSJed Brown if (isolver == 0) { 268c4762a1bSJed Brown Mat spRHS, spRHST, RHST; 269c4762a1bSJed Brown 2709566063dSJacob Faibussowitsch PetscCall(MatTranspose(RHS, MAT_INITIAL_MATRIX, &RHST)); 2719566063dSJacob Faibussowitsch PetscCall(MatConvert(RHST, MATSEQAIJ, MAT_INITIAL_MATRIX, &spRHST)); 2729566063dSJacob Faibussowitsch PetscCall(MatCreateTranspose(spRHST, &spRHS)); 273c4762a1bSJed Brown for (nsolve = 0; nsolve < 2; nsolve++) { 2749566063dSJacob Faibussowitsch PetscCall(MatMatSolve(F, spRHS, X)); 275c4762a1bSJed Brown 276c4762a1bSJed Brown /* Check the error */ 2779566063dSJacob Faibussowitsch PetscCall(MatAXPY(X, -1.0, C, SAME_NONZERO_PATTERN)); 2789566063dSJacob Faibussowitsch PetscCall(MatNorm(X, NORM_FROBENIUS, &norm)); 27948a46eb9SPierre 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)); 280c4762a1bSJed Brown } 2819566063dSJacob Faibussowitsch PetscCall(MatDestroy(&spRHST)); 2829566063dSJacob Faibussowitsch PetscCall(MatDestroy(&spRHS)); 2839566063dSJacob Faibussowitsch PetscCall(MatDestroy(&RHST)); 284c4762a1bSJed Brown } 2859566063dSJacob Faibussowitsch PetscCall(MatDestroy(&S)); 2869566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xschur)); 2879566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bschur)); 2889566063dSJacob Faibussowitsch PetscCall(VecDestroy(&uschur)); 289c4762a1bSJed Brown } 290c4762a1bSJed Brown /* Free data structures */ 2919566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 2929566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 2939566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F)); 2949566063dSJacob Faibussowitsch PetscCall(MatDestroy(&X)); 2959566063dSJacob Faibussowitsch PetscCall(MatDestroy(&RHS)); 2969566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand)); 2979566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 2989566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 2999566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 3009566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 301b122ec5aSJacob Faibussowitsch return 0; 302c4762a1bSJed Brown } 303c4762a1bSJed Brown 304c4762a1bSJed Brown /*TEST 305c4762a1bSJed Brown 306c4762a1bSJed Brown testset: 30738cf7977SStefano Zampini requires: mkl_pardiso double !complex 308c4762a1bSJed Brown args: -solver 1 309c4762a1bSJed Brown 310c4762a1bSJed Brown test: 311c4762a1bSJed Brown suffix: mkl_pardiso 312c4762a1bSJed Brown test: 313c4762a1bSJed Brown requires: cuda 314c4762a1bSJed Brown suffix: mkl_pardiso_cuda 315c4762a1bSJed Brown args: -cuda_solve 316c4762a1bSJed Brown output_file: output/ex192_mkl_pardiso.out 317c4762a1bSJed Brown test: 318c4762a1bSJed Brown suffix: mkl_pardiso_1 319c4762a1bSJed Brown args: -symmetric_solve 320c4762a1bSJed Brown output_file: output/ex192_mkl_pardiso_1.out 321c4762a1bSJed Brown test: 322c4762a1bSJed Brown requires: cuda 323c4762a1bSJed Brown suffix: mkl_pardiso_cuda_1 324c4762a1bSJed Brown args: -symmetric_solve -cuda_solve 325c4762a1bSJed Brown output_file: output/ex192_mkl_pardiso_1.out 326c4762a1bSJed Brown test: 327c4762a1bSJed Brown suffix: mkl_pardiso_3 328c4762a1bSJed Brown args: -symmetric_solve -hermitian_solve 329c4762a1bSJed Brown output_file: output/ex192_mkl_pardiso_3.out 330c4762a1bSJed Brown test: 331dfd57a17SPierre Jolivet requires: cuda defined(PETSC_HAVE_CUSOLVERDNDPOTRI) 332c4762a1bSJed Brown suffix: mkl_pardiso_cuda_3 333c4762a1bSJed Brown args: -symmetric_solve -hermitian_solve -cuda_solve 334c4762a1bSJed Brown output_file: output/ex192_mkl_pardiso_3.out 335c4762a1bSJed Brown 336c4762a1bSJed Brown testset: 337c4762a1bSJed Brown requires: mumps double !complex 338c4762a1bSJed Brown args: -solver 0 339c4762a1bSJed Brown 340c4762a1bSJed Brown test: 341c4762a1bSJed Brown suffix: mumps 342c4762a1bSJed Brown test: 343c4762a1bSJed Brown requires: cuda 344c4762a1bSJed Brown suffix: mumps_cuda 345c4762a1bSJed Brown args: -cuda_solve 346c4762a1bSJed Brown output_file: output/ex192_mumps.out 347c4762a1bSJed Brown test: 348c4762a1bSJed Brown suffix: mumps_2 349c4762a1bSJed Brown args: -symmetric_solve 350c4762a1bSJed Brown output_file: output/ex192_mumps_2.out 351c4762a1bSJed Brown test: 352c4762a1bSJed Brown requires: cuda 353c4762a1bSJed Brown suffix: mumps_cuda_2 354c4762a1bSJed Brown args: -symmetric_solve -cuda_solve 355c4762a1bSJed Brown output_file: output/ex192_mumps_2.out 356c4762a1bSJed Brown test: 357c4762a1bSJed Brown suffix: mumps_3 358c4762a1bSJed Brown args: -symmetric_solve -hermitian_solve 359c4762a1bSJed Brown output_file: output/ex192_mumps_3.out 360c4762a1bSJed Brown test: 361dfd57a17SPierre Jolivet requires: cuda defined(PETSC_HAVE_CUSOLVERDNDPOTRI) 362c4762a1bSJed Brown suffix: mumps_cuda_3 363c4762a1bSJed Brown args: -symmetric_solve -hermitian_solve -cuda_solve 364c4762a1bSJed Brown output_file: output/ex192_mumps_3.out 365c4762a1bSJed Brown 366c4762a1bSJed Brown TEST*/ 367