1c4762a1bSJed Brown static char help[] = "Tests ILU and ICC factorization with and without matrix ordering on seqaij format, and illustrates drawing of matrix sparsity structure with MatView().\n\ 2c4762a1bSJed Brown Input parameters are:\n\ 3c4762a1bSJed Brown -lf <level> : level of fill for ILU (default is 0)\n\ 4c4762a1bSJed Brown -lu : use full LU or Cholesky factorization\n\ 5c4762a1bSJed Brown -m <value>,-n <value> : grid dimensions\n\ 6c4762a1bSJed Brown Note that most users should employ the KSP interface to the\n\ 7c4762a1bSJed Brown linear solvers instead of using the factorization routines\n\ 8c4762a1bSJed Brown directly.\n\n"; 9c4762a1bSJed Brown 10c4762a1bSJed Brown #include <petscmat.h> 11c4762a1bSJed Brown 12d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 13d71ae5a4SJacob Faibussowitsch { 14c4762a1bSJed Brown Mat C, A; 15c4762a1bSJed Brown PetscInt i, j, m = 5, n = 5, Ii, J, lf = 0; 16bd5dbebeSNils Friess PetscBool LU = PETSC_FALSE, CHOLESKY = PETSC_FALSE, TRIANGULAR = PETSC_FALSE, MATDSPL = PETSC_FALSE, flg, matordering, use_mkl_pardiso = PETSC_FALSE; 17c4762a1bSJed Brown PetscScalar v; 18c4762a1bSJed Brown IS row, col; 19c4762a1bSJed Brown PetscViewer viewer1, viewer2; 20c4762a1bSJed Brown MatFactorInfo info; 21c4762a1bSJed Brown Vec x, y, b, ytmp; 22c4762a1bSJed Brown PetscReal norm2, norm2_inplace, tol = 100. * PETSC_MACHINE_EPSILON; 23c4762a1bSJed Brown PetscRandom rdm; 24c4762a1bSJed Brown PetscMPIInt size; 25bd5dbebeSNils Friess char pack[PETSC_MAX_PATH_LEN]; 26c4762a1bSJed Brown 27327415f7SBarry Smith PetscFunctionBeginUser; 28*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help)); 299566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 30be096a46SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 319566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL)); 329566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL)); 339566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-lf", &lf, NULL)); 34c4762a1bSJed Brown 359566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawOpen(PETSC_COMM_SELF, 0, 0, 0, 0, 400, 400, &viewer1)); 369566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawOpen(PETSC_COMM_SELF, 0, 0, 400, 0, 400, 400, &viewer2)); 37c4762a1bSJed Brown 389566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &C)); 399566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, m * n, m * n, m * n, m * n)); 409566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C)); 419566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 42c4762a1bSJed Brown 43c4762a1bSJed Brown /* Create matrix C in seqaij format and sC in seqsbaij. (This is five-point stencil with some extra elements) */ 44c4762a1bSJed Brown for (i = 0; i < m; i++) { 45c4762a1bSJed Brown for (j = 0; j < n; j++) { 469371c9d4SSatish Balay v = -1.0; 479371c9d4SSatish Balay Ii = j + n * i; 489371c9d4SSatish Balay J = Ii - n; 499371c9d4SSatish Balay if (J >= 0) PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 509371c9d4SSatish Balay J = Ii + n; 519371c9d4SSatish Balay if (J < m * n) PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 529371c9d4SSatish Balay J = Ii - 1; 539371c9d4SSatish Balay if (J >= 0) PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 549371c9d4SSatish Balay J = Ii + 1; 559371c9d4SSatish Balay if (J < m * n) PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 569371c9d4SSatish Balay v = 4.0; 579371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, INSERT_VALUES)); 58c4762a1bSJed Brown } 59c4762a1bSJed Brown } 609566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 619566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 62c4762a1bSJed Brown 639566063dSJacob Faibussowitsch PetscCall(MatIsSymmetric(C, 0.0, &flg)); 6428b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_SUP, "C is non-symmetric"); 65c4762a1bSJed Brown 66bd5dbebeSNils Friess PetscCall(MatSetOption(C, MAT_SPD, PETSC_TRUE)); 67bd5dbebeSNils Friess 68c4762a1bSJed Brown /* Create vectors for error checking */ 699566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(C, &x, &b)); 709566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &y)); 719566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &ytmp)); 729566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rdm)); 739566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rdm)); 749566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rdm)); 759566063dSJacob Faibussowitsch PetscCall(MatMult(C, x, b)); 76c4762a1bSJed Brown 779566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-mat_ordering", &matordering)); 78c4762a1bSJed Brown if (matordering) { 799566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(C, MATORDERINGRCM, &row, &col)); 80c4762a1bSJed Brown } else { 819566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(C, MATORDERINGNATURAL, &row, &col)); 82c4762a1bSJed Brown } 83c4762a1bSJed Brown 849566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-display_matrices", &MATDSPL)); 85c4762a1bSJed Brown if (MATDSPL) { 86c4762a1bSJed Brown printf("original matrix:\n"); 879566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_SELF, PETSC_VIEWER_ASCII_INFO)); 889566063dSJacob Faibussowitsch PetscCall(MatView(C, PETSC_VIEWER_STDOUT_SELF)); 899566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_SELF)); 909566063dSJacob Faibussowitsch PetscCall(MatView(C, PETSC_VIEWER_STDOUT_SELF)); 919566063dSJacob Faibussowitsch PetscCall(MatView(C, viewer1)); 92c4762a1bSJed Brown } 93c4762a1bSJed Brown 94c4762a1bSJed Brown /* Compute LU or ILU factor A */ 959566063dSJacob Faibussowitsch PetscCall(MatFactorInfoInitialize(&info)); 96c4762a1bSJed Brown 97c4762a1bSJed Brown info.fill = 1.0; 98c4762a1bSJed Brown info.diagonal_fill = 0; 99c4762a1bSJed Brown info.zeropivot = 0.0; 100c4762a1bSJed Brown 101bd5dbebeSNils Friess PetscCall(PetscOptionsGetString(NULL, NULL, "-mat_solver_type", pack, sizeof(pack), NULL)); 102bd5dbebeSNils Friess #if defined(PETSC_HAVE_MKL_PARDISO) 103bd5dbebeSNils Friess PetscCall(PetscStrcmp(MATSOLVERMKL_PARDISO, pack, &use_mkl_pardiso)); 104bd5dbebeSNils Friess #endif 105bd5dbebeSNils Friess 1069566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-lu", &LU)); 107c4762a1bSJed Brown if (LU) { 108c4762a1bSJed Brown printf("Test LU...\n"); 109bd5dbebeSNils Friess if (use_mkl_pardiso) PetscCall(MatGetFactor(C, MATSOLVERMKL_PARDISO, MAT_FACTOR_LU, &A)); 110bd5dbebeSNils Friess else PetscCall(MatGetFactor(C, MATSOLVERPETSC, MAT_FACTOR_LU, &A)); 1119566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(A, C, row, col, &info)); 112c4762a1bSJed Brown } else { 113c4762a1bSJed Brown printf("Test ILU...\n"); 114c4762a1bSJed Brown info.levels = lf; 115c4762a1bSJed Brown 1169566063dSJacob Faibussowitsch PetscCall(MatGetFactor(C, MATSOLVERPETSC, MAT_FACTOR_ILU, &A)); 1179566063dSJacob Faibussowitsch PetscCall(MatILUFactorSymbolic(A, C, row, col, &info)); 118c4762a1bSJed Brown } 1199566063dSJacob Faibussowitsch PetscCall(MatLUFactorNumeric(A, C, &info)); 120c4762a1bSJed Brown 121bd5dbebeSNils Friess /* test MatForwardSolve() and MatBackwardSolve() with MKL Pardiso*/ 122bd5dbebeSNils Friess if (LU && use_mkl_pardiso) { 123bd5dbebeSNils Friess PetscCall(PetscOptionsHasName(NULL, NULL, "-triangular_solve", &TRIANGULAR)); 124bd5dbebeSNils Friess if (TRIANGULAR) { 125bd5dbebeSNils Friess printf("Test MatForwardSolve...\n"); 126bd5dbebeSNils Friess PetscCall(MatForwardSolve(A, b, ytmp)); 127bd5dbebeSNils Friess printf("Test MatBackwardSolve...\n"); 128bd5dbebeSNils Friess PetscCall(MatBackwardSolve(A, ytmp, y)); 129bd5dbebeSNils Friess PetscCall(VecAXPY(y, -1.0, x)); 130bd5dbebeSNils Friess PetscCall(VecNorm(y, NORM_2, &norm2)); 131bd5dbebeSNils Friess if (norm2 > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "MatForwardSolve and BackwardSolve: Norm of error=%g\n", (double)norm2)); 132bd5dbebeSNils Friess } 133bd5dbebeSNils Friess } 134bd5dbebeSNils Friess 135c4762a1bSJed Brown /* Solve A*y = b, then check the error */ 1369566063dSJacob Faibussowitsch PetscCall(MatSolve(A, b, y)); 1379566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, x)); 1389566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm2)); 1399566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 140c4762a1bSJed Brown 141c4762a1bSJed Brown /* Test in-place ILU(0) and compare it with the out-place ILU(0) */ 142c4762a1bSJed Brown if (!LU && lf == 0) { 1439566063dSJacob Faibussowitsch PetscCall(MatDuplicate(C, MAT_COPY_VALUES, &A)); 1449566063dSJacob Faibussowitsch PetscCall(MatILUFactor(A, row, col, &info)); 145c4762a1bSJed Brown /* 146c4762a1bSJed Brown printf("In-place factored matrix:\n"); 1479566063dSJacob Faibussowitsch PetscCall(MatView(C,PETSC_VIEWER_STDOUT_SELF)); 148c4762a1bSJed Brown */ 1499566063dSJacob Faibussowitsch PetscCall(MatSolve(A, b, y)); 1509566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, x)); 1519566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm2_inplace)); 15208401ef6SPierre Jolivet PetscCheck(PetscAbs(norm2 - norm2_inplace) <= tol, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ILU(0) %g and in-place ILU(0) %g give different residuals", (double)norm2, (double)norm2_inplace); 1539566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 154c4762a1bSJed Brown } 155c4762a1bSJed Brown 156c4762a1bSJed Brown /* Test Cholesky and ICC on seqaij matrix with matrix reordering on aij matrix C */ 157bd5dbebeSNils Friess PetscCall(PetscOptionsHasName(NULL, NULL, "-chol", &CHOLESKY)); 158c4762a1bSJed Brown if (CHOLESKY) { 159c4762a1bSJed Brown printf("Test Cholesky...\n"); 160c4762a1bSJed Brown lf = -1; 161bd5dbebeSNils Friess if (use_mkl_pardiso) PetscCall(MatGetFactor(C, MATSOLVERMKL_PARDISO, MAT_FACTOR_CHOLESKY, &A)); 162bd5dbebeSNils Friess else PetscCall(MatGetFactor(C, MATSOLVERPETSC, MAT_FACTOR_CHOLESKY, &A)); 1639566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorSymbolic(A, C, row, &info)); 164c4762a1bSJed Brown } else { 165c4762a1bSJed Brown printf("Test ICC...\n"); 166c4762a1bSJed Brown info.levels = lf; 167c4762a1bSJed Brown info.fill = 1.0; 168c4762a1bSJed Brown info.diagonal_fill = 0; 169c4762a1bSJed Brown info.zeropivot = 0.0; 170c4762a1bSJed Brown 1719566063dSJacob Faibussowitsch PetscCall(MatGetFactor(C, MATSOLVERPETSC, MAT_FACTOR_ICC, &A)); 1729566063dSJacob Faibussowitsch PetscCall(MatICCFactorSymbolic(A, C, row, &info)); 173c4762a1bSJed Brown } 1749566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorNumeric(A, C, &info)); 175c4762a1bSJed Brown 176c4762a1bSJed Brown /* test MatForwardSolve() and MatBackwardSolve() with matrix reordering on aij matrix C */ 177bd5dbebeSNils Friess if (CHOLESKY) { 1789566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-triangular_solve", &TRIANGULAR)); 179c4762a1bSJed Brown if (TRIANGULAR) { 180c4762a1bSJed Brown printf("Test MatForwardSolve...\n"); 1819566063dSJacob Faibussowitsch PetscCall(MatForwardSolve(A, b, ytmp)); 182c4762a1bSJed Brown printf("Test MatBackwardSolve...\n"); 1839566063dSJacob Faibussowitsch PetscCall(MatBackwardSolve(A, ytmp, y)); 1849566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, x)); 1859566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm2)); 18648a46eb9SPierre Jolivet if (norm2 > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "MatForwardSolve and BackwardSolve: Norm of error=%g\n", (double)norm2)); 187c4762a1bSJed Brown } 188c4762a1bSJed Brown } 189c4762a1bSJed Brown 1909566063dSJacob Faibussowitsch PetscCall(MatSolve(A, b, y)); 1919566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1929566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, x)); 1939566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm2)); 19448a46eb9SPierre Jolivet if (lf == -1 && norm2 > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, " reordered SEQAIJ: Cholesky/ICC levels %" PetscInt_FMT ", residual %g\n", lf, (double)norm2)); 195c4762a1bSJed Brown 196c4762a1bSJed Brown /* Test in-place ICC(0) and compare it with the out-place ICC(0) */ 197c4762a1bSJed Brown if (!CHOLESKY && lf == 0 && !matordering) { 1989566063dSJacob Faibussowitsch PetscCall(MatConvert(C, MATSBAIJ, MAT_INITIAL_MATRIX, &A)); 1999566063dSJacob Faibussowitsch PetscCall(MatICCFactor(A, row, &info)); 200c4762a1bSJed Brown /* 201c4762a1bSJed Brown printf("In-place factored matrix:\n"); 2029566063dSJacob Faibussowitsch PetscCall(MatView(A,PETSC_VIEWER_STDOUT_SELF)); 203c4762a1bSJed Brown */ 2049566063dSJacob Faibussowitsch PetscCall(MatSolve(A, b, y)); 2059566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, x)); 2069566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm2_inplace)); 20708401ef6SPierre Jolivet PetscCheck(PetscAbs(norm2 - norm2_inplace) <= tol, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ICC(0) %g and in-place ICC(0) %g give different residuals", (double)norm2, (double)norm2_inplace); 2089566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 209c4762a1bSJed Brown } 210c4762a1bSJed Brown 211c4762a1bSJed Brown /* Free data structures */ 2129566063dSJacob Faibussowitsch PetscCall(ISDestroy(&row)); 2139566063dSJacob Faibussowitsch PetscCall(ISDestroy(&col)); 2149566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 2159566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer1)); 2169566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer2)); 2179566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm)); 2189566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 2199566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 2209566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ytmp)); 2219566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 2229566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 223b122ec5aSJacob Faibussowitsch return 0; 224c4762a1bSJed Brown } 225c4762a1bSJed Brown 226c4762a1bSJed Brown /*TEST 227c4762a1bSJed Brown 228c4762a1bSJed Brown test: 229c4762a1bSJed Brown args: -mat_ordering -display_matrices -nox 2308cc725e6SPierre Jolivet filter: grep -v " MPI process" 231c4762a1bSJed Brown 232c4762a1bSJed Brown test: 233c4762a1bSJed Brown suffix: 2 234bd5dbebeSNils Friess args: -mat_ordering -display_matrices -nox -lu -chol 235c4762a1bSJed Brown 236c4762a1bSJed Brown test: 237c4762a1bSJed Brown suffix: 3 238bd5dbebeSNils Friess args: -mat_ordering -lu -chol -triangular_solve 239c4762a1bSJed Brown 240c4762a1bSJed Brown test: 241c4762a1bSJed Brown suffix: 4 242c4762a1bSJed Brown 243c4762a1bSJed Brown test: 244c4762a1bSJed Brown suffix: 5 245bd5dbebeSNils Friess args: -lu -chol 246c4762a1bSJed Brown 247c4762a1bSJed Brown test: 248c4762a1bSJed Brown suffix: 6 249bd5dbebeSNils Friess args: -lu -chol -triangular_solve 250c4762a1bSJed Brown output_file: output/ex30_3.out 251c4762a1bSJed Brown 252bd5dbebeSNils Friess test: 253bd5dbebeSNils Friess suffix: 7 254bd5dbebeSNils Friess requires: mkl_pardiso 255bd5dbebeSNils Friess args: -lu -chol -mat_solver_type mkl_pardiso 256bd5dbebeSNils Friess output_file: output/ex30_5.out 257bd5dbebeSNils Friess 258bd5dbebeSNils Friess test: 259bd5dbebeSNils Friess suffix: 8 260bd5dbebeSNils Friess requires: mkl_pardiso 261bd5dbebeSNils Friess args: -lu -mat_solver_type mkl_pardiso -triangular_solve 262bd5dbebeSNils Friess 263c4762a1bSJed Brown TEST*/ 264