Actual source code: ex17.c
1: static char help[] = "Tests the use of MatSolveTranspose().\n\n";
3: #include <petscmat.h>
5: int main(int argc, char **args)
6: {
7: Mat C, A;
8: PetscInt i, j, m = 5, n = 5, Ii, J;
9: PetscScalar v, five = 5.0, one = 1.0;
10: IS isrow, row, col;
11: Vec x, u, b;
12: PetscReal norm;
13: MatFactorInfo info;
15: PetscFunctionBeginUser;
16: PetscCall(PetscInitialize(&argc, &args, NULL, help));
17: PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL));
18: PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL));
20: PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, m * n, m * n, 5, NULL, &C));
21: PetscCall(MatSetUp(C));
23: /* create the matrix for the five point stencil, YET AGAIN*/
24: for (i = 0; i < m; i++) {
25: for (j = 0; j < n; j++) {
26: v = -1.0;
27: Ii = j + n * i;
28: if (i > 0) {
29: J = Ii - n;
30: PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
31: }
32: if (i < m - 1) {
33: J = Ii + n;
34: PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
35: }
36: if (j > 0) {
37: J = Ii - 1;
38: PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
39: }
40: if (j < n - 1) {
41: J = Ii + 1;
42: PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
43: }
44: v = 4.0;
45: PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, INSERT_VALUES));
46: }
47: }
48: PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
49: PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
51: PetscCall(ISCreateStride(PETSC_COMM_SELF, (m * n) / 2, 0, 2, &isrow));
52: PetscCall(MatZeroRowsIS(C, isrow, five, 0, 0));
54: PetscCall(VecCreateSeq(PETSC_COMM_SELF, m * n, &u));
55: PetscCall(VecDuplicate(u, &x));
56: PetscCall(VecDuplicate(u, &b));
57: PetscCall(VecSet(u, one));
59: PetscCall(MatMultTranspose(C, u, b));
61: /* Set default ordering to be Quotient Minimum Degree; also read
62: orderings from the options database */
63: PetscCall(MatGetOrdering(C, MATORDERINGQMD, &row, &col));
65: PetscCall(MatFactorInfoInitialize(&info));
66: PetscCall(MatGetFactor(C, MATSOLVERPETSC, MAT_FACTOR_LU, &A));
67: PetscCall(MatLUFactorSymbolic(A, C, row, col, &info));
68: PetscCall(MatLUFactorNumeric(A, C, &info));
69: PetscCall(MatSolveTranspose(A, b, x));
71: PetscCall(ISView(row, PETSC_VIEWER_STDOUT_SELF));
72: PetscCall(VecAXPY(x, -1.0, u));
73: PetscCall(VecNorm(x, NORM_2, &norm));
74: PetscCall(PetscPrintf(PETSC_COMM_SELF, "Norm of error %g\n", (double)norm));
76: PetscCall(ISDestroy(&row));
77: PetscCall(ISDestroy(&col));
78: PetscCall(ISDestroy(&isrow));
79: PetscCall(VecDestroy(&u));
80: PetscCall(VecDestroy(&x));
81: PetscCall(VecDestroy(&b));
82: PetscCall(MatDestroy(&C));
83: PetscCall(MatDestroy(&A));
84: PetscCall(PetscFinalize());
85: return 0;
86: }
88: /*TEST
90: test:
92: TEST*/