1c4762a1bSJed Brown 2c4762a1bSJed 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\ 3c4762a1bSJed Brown Input parameters are:\n\ 4c4762a1bSJed Brown -lf <level> : level of fill for ILU (default is 0)\n\ 5c4762a1bSJed Brown -lu : use full LU or Cholesky factorization\n\ 6c4762a1bSJed Brown -m <value>,-n <value> : grid dimensions\n\ 7c4762a1bSJed Brown Note that most users should employ the KSP interface to the\n\ 8c4762a1bSJed Brown linear solvers instead of using the factorization routines\n\ 9c4762a1bSJed Brown directly.\n\n"; 10c4762a1bSJed Brown 11c4762a1bSJed Brown #include <petscmat.h> 12c4762a1bSJed Brown 13*d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 14*d71ae5a4SJacob Faibussowitsch { 15c4762a1bSJed Brown Mat C, A; 16c4762a1bSJed Brown PetscInt i, j, m = 5, n = 5, Ii, J, lf = 0; 17c4762a1bSJed Brown PetscBool LU = PETSC_FALSE, CHOLESKY, TRIANGULAR = PETSC_FALSE, MATDSPL = PETSC_FALSE, flg, matordering; 18c4762a1bSJed Brown PetscScalar v; 19c4762a1bSJed Brown IS row, col; 20c4762a1bSJed Brown PetscViewer viewer1, viewer2; 21c4762a1bSJed Brown MatFactorInfo info; 22c4762a1bSJed Brown Vec x, y, b, ytmp; 23c4762a1bSJed Brown PetscReal norm2, norm2_inplace, tol = 100. * PETSC_MACHINE_EPSILON; 24c4762a1bSJed Brown PetscRandom rdm; 25c4762a1bSJed Brown PetscMPIInt size; 26c4762a1bSJed Brown 27327415f7SBarry Smith PetscFunctionBeginUser; 289566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, 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 66c4762a1bSJed Brown /* Create vectors for error checking */ 679566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(C, &x, &b)); 689566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &y)); 699566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &ytmp)); 709566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rdm)); 719566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rdm)); 729566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rdm)); 739566063dSJacob Faibussowitsch PetscCall(MatMult(C, x, b)); 74c4762a1bSJed Brown 759566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-mat_ordering", &matordering)); 76c4762a1bSJed Brown if (matordering) { 779566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(C, MATORDERINGRCM, &row, &col)); 78c4762a1bSJed Brown } else { 799566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(C, MATORDERINGNATURAL, &row, &col)); 80c4762a1bSJed Brown } 81c4762a1bSJed Brown 829566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-display_matrices", &MATDSPL)); 83c4762a1bSJed Brown if (MATDSPL) { 84c4762a1bSJed Brown printf("original matrix:\n"); 859566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_SELF, PETSC_VIEWER_ASCII_INFO)); 869566063dSJacob Faibussowitsch PetscCall(MatView(C, PETSC_VIEWER_STDOUT_SELF)); 879566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_SELF)); 889566063dSJacob Faibussowitsch PetscCall(MatView(C, PETSC_VIEWER_STDOUT_SELF)); 899566063dSJacob Faibussowitsch PetscCall(MatView(C, viewer1)); 90c4762a1bSJed Brown } 91c4762a1bSJed Brown 92c4762a1bSJed Brown /* Compute LU or ILU factor A */ 939566063dSJacob Faibussowitsch PetscCall(MatFactorInfoInitialize(&info)); 94c4762a1bSJed Brown 95c4762a1bSJed Brown info.fill = 1.0; 96c4762a1bSJed Brown info.diagonal_fill = 0; 97c4762a1bSJed Brown info.zeropivot = 0.0; 98c4762a1bSJed Brown 999566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-lu", &LU)); 100c4762a1bSJed Brown if (LU) { 101c4762a1bSJed Brown printf("Test LU...\n"); 1029566063dSJacob Faibussowitsch PetscCall(MatGetFactor(C, MATSOLVERPETSC, MAT_FACTOR_LU, &A)); 1039566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(A, C, row, col, &info)); 104c4762a1bSJed Brown } else { 105c4762a1bSJed Brown printf("Test ILU...\n"); 106c4762a1bSJed Brown info.levels = lf; 107c4762a1bSJed Brown 1089566063dSJacob Faibussowitsch PetscCall(MatGetFactor(C, MATSOLVERPETSC, MAT_FACTOR_ILU, &A)); 1099566063dSJacob Faibussowitsch PetscCall(MatILUFactorSymbolic(A, C, row, col, &info)); 110c4762a1bSJed Brown } 1119566063dSJacob Faibussowitsch PetscCall(MatLUFactorNumeric(A, C, &info)); 112c4762a1bSJed Brown 113c4762a1bSJed Brown /* Solve A*y = b, then check the error */ 1149566063dSJacob Faibussowitsch PetscCall(MatSolve(A, b, y)); 1159566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, x)); 1169566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm2)); 1179566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 118c4762a1bSJed Brown 119c4762a1bSJed Brown /* Test in-place ILU(0) and compare it with the out-place ILU(0) */ 120c4762a1bSJed Brown if (!LU && lf == 0) { 1219566063dSJacob Faibussowitsch PetscCall(MatDuplicate(C, MAT_COPY_VALUES, &A)); 1229566063dSJacob Faibussowitsch PetscCall(MatILUFactor(A, row, col, &info)); 123c4762a1bSJed Brown /* 124c4762a1bSJed Brown printf("In-place factored matrix:\n"); 1259566063dSJacob Faibussowitsch PetscCall(MatView(C,PETSC_VIEWER_STDOUT_SELF)); 126c4762a1bSJed Brown */ 1279566063dSJacob Faibussowitsch PetscCall(MatSolve(A, b, y)); 1289566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, x)); 1299566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm2_inplace)); 13008401ef6SPierre 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); 1319566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 132c4762a1bSJed Brown } 133c4762a1bSJed Brown 134c4762a1bSJed Brown /* Test Cholesky and ICC on seqaij matrix with matrix reordering on aij matrix C */ 135c4762a1bSJed Brown CHOLESKY = LU; 136c4762a1bSJed Brown if (CHOLESKY) { 137c4762a1bSJed Brown printf("Test Cholesky...\n"); 138c4762a1bSJed Brown lf = -1; 1399566063dSJacob Faibussowitsch PetscCall(MatGetFactor(C, MATSOLVERPETSC, MAT_FACTOR_CHOLESKY, &A)); 1409566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorSymbolic(A, C, row, &info)); 141c4762a1bSJed Brown } else { 142c4762a1bSJed Brown printf("Test ICC...\n"); 143c4762a1bSJed Brown info.levels = lf; 144c4762a1bSJed Brown info.fill = 1.0; 145c4762a1bSJed Brown info.diagonal_fill = 0; 146c4762a1bSJed Brown info.zeropivot = 0.0; 147c4762a1bSJed Brown 1489566063dSJacob Faibussowitsch PetscCall(MatGetFactor(C, MATSOLVERPETSC, MAT_FACTOR_ICC, &A)); 1499566063dSJacob Faibussowitsch PetscCall(MatICCFactorSymbolic(A, C, row, &info)); 150c4762a1bSJed Brown } 1519566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorNumeric(A, C, &info)); 152c4762a1bSJed Brown 153c4762a1bSJed Brown /* test MatForwardSolve() and MatBackwardSolve() with matrix reordering on aij matrix C */ 154c4762a1bSJed Brown if (lf == -1) { 1559566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-triangular_solve", &TRIANGULAR)); 156c4762a1bSJed Brown if (TRIANGULAR) { 157c4762a1bSJed Brown printf("Test MatForwardSolve...\n"); 1589566063dSJacob Faibussowitsch PetscCall(MatForwardSolve(A, b, ytmp)); 159c4762a1bSJed Brown printf("Test MatBackwardSolve...\n"); 1609566063dSJacob Faibussowitsch PetscCall(MatBackwardSolve(A, ytmp, y)); 1619566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, x)); 1629566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm2)); 16348a46eb9SPierre Jolivet if (norm2 > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "MatForwardSolve and BackwardSolve: Norm of error=%g\n", (double)norm2)); 164c4762a1bSJed Brown } 165c4762a1bSJed Brown } 166c4762a1bSJed Brown 1679566063dSJacob Faibussowitsch PetscCall(MatSolve(A, b, y)); 1689566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1699566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, x)); 1709566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm2)); 17148a46eb9SPierre Jolivet if (lf == -1 && norm2 > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, " reordered SEQAIJ: Cholesky/ICC levels %" PetscInt_FMT ", residual %g\n", lf, (double)norm2)); 172c4762a1bSJed Brown 173c4762a1bSJed Brown /* Test in-place ICC(0) and compare it with the out-place ICC(0) */ 174c4762a1bSJed Brown if (!CHOLESKY && lf == 0 && !matordering) { 1759566063dSJacob Faibussowitsch PetscCall(MatConvert(C, MATSBAIJ, MAT_INITIAL_MATRIX, &A)); 1769566063dSJacob Faibussowitsch PetscCall(MatICCFactor(A, row, &info)); 177c4762a1bSJed Brown /* 178c4762a1bSJed Brown printf("In-place factored matrix:\n"); 1799566063dSJacob Faibussowitsch PetscCall(MatView(A,PETSC_VIEWER_STDOUT_SELF)); 180c4762a1bSJed Brown */ 1819566063dSJacob Faibussowitsch PetscCall(MatSolve(A, b, y)); 1829566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, x)); 1839566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm2_inplace)); 18408401ef6SPierre 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); 1859566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 186c4762a1bSJed Brown } 187c4762a1bSJed Brown 188c4762a1bSJed Brown /* Free data structures */ 1899566063dSJacob Faibussowitsch PetscCall(ISDestroy(&row)); 1909566063dSJacob Faibussowitsch PetscCall(ISDestroy(&col)); 1919566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 1929566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer1)); 1939566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer2)); 1949566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm)); 1959566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1969566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 1979566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ytmp)); 1989566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 1999566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 200b122ec5aSJacob Faibussowitsch return 0; 201c4762a1bSJed Brown } 202c4762a1bSJed Brown 203c4762a1bSJed Brown /*TEST 204c4762a1bSJed Brown 205c4762a1bSJed Brown test: 206c4762a1bSJed Brown args: -mat_ordering -display_matrices -nox 2078cc725e6SPierre Jolivet filter: grep -v " MPI process" 208c4762a1bSJed Brown 209c4762a1bSJed Brown test: 210c4762a1bSJed Brown suffix: 2 211c4762a1bSJed Brown args: -mat_ordering -display_matrices -nox -lu 212c4762a1bSJed Brown 213c4762a1bSJed Brown test: 214c4762a1bSJed Brown suffix: 3 215c4762a1bSJed Brown args: -mat_ordering -lu -triangular_solve 216c4762a1bSJed Brown 217c4762a1bSJed Brown test: 218c4762a1bSJed Brown suffix: 4 219c4762a1bSJed Brown 220c4762a1bSJed Brown test: 221c4762a1bSJed Brown suffix: 5 222c4762a1bSJed Brown args: -lu 223c4762a1bSJed Brown 224c4762a1bSJed Brown test: 225c4762a1bSJed Brown suffix: 6 226c4762a1bSJed Brown args: -lu -triangular_solve 227c4762a1bSJed Brown output_file: output/ex30_3.out 228c4762a1bSJed Brown 229c4762a1bSJed Brown TEST*/ 230