xref: /petsc/src/mat/impls/baij/seq/baijsolvtrannat5.c (revision 58d68138c660dfb4e9f5b03334792cd4f2ffd7cc)
1 #include <../src/mat/impls/baij/seq/baij.h>
2 
3 PetscErrorCode MatSolveTranspose_SeqBAIJ_5_NaturalOrdering_inplace(Mat A, Vec bb, Vec xx) {
4   Mat_SeqBAIJ     *a    = (Mat_SeqBAIJ *)A->data;
5   const PetscInt  *diag = a->diag, n = a->mbs, *vi, *ai = a->i, *aj = a->j;
6   PetscInt         i, nz, idx, idt, oidx;
7   const MatScalar *aa = a->a, *v;
8   PetscScalar      s1, s2, s3, s4, s5, x1, x2, x3, x4, x5, *x;
9 
10   PetscFunctionBegin;
11   PetscCall(VecCopy(bb, xx));
12   PetscCall(VecGetArray(xx, &x));
13 
14   /* forward solve the U^T */
15   idx = 0;
16   for (i = 0; i < n; i++) {
17     v  = aa + 25 * diag[i];
18     /* multiply by the inverse of the block diagonal */
19     x1 = x[idx];
20     x2 = x[1 + idx];
21     x3 = x[2 + idx];
22     x4 = x[3 + idx];
23     x5 = x[4 + idx];
24     s1 = v[0] * x1 + v[1] * x2 + v[2] * x3 + v[3] * x4 + v[4] * x5;
25     s2 = v[5] * x1 + v[6] * x2 + v[7] * x3 + v[8] * x4 + v[9] * x5;
26     s3 = v[10] * x1 + v[11] * x2 + v[12] * x3 + v[13] * x4 + v[14] * x5;
27     s4 = v[15] * x1 + v[16] * x2 + v[17] * x3 + v[18] * x4 + v[19] * x5;
28     s5 = v[20] * x1 + v[21] * x2 + v[22] * x3 + v[23] * x4 + v[24] * x5;
29     v += 25;
30 
31     vi = aj + diag[i] + 1;
32     nz = ai[i + 1] - diag[i] - 1;
33     while (nz--) {
34       oidx = 5 * (*vi++);
35       x[oidx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4 + v[4] * s5;
36       x[oidx + 1] -= v[5] * s1 + v[6] * s2 + v[7] * s3 + v[8] * s4 + v[9] * s5;
37       x[oidx + 2] -= v[10] * s1 + v[11] * s2 + v[12] * s3 + v[13] * s4 + v[14] * s5;
38       x[oidx + 3] -= v[15] * s1 + v[16] * s2 + v[17] * s3 + v[18] * s4 + v[19] * s5;
39       x[oidx + 4] -= v[20] * s1 + v[21] * s2 + v[22] * s3 + v[23] * s4 + v[24] * s5;
40       v += 25;
41     }
42     x[idx]     = s1;
43     x[1 + idx] = s2;
44     x[2 + idx] = s3;
45     x[3 + idx] = s4;
46     x[4 + idx] = s5;
47     idx += 5;
48   }
49   /* backward solve the L^T */
50   for (i = n - 1; i >= 0; i--) {
51     v   = aa + 25 * diag[i] - 25;
52     vi  = aj + diag[i] - 1;
53     nz  = diag[i] - ai[i];
54     idt = 5 * i;
55     s1  = x[idt];
56     s2  = x[1 + idt];
57     s3  = x[2 + idt];
58     s4  = x[3 + idt];
59     s5  = x[4 + idt];
60     while (nz--) {
61       idx = 5 * (*vi--);
62       x[idx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4 + v[4] * s5;
63       x[idx + 1] -= v[5] * s1 + v[6] * s2 + v[7] * s3 + v[8] * s4 + v[9] * s5;
64       x[idx + 2] -= v[10] * s1 + v[11] * s2 + v[12] * s3 + v[13] * s4 + v[14] * s5;
65       x[idx + 3] -= v[15] * s1 + v[16] * s2 + v[17] * s3 + v[18] * s4 + v[19] * s5;
66       x[idx + 4] -= v[20] * s1 + v[21] * s2 + v[22] * s3 + v[23] * s4 + v[24] * s5;
67       v -= 25;
68     }
69   }
70   PetscCall(VecRestoreArray(xx, &x));
71   PetscCall(PetscLogFlops(2.0 * 25 * (a->nz) - 5.0 * A->cmap->n));
72   PetscFunctionReturn(0);
73 }
74 
75 PetscErrorCode MatSolveTranspose_SeqBAIJ_5_NaturalOrdering(Mat A, Vec bb, Vec xx) {
76   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *)A->data;
77   const PetscInt   n = a->mbs, *vi, *ai = a->i, *aj = a->j, *diag = a->diag;
78   PetscInt         nz, idx, idt, j, i, oidx;
79   const PetscInt   bs = A->rmap->bs, bs2 = a->bs2;
80   const MatScalar *aa = a->a, *v;
81   PetscScalar      s1, s2, s3, s4, s5, x1, x2, x3, x4, x5, *x;
82 
83   PetscFunctionBegin;
84   PetscCall(VecCopy(bb, xx));
85   PetscCall(VecGetArray(xx, &x));
86 
87   /* forward solve the U^T */
88   idx = 0;
89   for (i = 0; i < n; i++) {
90     v  = aa + bs2 * diag[i];
91     /* multiply by the inverse of the block diagonal */
92     x1 = x[idx];
93     x2 = x[1 + idx];
94     x3 = x[2 + idx];
95     x4 = x[3 + idx];
96     x5 = x[4 + idx];
97     s1 = v[0] * x1 + v[1] * x2 + v[2] * x3 + v[3] * x4 + v[4] * x5;
98     s2 = v[5] * x1 + v[6] * x2 + v[7] * x3 + v[8] * x4 + v[9] * x5;
99     s3 = v[10] * x1 + v[11] * x2 + v[12] * x3 + v[13] * x4 + v[14] * x5;
100     s4 = v[15] * x1 + v[16] * x2 + v[17] * x3 + v[18] * x4 + v[19] * x5;
101     s5 = v[20] * x1 + v[21] * x2 + v[22] * x3 + v[23] * x4 + v[24] * x5;
102     v -= bs2;
103 
104     vi = aj + diag[i] - 1;
105     nz = diag[i] - diag[i + 1] - 1;
106     for (j = 0; j > -nz; j--) {
107       oidx = bs * vi[j];
108       x[oidx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4 + v[4] * s5;
109       x[oidx + 1] -= v[5] * s1 + v[6] * s2 + v[7] * s3 + v[8] * s4 + v[9] * s5;
110       x[oidx + 2] -= v[10] * s1 + v[11] * s2 + v[12] * s3 + v[13] * s4 + v[14] * s5;
111       x[oidx + 3] -= v[15] * s1 + v[16] * s2 + v[17] * s3 + v[18] * s4 + v[19] * s5;
112       x[oidx + 4] -= v[20] * s1 + v[21] * s2 + v[22] * s3 + v[23] * s4 + v[24] * s5;
113       v -= bs2;
114     }
115     x[idx]     = s1;
116     x[1 + idx] = s2;
117     x[2 + idx] = s3;
118     x[3 + idx] = s4;
119     x[4 + idx] = s5;
120     idx += bs;
121   }
122   /* backward solve the L^T */
123   for (i = n - 1; i >= 0; i--) {
124     v   = aa + bs2 * ai[i];
125     vi  = aj + ai[i];
126     nz  = ai[i + 1] - ai[i];
127     idt = bs * i;
128     s1  = x[idt];
129     s2  = x[1 + idt];
130     s3  = x[2 + idt];
131     s4  = x[3 + idt];
132     s5  = x[4 + idt];
133     for (j = 0; j < nz; j++) {
134       idx = bs * vi[j];
135       x[idx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4 + v[4] * s5;
136       x[idx + 1] -= v[5] * s1 + v[6] * s2 + v[7] * s3 + v[8] * s4 + v[9] * s5;
137       x[idx + 2] -= v[10] * s1 + v[11] * s2 + v[12] * s3 + v[13] * s4 + v[14] * s5;
138       x[idx + 3] -= v[15] * s1 + v[16] * s2 + v[17] * s3 + v[18] * s4 + v[19] * s5;
139       x[idx + 4] -= v[20] * s1 + v[21] * s2 + v[22] * s3 + v[23] * s4 + v[24] * s5;
140       v += bs2;
141     }
142   }
143   PetscCall(VecRestoreArray(xx, &x));
144   PetscCall(PetscLogFlops(2.0 * bs2 * (a->nz) - bs * A->cmap->n));
145   PetscFunctionReturn(0);
146 }
147