1c4762a1bSJed Brown static char help[] = "Tests matrix factorization. Note that most users should\n\ 2c4762a1bSJed Brown employ the KSP interface to the linear solvers instead of using the factorization\n\ 3c4762a1bSJed Brown routines directly.\n\n"; 4c4762a1bSJed Brown 5c4762a1bSJed Brown #include <petscmat.h> 6c4762a1bSJed Brown 7d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 8d71ae5a4SJacob Faibussowitsch { 9c4762a1bSJed Brown Mat C, LU; 10c4762a1bSJed Brown MatInfo info; 111a6d72e3SSatish Balay PetscInt i, j, m, n, Ii, J; 12c4762a1bSJed Brown PetscScalar v, one = 1.0; 13c4762a1bSJed Brown IS perm, iperm; 14c4762a1bSJed Brown Vec x, u, b; 151a6d72e3SSatish Balay PetscReal norm, fill; 16c4762a1bSJed Brown MatFactorInfo luinfo; 17c4762a1bSJed Brown 18327415f7SBarry Smith PetscFunctionBeginUser; 19*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help)); 201a6d72e3SSatish Balay 21d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Mat test ex7 options", "Mat"); 229371c9d4SSatish Balay m = 3; 239371c9d4SSatish Balay n = 3; 249371c9d4SSatish Balay fill = 2.0; 259566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-m", "Number of rows in grid", NULL, m, &m, NULL)); 269566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-n", "Number of columns in grid", NULL, n, &n, NULL)); 279566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-fill", "Expected fill ratio for factorization", NULL, fill, &fill, NULL)); 281a6d72e3SSatish Balay 29d0609cedSBarry Smith PetscOptionsEnd(); 301a6d72e3SSatish Balay 31c4762a1bSJed Brown /* Create the matrix for the five point stencil, YET AGAIN */ 329566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &C)); 339566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, m * n, m * n)); 349566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C)); 359566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 36c4762a1bSJed Brown for (i = 0; i < m; i++) { 37c4762a1bSJed Brown for (j = 0; j < n; j++) { 389371c9d4SSatish Balay v = -1.0; 399371c9d4SSatish Balay Ii = j + n * i; 409371c9d4SSatish Balay if (i > 0) { 419371c9d4SSatish Balay J = Ii - n; 429371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 439371c9d4SSatish Balay } 449371c9d4SSatish Balay if (i < m - 1) { 459371c9d4SSatish Balay J = Ii + n; 469371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 479371c9d4SSatish Balay } 489371c9d4SSatish Balay if (j > 0) { 499371c9d4SSatish Balay J = Ii - 1; 509371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 519371c9d4SSatish Balay } 529371c9d4SSatish Balay if (j < n - 1) { 539371c9d4SSatish Balay J = Ii + 1; 549371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 559371c9d4SSatish Balay } 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)); 629566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(C, MATORDERINGRCM, &perm, &iperm)); 639566063dSJacob Faibussowitsch PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD)); 649566063dSJacob Faibussowitsch PetscCall(ISView(perm, PETSC_VIEWER_STDOUT_SELF)); 65c4762a1bSJed Brown 669566063dSJacob Faibussowitsch PetscCall(MatFactorInfoInitialize(&luinfo)); 67c4762a1bSJed Brown 681a6d72e3SSatish Balay luinfo.fill = fill; 69c4762a1bSJed Brown luinfo.dtcol = 0.0; 70c4762a1bSJed Brown luinfo.zeropivot = 1.e-14; 71c4762a1bSJed Brown luinfo.pivotinblocks = 1.0; 72c4762a1bSJed Brown 739566063dSJacob Faibussowitsch PetscCall(MatGetFactor(C, MATSOLVERPETSC, MAT_FACTOR_LU, &LU)); 749566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(LU, C, perm, iperm, &luinfo)); 759566063dSJacob Faibussowitsch PetscCall(MatLUFactorNumeric(LU, C, &luinfo)); 76c4762a1bSJed Brown 779566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, m * n, &u)); 789566063dSJacob Faibussowitsch PetscCall(VecSet(u, one)); 799566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &x)); 809566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &b)); 81c4762a1bSJed Brown 829566063dSJacob Faibussowitsch PetscCall(MatMult(C, u, b)); 839566063dSJacob Faibussowitsch PetscCall(MatSolve(LU, b, x)); 84c4762a1bSJed Brown 859566063dSJacob Faibussowitsch PetscCall(VecView(b, PETSC_VIEWER_STDOUT_SELF)); 869566063dSJacob Faibussowitsch PetscCall(VecView(x, PETSC_VIEWER_STDOUT_SELF)); 87c4762a1bSJed Brown 889566063dSJacob Faibussowitsch PetscCall(VecAXPY(x, -1.0, u)); 899566063dSJacob Faibussowitsch PetscCall(VecNorm(x, NORM_2, &norm)); 909566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Norm of error %g\n", (double)norm)); 91c4762a1bSJed Brown 929566063dSJacob Faibussowitsch PetscCall(MatGetInfo(C, MAT_LOCAL, &info)); 939566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "original matrix nonzeros = %" PetscInt_FMT "\n", (PetscInt)info.nz_used)); 949566063dSJacob Faibussowitsch PetscCall(MatGetInfo(LU, MAT_LOCAL, &info)); 959566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "factored matrix nonzeros = %" PetscInt_FMT "\n", (PetscInt)info.nz_used)); 96c4762a1bSJed Brown 979566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 989566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 999566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1009566063dSJacob Faibussowitsch PetscCall(ISDestroy(&perm)); 1019566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iperm)); 1029566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 1039566063dSJacob Faibussowitsch PetscCall(MatDestroy(&LU)); 104c4762a1bSJed Brown 1059566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 106b122ec5aSJacob Faibussowitsch return 0; 107c4762a1bSJed Brown } 108c4762a1bSJed Brown 109c4762a1bSJed Brown /*TEST 110c4762a1bSJed Brown 111c4762a1bSJed Brown test: 1121a6d72e3SSatish Balay suffix: 1 1138cc725e6SPierre Jolivet filter: grep -v " MPI process" 1141a6d72e3SSatish Balay 1151a6d72e3SSatish Balay test: 1161a6d72e3SSatish Balay suffix: 2 1171a6d72e3SSatish Balay args: -m 1 -n 1 -fill 0.49 1188cc725e6SPierre Jolivet filter: grep -v " MPI process" 119c4762a1bSJed Brown 120c4762a1bSJed Brown TEST*/ 121