xref: /petsc/src/mat/impls/baij/seq/baijsolvtrann.c (revision d5b43468fb8780a8feea140ccd6fa3e6a50411cc)
1 #include <../src/mat/impls/baij/seq/baij.h>
2 #include <petsc/private/kernels/blockinvert.h>
3 
4 /* ----------------------------------------------------------- */
5 PetscErrorCode MatSolveTranspose_SeqBAIJ_N_inplace(Mat A, Vec bb, Vec xx)
6 {
7   Mat_SeqBAIJ       *a     = (Mat_SeqBAIJ *)A->data;
8   IS                 iscol = a->col, isrow = a->row;
9   const PetscInt    *r, *c, *rout, *cout, *ai = a->i, *aj = a->j, *vi;
10   PetscInt           i, nz, j;
11   const PetscInt     n = a->mbs, bs = A->rmap->bs, bs2 = a->bs2;
12   const MatScalar   *aa = a->a, *v;
13   PetscScalar       *x, *t, *ls;
14   const PetscScalar *b;
15 
16   PetscFunctionBegin;
17   PetscCall(VecGetArrayRead(bb, &b));
18   PetscCall(VecGetArray(xx, &x));
19   t = a->solve_work;
20 
21   PetscCall(ISGetIndices(isrow, &rout));
22   r = rout;
23   PetscCall(ISGetIndices(iscol, &cout));
24   c = cout;
25 
26   /* copy the b into temp work space according to permutation */
27   for (i = 0; i < n; i++) {
28     for (j = 0; j < bs; j++) t[i * bs + j] = b[c[i] * bs + j];
29   }
30 
31   /* forward solve the upper triangular transpose */
32   ls = a->solve_work + A->cmap->n;
33   for (i = 0; i < n; i++) {
34     PetscCall(PetscArraycpy(ls, t + i * bs, bs));
35     PetscKernel_w_gets_transA_times_v(bs, ls, aa + bs2 * a->diag[i], t + i * bs);
36     v  = aa + bs2 * (a->diag[i] + 1);
37     vi = aj + a->diag[i] + 1;
38     nz = ai[i + 1] - a->diag[i] - 1;
39     while (nz--) {
40       PetscKernel_v_gets_v_minus_transA_times_w(bs, t + bs * (*vi++), v, t + i * bs);
41       v += bs2;
42     }
43   }
44 
45   /* backward solve the lower triangular transpose */
46   for (i = n - 1; i >= 0; i--) {
47     v  = aa + bs2 * ai[i];
48     vi = aj + ai[i];
49     nz = a->diag[i] - ai[i];
50     while (nz--) {
51       PetscKernel_v_gets_v_minus_transA_times_w(bs, t + bs * (*vi++), v, t + i * bs);
52       v += bs2;
53     }
54   }
55 
56   /* copy t into x according to permutation */
57   for (i = 0; i < n; i++) {
58     for (j = 0; j < bs; j++) x[bs * r[i] + j] = t[bs * i + j];
59   }
60 
61   PetscCall(ISRestoreIndices(isrow, &rout));
62   PetscCall(ISRestoreIndices(iscol, &cout));
63   PetscCall(VecRestoreArrayRead(bb, &b));
64   PetscCall(VecRestoreArray(xx, &x));
65   PetscCall(PetscLogFlops(2.0 * (a->bs2) * (a->nz) - A->rmap->bs * A->cmap->n));
66   PetscFunctionReturn(0);
67 }
68 
69 PetscErrorCode MatSolveTranspose_SeqBAIJ_N(Mat A, Vec bb, Vec xx)
70 {
71   Mat_SeqBAIJ       *a     = (Mat_SeqBAIJ *)A->data;
72   IS                 iscol = a->col, isrow = a->row;
73   const PetscInt    *r, *c, *rout, *cout;
74   const PetscInt     n = a->mbs, *ai = a->i, *aj = a->j, *vi, *diag = a->diag;
75   PetscInt           i, j, nz;
76   const PetscInt     bs = A->rmap->bs, bs2 = a->bs2;
77   const MatScalar   *aa = a->a, *v;
78   PetscScalar       *x, *t, *ls;
79   const PetscScalar *b;
80 
81   PetscFunctionBegin;
82   PetscCall(VecGetArrayRead(bb, &b));
83   PetscCall(VecGetArray(xx, &x));
84   t = a->solve_work;
85 
86   PetscCall(ISGetIndices(isrow, &rout));
87   r = rout;
88   PetscCall(ISGetIndices(iscol, &cout));
89   c = cout;
90 
91   /* copy the b into temp work space according to permutation */
92   for (i = 0; i < n; i++) {
93     for (j = 0; j < bs; j++) t[i * bs + j] = b[c[i] * bs + j];
94   }
95 
96   /* forward solve the upper triangular transpose */
97   ls = a->solve_work + A->cmap->n;
98   for (i = 0; i < n; i++) {
99     PetscCall(PetscArraycpy(ls, t + i * bs, bs));
100     PetscKernel_w_gets_transA_times_v(bs, ls, aa + bs2 * diag[i], t + i * bs);
101     v  = aa + bs2 * (diag[i] - 1);
102     vi = aj + diag[i] - 1;
103     nz = diag[i] - diag[i + 1] - 1;
104     for (j = 0; j > -nz; j--) {
105       PetscKernel_v_gets_v_minus_transA_times_w(bs, t + bs * (vi[j]), v, t + i * bs);
106       v -= bs2;
107     }
108   }
109 
110   /* backward solve the lower triangular transpose */
111   for (i = n - 1; i >= 0; i--) {
112     v  = aa + bs2 * ai[i];
113     vi = aj + ai[i];
114     nz = ai[i + 1] - ai[i];
115     for (j = 0; j < nz; j++) {
116       PetscKernel_v_gets_v_minus_transA_times_w(bs, t + bs * (vi[j]), v, t + i * bs);
117       v += bs2;
118     }
119   }
120 
121   /* copy t into x according to permutation */
122   for (i = 0; i < n; i++) {
123     for (j = 0; j < bs; j++) x[bs * r[i] + j] = t[bs * i + j];
124   }
125 
126   PetscCall(ISRestoreIndices(isrow, &rout));
127   PetscCall(ISRestoreIndices(iscol, &cout));
128   PetscCall(VecRestoreArrayRead(bb, &b));
129   PetscCall(VecRestoreArray(xx, &x));
130   PetscCall(PetscLogFlops(2.0 * (a->bs2) * (a->nz) - A->rmap->bs * A->cmap->n));
131   PetscFunctionReturn(0);
132 }
133