1c4762a1bSJed Brown static char help[] = "Tests the use of MatSolveTranspose().\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscmat.h> 4c4762a1bSJed Brown 5d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 6d71ae5a4SJacob Faibussowitsch { 7c4762a1bSJed Brown Mat C, A; 8c4762a1bSJed Brown PetscInt i, j, m = 5, n = 5, Ii, J; 9c4762a1bSJed Brown PetscScalar v, five = 5.0, one = 1.0; 10c4762a1bSJed Brown IS isrow, row, col; 11c4762a1bSJed Brown Vec x, u, b; 12c4762a1bSJed Brown PetscReal norm; 13c4762a1bSJed Brown MatFactorInfo info; 14c4762a1bSJed Brown 15327415f7SBarry Smith PetscFunctionBeginUser; 16*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help)); 179566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL)); 189566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL)); 19c4762a1bSJed Brown 209566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, m * n, m * n, 5, NULL, &C)); 219566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 22c4762a1bSJed Brown 23c4762a1bSJed Brown /* create the matrix for the five point stencil, YET AGAIN*/ 24c4762a1bSJed Brown for (i = 0; i < m; i++) { 25c4762a1bSJed Brown for (j = 0; j < n; j++) { 269371c9d4SSatish Balay v = -1.0; 279371c9d4SSatish Balay Ii = j + n * i; 289371c9d4SSatish Balay if (i > 0) { 299371c9d4SSatish Balay J = Ii - n; 309371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 319371c9d4SSatish Balay } 329371c9d4SSatish Balay if (i < m - 1) { 339371c9d4SSatish Balay J = Ii + n; 349371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 359371c9d4SSatish Balay } 369371c9d4SSatish Balay if (j > 0) { 379371c9d4SSatish Balay J = Ii - 1; 389371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 399371c9d4SSatish Balay } 409371c9d4SSatish Balay if (j < n - 1) { 419371c9d4SSatish Balay J = Ii + 1; 429371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 439371c9d4SSatish Balay } 449371c9d4SSatish Balay v = 4.0; 459371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, INSERT_VALUES)); 46c4762a1bSJed Brown } 47c4762a1bSJed Brown } 489566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 499566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 50c4762a1bSJed Brown 519566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, (m * n) / 2, 0, 2, &isrow)); 529566063dSJacob Faibussowitsch PetscCall(MatZeroRowsIS(C, isrow, five, 0, 0)); 53c4762a1bSJed Brown 549566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, m * n, &u)); 559566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &x)); 569566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &b)); 579566063dSJacob Faibussowitsch PetscCall(VecSet(u, one)); 58c4762a1bSJed Brown 599566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(C, u, b)); 60c4762a1bSJed Brown 61c4762a1bSJed Brown /* Set default ordering to be Quotient Minimum Degree; also read 62c4762a1bSJed Brown orderings from the options database */ 639566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(C, MATORDERINGQMD, &row, &col)); 64c4762a1bSJed Brown 659566063dSJacob Faibussowitsch PetscCall(MatFactorInfoInitialize(&info)); 669566063dSJacob Faibussowitsch PetscCall(MatGetFactor(C, MATSOLVERPETSC, MAT_FACTOR_LU, &A)); 679566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(A, C, row, col, &info)); 689566063dSJacob Faibussowitsch PetscCall(MatLUFactorNumeric(A, C, &info)); 699566063dSJacob Faibussowitsch PetscCall(MatSolveTranspose(A, b, x)); 70c4762a1bSJed Brown 719566063dSJacob Faibussowitsch PetscCall(ISView(row, PETSC_VIEWER_STDOUT_SELF)); 729566063dSJacob Faibussowitsch PetscCall(VecAXPY(x, -1.0, u)); 739566063dSJacob Faibussowitsch PetscCall(VecNorm(x, NORM_2, &norm)); 749566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Norm of error %g\n", (double)norm)); 75c4762a1bSJed Brown 769566063dSJacob Faibussowitsch PetscCall(ISDestroy(&row)); 779566063dSJacob Faibussowitsch PetscCall(ISDestroy(&col)); 789566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrow)); 799566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 819566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 829566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 839566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 849566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 85b122ec5aSJacob Faibussowitsch return 0; 86c4762a1bSJed Brown } 87c4762a1bSJed Brown 88c4762a1bSJed Brown /*TEST 89c4762a1bSJed Brown 90c4762a1bSJed Brown test: 91c4762a1bSJed Brown 92c4762a1bSJed Brown TEST*/ 93