1c4762a1bSJed Brown static char help[] = "Tests ILU and ICC factorization with and without matrix ordering on seqsbaij format. Modified from ex30.c\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, sC, sA; 15c4762a1bSJed Brown PetscInt i, j, m = 5, n = 5, Ii, J, lf = 0; 16c4762a1bSJed Brown PetscBool CHOLESKY = PETSC_FALSE, TRIANGULAR = PETSC_FALSE, flg; 17c4762a1bSJed Brown PetscScalar v; 18c4762a1bSJed Brown IS row, col; 19c4762a1bSJed Brown MatFactorInfo info; 20c4762a1bSJed Brown Vec x, y, b, ytmp; 21c4762a1bSJed Brown PetscReal norm2, tol = 100 * PETSC_MACHINE_EPSILON; 22c4762a1bSJed Brown PetscRandom rdm; 23c4762a1bSJed Brown PetscMPIInt size; 24c4762a1bSJed Brown 25327415f7SBarry Smith PetscFunctionBeginUser; 26*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help)); 279566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 28be096a46SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 299566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL)); 309566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL)); 319566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-lf", &lf, NULL)); 32c4762a1bSJed Brown 339566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &C)); 349566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, m * n, m * n, m * n, m * n)); 359566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C)); 369566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 37c4762a1bSJed Brown 38c4762a1bSJed Brown /* Create matrix C in seqaij format and sC in seqsbaij. (This is five-point stencil with some extra elements) */ 39c4762a1bSJed Brown for (i = 0; i < m; i++) { 40c4762a1bSJed Brown for (j = 0; j < n; j++) { 419371c9d4SSatish Balay v = -1.0; 429371c9d4SSatish Balay Ii = j + n * i; 439371c9d4SSatish Balay J = Ii - n; 449371c9d4SSatish Balay if (J >= 0) PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 459371c9d4SSatish Balay J = Ii + n; 469371c9d4SSatish Balay if (J < m * n) PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 479371c9d4SSatish Balay J = Ii - 1; 489371c9d4SSatish Balay if (J >= 0) PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 499371c9d4SSatish Balay J = Ii + 1; 509371c9d4SSatish Balay if (J < m * n) PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 519371c9d4SSatish Balay v = 4.0; 529371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, INSERT_VALUES)); 53c4762a1bSJed Brown } 54c4762a1bSJed Brown } 559566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 569566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 57c4762a1bSJed Brown 589566063dSJacob Faibussowitsch PetscCall(MatIsSymmetric(C, 0.0, &flg)); 5928b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_SUP, "C is non-symmetric"); 609566063dSJacob Faibussowitsch PetscCall(MatConvert(C, MATSEQSBAIJ, MAT_INITIAL_MATRIX, &sC)); 61c4762a1bSJed Brown 62c4762a1bSJed Brown /* Create vectors for error checking */ 639566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(C, &x, &b)); 649566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &y)); 659566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &ytmp)); 669566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rdm)); 679566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rdm)); 689566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rdm)); 699566063dSJacob Faibussowitsch PetscCall(MatMult(C, x, b)); 70c4762a1bSJed Brown 719566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(C, MATORDERINGNATURAL, &row, &col)); 72c4762a1bSJed Brown 73c4762a1bSJed Brown /* Compute CHOLESKY or ICC factor sA */ 749566063dSJacob Faibussowitsch PetscCall(MatFactorInfoInitialize(&info)); 75c4762a1bSJed Brown 76c4762a1bSJed Brown info.fill = 1.0; 77c4762a1bSJed Brown info.diagonal_fill = 0; 78c4762a1bSJed Brown info.zeropivot = 0.0; 79c4762a1bSJed Brown 809566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-cholesky", &CHOLESKY)); 81c4762a1bSJed Brown if (CHOLESKY) { 829566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Test CHOLESKY...\n")); 839566063dSJacob Faibussowitsch PetscCall(MatGetFactor(sC, MATSOLVERPETSC, MAT_FACTOR_CHOLESKY, &sA)); 849566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorSymbolic(sA, sC, row, &info)); 85c4762a1bSJed Brown } else { 869566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Test ICC...\n")); 87c4762a1bSJed Brown info.levels = lf; 88c4762a1bSJed Brown 899566063dSJacob Faibussowitsch PetscCall(MatGetFactor(sC, MATSOLVERPETSC, MAT_FACTOR_ICC, &sA)); 909566063dSJacob Faibussowitsch PetscCall(MatICCFactorSymbolic(sA, sC, row, &info)); 91c4762a1bSJed Brown } 929566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorNumeric(sA, sC, &info)); 93c4762a1bSJed Brown 94c4762a1bSJed Brown /* test MatForwardSolve() and MatBackwardSolve() with matrix reordering on aij matrix C */ 95c4762a1bSJed Brown if (CHOLESKY) { 969566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-triangular_solve", &TRIANGULAR)); 97c4762a1bSJed Brown if (TRIANGULAR) { 989566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Test MatForwardSolve...\n")); 999566063dSJacob Faibussowitsch PetscCall(MatForwardSolve(sA, b, ytmp)); 1009566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Test MatBackwardSolve...\n")); 1019566063dSJacob Faibussowitsch PetscCall(MatBackwardSolve(sA, ytmp, y)); 1029566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, x)); 1039566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm2)); 10448a46eb9SPierre Jolivet if (norm2 > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "MatForwardSolve and BackwardSolve: Norm of error=%g\n", (double)norm2)); 105c4762a1bSJed Brown } 106c4762a1bSJed Brown } 107c4762a1bSJed Brown 1089566063dSJacob Faibussowitsch PetscCall(MatSolve(sA, b, y)); 1099566063dSJacob Faibussowitsch PetscCall(MatDestroy(&sC)); 1109566063dSJacob Faibussowitsch PetscCall(MatDestroy(&sA)); 1119566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, x)); 1129566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm2)); 11348a46eb9SPierre Jolivet if (lf == -1 && norm2 > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, " reordered SEQAIJ: Cholesky/ICC levels %" PetscInt_FMT ", residual %g\n", lf, (double)norm2)); 114c4762a1bSJed Brown 115c4762a1bSJed Brown /* Free data structures */ 1169566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 1179566063dSJacob Faibussowitsch PetscCall(ISDestroy(&row)); 1189566063dSJacob Faibussowitsch PetscCall(ISDestroy(&col)); 1199566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm)); 1209566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1219566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 1229566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ytmp)); 1239566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 1249566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 125b122ec5aSJacob Faibussowitsch return 0; 126c4762a1bSJed Brown } 127c4762a1bSJed Brown 128c4762a1bSJed Brown /*TEST 129c4762a1bSJed Brown 130c4762a1bSJed Brown test: 131c4762a1bSJed Brown output_file: output/ex128.out 132c4762a1bSJed Brown 133c4762a1bSJed Brown test: 134c4762a1bSJed Brown suffix: 2 135c4762a1bSJed Brown args: -cholesky -triangular_solve 136c4762a1bSJed Brown 137c4762a1bSJed Brown TEST*/ 138