xref: /petsc/src/mat/impls/baij/seq/baijfact7.c (revision d9acb416d05abeed0a33bde3a81aeb2ea0364f6a)
1 /*
2     Factorization code for BAIJ format.
3 */
4 #include <../src/mat/impls/baij/seq/baij.h>
5 #include <petsc/private/kernels/blockinvert.h>
6 
7 /*
8       Version for when blocks are 6 by 6
9 */
10 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_6_inplace(Mat C, Mat A, const MatFactorInfo *info)
11 {
12   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data;
13   IS              isrow = b->row, isicol = b->icol;
14   const PetscInt *ajtmpold, *ajtmp, *diag_offset = b->diag, *r, *ic, *bi = b->i, *bj = b->j, *ai = a->i, *aj = a->j, *pj;
15   PetscInt        nz, row, i, j, n = a->mbs, idx;
16   MatScalar      *pv, *v, *rtmp, *pc, *w, *x;
17   MatScalar       p1, p2, p3, p4, m1, m2, m3, m4, m5, m6, m7, m8, m9, x1, x2, x3, x4;
18   MatScalar       p5, p6, p7, p8, p9, x5, x6, x7, x8, x9, x10, x11, x12, x13, x14, x15, x16;
19   MatScalar       x17, x18, x19, x20, x21, x22, x23, x24, x25, p10, p11, p12, p13, p14;
20   MatScalar       p15, p16, p17, p18, p19, p20, p21, p22, p23, p24, p25, m10, m11, m12;
21   MatScalar       m13, m14, m15, m16, m17, m18, m19, m20, m21, m22, m23, m24, m25;
22   MatScalar       p26, p27, p28, p29, p30, p31, p32, p33, p34, p35, p36;
23   MatScalar       x26, x27, x28, x29, x30, x31, x32, x33, x34, x35, x36;
24   MatScalar       m26, m27, m28, m29, m30, m31, m32, m33, m34, m35, m36;
25   MatScalar      *ba = b->a, *aa = a->a;
26   PetscReal       shift = info->shiftamount;
27   PetscBool       allowzeropivot, zeropivotdetected;
28 
29   PetscFunctionBegin;
30   allowzeropivot = PetscNot(A->erroriffailure);
31   PetscCall(ISGetIndices(isrow, &r));
32   PetscCall(ISGetIndices(isicol, &ic));
33   PetscCall(PetscMalloc1(36 * (n + 1), &rtmp));
34 
35   for (i = 0; i < n; i++) {
36     nz    = bi[i + 1] - bi[i];
37     ajtmp = bj + bi[i];
38     for (j = 0; j < nz; j++) {
39       x    = rtmp + 36 * ajtmp[j];
40       x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = x[9] = 0.0;
41       x[10] = x[11] = x[12] = x[13] = x[14] = x[15] = x[16] = x[17] = 0.0;
42       x[18] = x[19] = x[20] = x[21] = x[22] = x[23] = x[24] = x[25] = 0.0;
43       x[26] = x[27] = x[28] = x[29] = x[30] = x[31] = x[32] = x[33] = 0.0;
44       x[34] = x[35] = 0.0;
45     }
46     /* load in initial (unfactored row) */
47     idx      = r[i];
48     nz       = ai[idx + 1] - ai[idx];
49     ajtmpold = aj + ai[idx];
50     v        = aa + 36 * ai[idx];
51     for (j = 0; j < nz; j++) {
52       x     = rtmp + 36 * ic[ajtmpold[j]];
53       x[0]  = v[0];
54       x[1]  = v[1];
55       x[2]  = v[2];
56       x[3]  = v[3];
57       x[4]  = v[4];
58       x[5]  = v[5];
59       x[6]  = v[6];
60       x[7]  = v[7];
61       x[8]  = v[8];
62       x[9]  = v[9];
63       x[10] = v[10];
64       x[11] = v[11];
65       x[12] = v[12];
66       x[13] = v[13];
67       x[14] = v[14];
68       x[15] = v[15];
69       x[16] = v[16];
70       x[17] = v[17];
71       x[18] = v[18];
72       x[19] = v[19];
73       x[20] = v[20];
74       x[21] = v[21];
75       x[22] = v[22];
76       x[23] = v[23];
77       x[24] = v[24];
78       x[25] = v[25];
79       x[26] = v[26];
80       x[27] = v[27];
81       x[28] = v[28];
82       x[29] = v[29];
83       x[30] = v[30];
84       x[31] = v[31];
85       x[32] = v[32];
86       x[33] = v[33];
87       x[34] = v[34];
88       x[35] = v[35];
89       v += 36;
90     }
91     row = *ajtmp++;
92     while (row < i) {
93       pc  = rtmp + 36 * row;
94       p1  = pc[0];
95       p2  = pc[1];
96       p3  = pc[2];
97       p4  = pc[3];
98       p5  = pc[4];
99       p6  = pc[5];
100       p7  = pc[6];
101       p8  = pc[7];
102       p9  = pc[8];
103       p10 = pc[9];
104       p11 = pc[10];
105       p12 = pc[11];
106       p13 = pc[12];
107       p14 = pc[13];
108       p15 = pc[14];
109       p16 = pc[15];
110       p17 = pc[16];
111       p18 = pc[17];
112       p19 = pc[18];
113       p20 = pc[19];
114       p21 = pc[20];
115       p22 = pc[21];
116       p23 = pc[22];
117       p24 = pc[23];
118       p25 = pc[24];
119       p26 = pc[25];
120       p27 = pc[26];
121       p28 = pc[27];
122       p29 = pc[28];
123       p30 = pc[29];
124       p31 = pc[30];
125       p32 = pc[31];
126       p33 = pc[32];
127       p34 = pc[33];
128       p35 = pc[34];
129       p36 = pc[35];
130       if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 || p5 != 0.0 || p6 != 0.0 || p7 != 0.0 || p8 != 0.0 || p9 != 0.0 || p10 != 0.0 || p11 != 0.0 || p12 != 0.0 || p13 != 0.0 || p14 != 0.0 || p15 != 0.0 || p16 != 0.0 || p17 != 0.0 || p18 != 0.0 || p19 != 0.0 || p20 != 0.0 || p21 != 0.0 || p22 != 0.0 || p23 != 0.0 || p24 != 0.0 || p25 != 0.0 || p26 != 0.0 || p27 != 0.0 || p28 != 0.0 || p29 != 0.0 || p30 != 0.0 || p31 != 0.0 || p32 != 0.0 || p33 != 0.0 || p34 != 0.0 || p35 != 0.0 || p36 != 0.0) {
131         pv    = ba + 36 * diag_offset[row];
132         pj    = bj + diag_offset[row] + 1;
133         x1    = pv[0];
134         x2    = pv[1];
135         x3    = pv[2];
136         x4    = pv[3];
137         x5    = pv[4];
138         x6    = pv[5];
139         x7    = pv[6];
140         x8    = pv[7];
141         x9    = pv[8];
142         x10   = pv[9];
143         x11   = pv[10];
144         x12   = pv[11];
145         x13   = pv[12];
146         x14   = pv[13];
147         x15   = pv[14];
148         x16   = pv[15];
149         x17   = pv[16];
150         x18   = pv[17];
151         x19   = pv[18];
152         x20   = pv[19];
153         x21   = pv[20];
154         x22   = pv[21];
155         x23   = pv[22];
156         x24   = pv[23];
157         x25   = pv[24];
158         x26   = pv[25];
159         x27   = pv[26];
160         x28   = pv[27];
161         x29   = pv[28];
162         x30   = pv[29];
163         x31   = pv[30];
164         x32   = pv[31];
165         x33   = pv[32];
166         x34   = pv[33];
167         x35   = pv[34];
168         x36   = pv[35];
169         pc[0] = m1 = p1 * x1 + p7 * x2 + p13 * x3 + p19 * x4 + p25 * x5 + p31 * x6;
170         pc[1] = m2 = p2 * x1 + p8 * x2 + p14 * x3 + p20 * x4 + p26 * x5 + p32 * x6;
171         pc[2] = m3 = p3 * x1 + p9 * x2 + p15 * x3 + p21 * x4 + p27 * x5 + p33 * x6;
172         pc[3] = m4 = p4 * x1 + p10 * x2 + p16 * x3 + p22 * x4 + p28 * x5 + p34 * x6;
173         pc[4] = m5 = p5 * x1 + p11 * x2 + p17 * x3 + p23 * x4 + p29 * x5 + p35 * x6;
174         pc[5] = m6 = p6 * x1 + p12 * x2 + p18 * x3 + p24 * x4 + p30 * x5 + p36 * x6;
175 
176         pc[6] = m7 = p1 * x7 + p7 * x8 + p13 * x9 + p19 * x10 + p25 * x11 + p31 * x12;
177         pc[7] = m8 = p2 * x7 + p8 * x8 + p14 * x9 + p20 * x10 + p26 * x11 + p32 * x12;
178         pc[8] = m9 = p3 * x7 + p9 * x8 + p15 * x9 + p21 * x10 + p27 * x11 + p33 * x12;
179         pc[9] = m10 = p4 * x7 + p10 * x8 + p16 * x9 + p22 * x10 + p28 * x11 + p34 * x12;
180         pc[10] = m11 = p5 * x7 + p11 * x8 + p17 * x9 + p23 * x10 + p29 * x11 + p35 * x12;
181         pc[11] = m12 = p6 * x7 + p12 * x8 + p18 * x9 + p24 * x10 + p30 * x11 + p36 * x12;
182 
183         pc[12] = m13 = p1 * x13 + p7 * x14 + p13 * x15 + p19 * x16 + p25 * x17 + p31 * x18;
184         pc[13] = m14 = p2 * x13 + p8 * x14 + p14 * x15 + p20 * x16 + p26 * x17 + p32 * x18;
185         pc[14] = m15 = p3 * x13 + p9 * x14 + p15 * x15 + p21 * x16 + p27 * x17 + p33 * x18;
186         pc[15] = m16 = p4 * x13 + p10 * x14 + p16 * x15 + p22 * x16 + p28 * x17 + p34 * x18;
187         pc[16] = m17 = p5 * x13 + p11 * x14 + p17 * x15 + p23 * x16 + p29 * x17 + p35 * x18;
188         pc[17] = m18 = p6 * x13 + p12 * x14 + p18 * x15 + p24 * x16 + p30 * x17 + p36 * x18;
189 
190         pc[18] = m19 = p1 * x19 + p7 * x20 + p13 * x21 + p19 * x22 + p25 * x23 + p31 * x24;
191         pc[19] = m20 = p2 * x19 + p8 * x20 + p14 * x21 + p20 * x22 + p26 * x23 + p32 * x24;
192         pc[20] = m21 = p3 * x19 + p9 * x20 + p15 * x21 + p21 * x22 + p27 * x23 + p33 * x24;
193         pc[21] = m22 = p4 * x19 + p10 * x20 + p16 * x21 + p22 * x22 + p28 * x23 + p34 * x24;
194         pc[22] = m23 = p5 * x19 + p11 * x20 + p17 * x21 + p23 * x22 + p29 * x23 + p35 * x24;
195         pc[23] = m24 = p6 * x19 + p12 * x20 + p18 * x21 + p24 * x22 + p30 * x23 + p36 * x24;
196 
197         pc[24] = m25 = p1 * x25 + p7 * x26 + p13 * x27 + p19 * x28 + p25 * x29 + p31 * x30;
198         pc[25] = m26 = p2 * x25 + p8 * x26 + p14 * x27 + p20 * x28 + p26 * x29 + p32 * x30;
199         pc[26] = m27 = p3 * x25 + p9 * x26 + p15 * x27 + p21 * x28 + p27 * x29 + p33 * x30;
200         pc[27] = m28 = p4 * x25 + p10 * x26 + p16 * x27 + p22 * x28 + p28 * x29 + p34 * x30;
201         pc[28] = m29 = p5 * x25 + p11 * x26 + p17 * x27 + p23 * x28 + p29 * x29 + p35 * x30;
202         pc[29] = m30 = p6 * x25 + p12 * x26 + p18 * x27 + p24 * x28 + p30 * x29 + p36 * x30;
203 
204         pc[30] = m31 = p1 * x31 + p7 * x32 + p13 * x33 + p19 * x34 + p25 * x35 + p31 * x36;
205         pc[31] = m32 = p2 * x31 + p8 * x32 + p14 * x33 + p20 * x34 + p26 * x35 + p32 * x36;
206         pc[32] = m33 = p3 * x31 + p9 * x32 + p15 * x33 + p21 * x34 + p27 * x35 + p33 * x36;
207         pc[33] = m34 = p4 * x31 + p10 * x32 + p16 * x33 + p22 * x34 + p28 * x35 + p34 * x36;
208         pc[34] = m35 = p5 * x31 + p11 * x32 + p17 * x33 + p23 * x34 + p29 * x35 + p35 * x36;
209         pc[35] = m36 = p6 * x31 + p12 * x32 + p18 * x33 + p24 * x34 + p30 * x35 + p36 * x36;
210 
211         nz = bi[row + 1] - diag_offset[row] - 1;
212         pv += 36;
213         for (j = 0; j < nz; j++) {
214           x1  = pv[0];
215           x2  = pv[1];
216           x3  = pv[2];
217           x4  = pv[3];
218           x5  = pv[4];
219           x6  = pv[5];
220           x7  = pv[6];
221           x8  = pv[7];
222           x9  = pv[8];
223           x10 = pv[9];
224           x11 = pv[10];
225           x12 = pv[11];
226           x13 = pv[12];
227           x14 = pv[13];
228           x15 = pv[14];
229           x16 = pv[15];
230           x17 = pv[16];
231           x18 = pv[17];
232           x19 = pv[18];
233           x20 = pv[19];
234           x21 = pv[20];
235           x22 = pv[21];
236           x23 = pv[22];
237           x24 = pv[23];
238           x25 = pv[24];
239           x26 = pv[25];
240           x27 = pv[26];
241           x28 = pv[27];
242           x29 = pv[28];
243           x30 = pv[29];
244           x31 = pv[30];
245           x32 = pv[31];
246           x33 = pv[32];
247           x34 = pv[33];
248           x35 = pv[34];
249           x36 = pv[35];
250           x   = rtmp + 36 * pj[j];
251           x[0] -= m1 * x1 + m7 * x2 + m13 * x3 + m19 * x4 + m25 * x5 + m31 * x6;
252           x[1] -= m2 * x1 + m8 * x2 + m14 * x3 + m20 * x4 + m26 * x5 + m32 * x6;
253           x[2] -= m3 * x1 + m9 * x2 + m15 * x3 + m21 * x4 + m27 * x5 + m33 * x6;
254           x[3] -= m4 * x1 + m10 * x2 + m16 * x3 + m22 * x4 + m28 * x5 + m34 * x6;
255           x[4] -= m5 * x1 + m11 * x2 + m17 * x3 + m23 * x4 + m29 * x5 + m35 * x6;
256           x[5] -= m6 * x1 + m12 * x2 + m18 * x3 + m24 * x4 + m30 * x5 + m36 * x6;
257 
258           x[6] -= m1 * x7 + m7 * x8 + m13 * x9 + m19 * x10 + m25 * x11 + m31 * x12;
259           x[7] -= m2 * x7 + m8 * x8 + m14 * x9 + m20 * x10 + m26 * x11 + m32 * x12;
260           x[8] -= m3 * x7 + m9 * x8 + m15 * x9 + m21 * x10 + m27 * x11 + m33 * x12;
261           x[9] -= m4 * x7 + m10 * x8 + m16 * x9 + m22 * x10 + m28 * x11 + m34 * x12;
262           x[10] -= m5 * x7 + m11 * x8 + m17 * x9 + m23 * x10 + m29 * x11 + m35 * x12;
263           x[11] -= m6 * x7 + m12 * x8 + m18 * x9 + m24 * x10 + m30 * x11 + m36 * x12;
264 
265           x[12] -= m1 * x13 + m7 * x14 + m13 * x15 + m19 * x16 + m25 * x17 + m31 * x18;
266           x[13] -= m2 * x13 + m8 * x14 + m14 * x15 + m20 * x16 + m26 * x17 + m32 * x18;
267           x[14] -= m3 * x13 + m9 * x14 + m15 * x15 + m21 * x16 + m27 * x17 + m33 * x18;
268           x[15] -= m4 * x13 + m10 * x14 + m16 * x15 + m22 * x16 + m28 * x17 + m34 * x18;
269           x[16] -= m5 * x13 + m11 * x14 + m17 * x15 + m23 * x16 + m29 * x17 + m35 * x18;
270           x[17] -= m6 * x13 + m12 * x14 + m18 * x15 + m24 * x16 + m30 * x17 + m36 * x18;
271 
272           x[18] -= m1 * x19 + m7 * x20 + m13 * x21 + m19 * x22 + m25 * x23 + m31 * x24;
273           x[19] -= m2 * x19 + m8 * x20 + m14 * x21 + m20 * x22 + m26 * x23 + m32 * x24;
274           x[20] -= m3 * x19 + m9 * x20 + m15 * x21 + m21 * x22 + m27 * x23 + m33 * x24;
275           x[21] -= m4 * x19 + m10 * x20 + m16 * x21 + m22 * x22 + m28 * x23 + m34 * x24;
276           x[22] -= m5 * x19 + m11 * x20 + m17 * x21 + m23 * x22 + m29 * x23 + m35 * x24;
277           x[23] -= m6 * x19 + m12 * x20 + m18 * x21 + m24 * x22 + m30 * x23 + m36 * x24;
278 
279           x[24] -= m1 * x25 + m7 * x26 + m13 * x27 + m19 * x28 + m25 * x29 + m31 * x30;
280           x[25] -= m2 * x25 + m8 * x26 + m14 * x27 + m20 * x28 + m26 * x29 + m32 * x30;
281           x[26] -= m3 * x25 + m9 * x26 + m15 * x27 + m21 * x28 + m27 * x29 + m33 * x30;
282           x[27] -= m4 * x25 + m10 * x26 + m16 * x27 + m22 * x28 + m28 * x29 + m34 * x30;
283           x[28] -= m5 * x25 + m11 * x26 + m17 * x27 + m23 * x28 + m29 * x29 + m35 * x30;
284           x[29] -= m6 * x25 + m12 * x26 + m18 * x27 + m24 * x28 + m30 * x29 + m36 * x30;
285 
286           x[30] -= m1 * x31 + m7 * x32 + m13 * x33 + m19 * x34 + m25 * x35 + m31 * x36;
287           x[31] -= m2 * x31 + m8 * x32 + m14 * x33 + m20 * x34 + m26 * x35 + m32 * x36;
288           x[32] -= m3 * x31 + m9 * x32 + m15 * x33 + m21 * x34 + m27 * x35 + m33 * x36;
289           x[33] -= m4 * x31 + m10 * x32 + m16 * x33 + m22 * x34 + m28 * x35 + m34 * x36;
290           x[34] -= m5 * x31 + m11 * x32 + m17 * x33 + m23 * x34 + m29 * x35 + m35 * x36;
291           x[35] -= m6 * x31 + m12 * x32 + m18 * x33 + m24 * x34 + m30 * x35 + m36 * x36;
292 
293           pv += 36;
294         }
295         PetscCall(PetscLogFlops(432.0 * nz + 396.0));
296       }
297       row = *ajtmp++;
298     }
299     /* finished row so stick it into b->a */
300     pv = ba + 36 * bi[i];
301     pj = bj + bi[i];
302     nz = bi[i + 1] - bi[i];
303     for (j = 0; j < nz; j++) {
304       x      = rtmp + 36 * pj[j];
305       pv[0]  = x[0];
306       pv[1]  = x[1];
307       pv[2]  = x[2];
308       pv[3]  = x[3];
309       pv[4]  = x[4];
310       pv[5]  = x[5];
311       pv[6]  = x[6];
312       pv[7]  = x[7];
313       pv[8]  = x[8];
314       pv[9]  = x[9];
315       pv[10] = x[10];
316       pv[11] = x[11];
317       pv[12] = x[12];
318       pv[13] = x[13];
319       pv[14] = x[14];
320       pv[15] = x[15];
321       pv[16] = x[16];
322       pv[17] = x[17];
323       pv[18] = x[18];
324       pv[19] = x[19];
325       pv[20] = x[20];
326       pv[21] = x[21];
327       pv[22] = x[22];
328       pv[23] = x[23];
329       pv[24] = x[24];
330       pv[25] = x[25];
331       pv[26] = x[26];
332       pv[27] = x[27];
333       pv[28] = x[28];
334       pv[29] = x[29];
335       pv[30] = x[30];
336       pv[31] = x[31];
337       pv[32] = x[32];
338       pv[33] = x[33];
339       pv[34] = x[34];
340       pv[35] = x[35];
341       pv += 36;
342     }
343     /* invert diagonal block */
344     w = ba + 36 * diag_offset[i];
345     PetscCall(PetscKernel_A_gets_inverse_A_6(w, shift, allowzeropivot, &zeropivotdetected));
346     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
347   }
348 
349   PetscCall(PetscFree(rtmp));
350   PetscCall(ISRestoreIndices(isicol, &ic));
351   PetscCall(ISRestoreIndices(isrow, &r));
352 
353   C->ops->solve          = MatSolve_SeqBAIJ_6_inplace;
354   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_6_inplace;
355   C->assembled           = PETSC_TRUE;
356 
357   PetscCall(PetscLogFlops(1.333333333333 * 6 * 6 * 6 * b->mbs)); /* from inverting diagonal blocks */
358   PetscFunctionReturn(PETSC_SUCCESS);
359 }
360 
361 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_6(Mat B, Mat A, const MatFactorInfo *info)
362 {
363   Mat             C = B;
364   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data;
365   IS              isrow = b->row, isicol = b->icol;
366   const PetscInt *r, *ic;
367   PetscInt        i, j, k, nz, nzL, row;
368   const PetscInt  n = a->mbs, *ai = a->i, *aj = a->j, *bi = b->i, *bj = b->j;
369   const PetscInt *ajtmp, *bjtmp, *bdiag = b->diag, *pj, bs2 = a->bs2;
370   MatScalar      *rtmp, *pc, *mwork, *v, *pv, *aa = a->a;
371   PetscInt        flg;
372   PetscReal       shift = info->shiftamount;
373   PetscBool       allowzeropivot, zeropivotdetected;
374 
375   PetscFunctionBegin;
376   allowzeropivot = PetscNot(A->erroriffailure);
377   PetscCall(ISGetIndices(isrow, &r));
378   PetscCall(ISGetIndices(isicol, &ic));
379 
380   /* generate work space needed by the factorization */
381   PetscCall(PetscMalloc2(bs2 * n, &rtmp, bs2, &mwork));
382   PetscCall(PetscArrayzero(rtmp, bs2 * n));
383 
384   for (i = 0; i < n; i++) {
385     /* zero rtmp */
386     /* L part */
387     nz    = bi[i + 1] - bi[i];
388     bjtmp = bj + bi[i];
389     for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2));
390 
391     /* U part */
392     nz    = bdiag[i] - bdiag[i + 1];
393     bjtmp = bj + bdiag[i + 1] + 1;
394     for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2));
395 
396     /* load in initial (unfactored row) */
397     nz    = ai[r[i] + 1] - ai[r[i]];
398     ajtmp = aj + ai[r[i]];
399     v     = aa + bs2 * ai[r[i]];
400     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(rtmp + bs2 * ic[ajtmp[j]], v + bs2 * j, bs2));
401 
402     /* elimination */
403     bjtmp = bj + bi[i];
404     nzL   = bi[i + 1] - bi[i];
405     for (k = 0; k < nzL; k++) {
406       row = bjtmp[k];
407       pc  = rtmp + bs2 * row;
408       for (flg = 0, j = 0; j < bs2; j++) {
409         if (pc[j] != 0.0) {
410           flg = 1;
411           break;
412         }
413       }
414       if (flg) {
415         pv = b->a + bs2 * bdiag[row];
416         /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
417         PetscCall(PetscKernel_A_gets_A_times_B_6(pc, pv, mwork));
418 
419         pj = b->j + bdiag[row + 1] + 1; /* beginning of U(row,:) */
420         pv = b->a + bs2 * (bdiag[row + 1] + 1);
421         nz = bdiag[row] - bdiag[row + 1] - 1; /* num of entries inU(row,:), excluding diag */
422         for (j = 0; j < nz; j++) {
423           /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
424           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
425           v = rtmp + bs2 * pj[j];
426           PetscCall(PetscKernel_A_gets_A_minus_B_times_C_6(v, pc, pv));
427           pv += bs2;
428         }
429         PetscCall(PetscLogFlops(432.0 * nz + 396)); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
430       }
431     }
432 
433     /* finished row so stick it into b->a */
434     /* L part */
435     pv = b->a + bs2 * bi[i];
436     pj = b->j + bi[i];
437     nz = bi[i + 1] - bi[i];
438     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2));
439 
440     /* Mark diagonal and invert diagonal for simpler triangular solves */
441     pv = b->a + bs2 * bdiag[i];
442     pj = b->j + bdiag[i];
443     PetscCall(PetscArraycpy(pv, rtmp + bs2 * pj[0], bs2));
444     PetscCall(PetscKernel_A_gets_inverse_A_6(pv, shift, allowzeropivot, &zeropivotdetected));
445     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
446 
447     /* U part */
448     pv = b->a + bs2 * (bdiag[i + 1] + 1);
449     pj = b->j + bdiag[i + 1] + 1;
450     nz = bdiag[i] - bdiag[i + 1] - 1;
451     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2));
452   }
453 
454   PetscCall(PetscFree2(rtmp, mwork));
455   PetscCall(ISRestoreIndices(isicol, &ic));
456   PetscCall(ISRestoreIndices(isrow, &r));
457 
458   C->ops->solve          = MatSolve_SeqBAIJ_6;
459   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_6;
460   C->assembled           = PETSC_TRUE;
461 
462   PetscCall(PetscLogFlops(1.333333333333 * 6 * 6 * 6 * n)); /* from inverting diagonal blocks */
463   PetscFunctionReturn(PETSC_SUCCESS);
464 }
465 
466 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering_inplace(Mat C, Mat A, const MatFactorInfo *info)
467 {
468   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data;
469   PetscInt     i, j, n = a->mbs, *bi = b->i, *bj = b->j;
470   PetscInt    *ajtmpold, *ajtmp, nz, row;
471   PetscInt    *diag_offset = b->diag, *ai = a->i, *aj = a->j, *pj;
472   MatScalar   *pv, *v, *rtmp, *pc, *w, *x;
473   MatScalar    x1, x2, x3, x4, x5, x6, x7, x8, x9, x10, x11, x12, x13, x14, x15;
474   MatScalar    x16, x17, x18, x19, x20, x21, x22, x23, x24, x25;
475   MatScalar    p1, p2, p3, p4, p5, p6, p7, p8, p9, p10, p11, p12, p13, p14, p15;
476   MatScalar    p16, p17, p18, p19, p20, p21, p22, p23, p24, p25;
477   MatScalar    m1, m2, m3, m4, m5, m6, m7, m8, m9, m10, m11, m12, m13, m14, m15;
478   MatScalar    m16, m17, m18, m19, m20, m21, m22, m23, m24, m25;
479   MatScalar    p26, p27, p28, p29, p30, p31, p32, p33, p34, p35, p36;
480   MatScalar    x26, x27, x28, x29, x30, x31, x32, x33, x34, x35, x36;
481   MatScalar    m26, m27, m28, m29, m30, m31, m32, m33, m34, m35, m36;
482   MatScalar   *ba = b->a, *aa = a->a;
483   PetscReal    shift = info->shiftamount;
484   PetscBool    allowzeropivot, zeropivotdetected;
485 
486   PetscFunctionBegin;
487   allowzeropivot = PetscNot(A->erroriffailure);
488   PetscCall(PetscMalloc1(36 * (n + 1), &rtmp));
489   for (i = 0; i < n; i++) {
490     nz    = bi[i + 1] - bi[i];
491     ajtmp = bj + bi[i];
492     for (j = 0; j < nz; j++) {
493       x    = rtmp + 36 * ajtmp[j];
494       x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = x[9] = 0.0;
495       x[10] = x[11] = x[12] = x[13] = x[14] = x[15] = x[16] = x[17] = 0.0;
496       x[18] = x[19] = x[20] = x[21] = x[22] = x[23] = x[24] = x[25] = 0.0;
497       x[26] = x[27] = x[28] = x[29] = x[30] = x[31] = x[32] = x[33] = 0.0;
498       x[34] = x[35] = 0.0;
499     }
500     /* load in initial (unfactored row) */
501     nz       = ai[i + 1] - ai[i];
502     ajtmpold = aj + ai[i];
503     v        = aa + 36 * ai[i];
504     for (j = 0; j < nz; j++) {
505       x     = rtmp + 36 * ajtmpold[j];
506       x[0]  = v[0];
507       x[1]  = v[1];
508       x[2]  = v[2];
509       x[3]  = v[3];
510       x[4]  = v[4];
511       x[5]  = v[5];
512       x[6]  = v[6];
513       x[7]  = v[7];
514       x[8]  = v[8];
515       x[9]  = v[9];
516       x[10] = v[10];
517       x[11] = v[11];
518       x[12] = v[12];
519       x[13] = v[13];
520       x[14] = v[14];
521       x[15] = v[15];
522       x[16] = v[16];
523       x[17] = v[17];
524       x[18] = v[18];
525       x[19] = v[19];
526       x[20] = v[20];
527       x[21] = v[21];
528       x[22] = v[22];
529       x[23] = v[23];
530       x[24] = v[24];
531       x[25] = v[25];
532       x[26] = v[26];
533       x[27] = v[27];
534       x[28] = v[28];
535       x[29] = v[29];
536       x[30] = v[30];
537       x[31] = v[31];
538       x[32] = v[32];
539       x[33] = v[33];
540       x[34] = v[34];
541       x[35] = v[35];
542       v += 36;
543     }
544     row = *ajtmp++;
545     while (row < i) {
546       pc  = rtmp + 36 * row;
547       p1  = pc[0];
548       p2  = pc[1];
549       p3  = pc[2];
550       p4  = pc[3];
551       p5  = pc[4];
552       p6  = pc[5];
553       p7  = pc[6];
554       p8  = pc[7];
555       p9  = pc[8];
556       p10 = pc[9];
557       p11 = pc[10];
558       p12 = pc[11];
559       p13 = pc[12];
560       p14 = pc[13];
561       p15 = pc[14];
562       p16 = pc[15];
563       p17 = pc[16];
564       p18 = pc[17];
565       p19 = pc[18];
566       p20 = pc[19];
567       p21 = pc[20];
568       p22 = pc[21];
569       p23 = pc[22];
570       p24 = pc[23];
571       p25 = pc[24];
572       p26 = pc[25];
573       p27 = pc[26];
574       p28 = pc[27];
575       p29 = pc[28];
576       p30 = pc[29];
577       p31 = pc[30];
578       p32 = pc[31];
579       p33 = pc[32];
580       p34 = pc[33];
581       p35 = pc[34];
582       p36 = pc[35];
583       if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 || p5 != 0.0 || p6 != 0.0 || p7 != 0.0 || p8 != 0.0 || p9 != 0.0 || p10 != 0.0 || p11 != 0.0 || p12 != 0.0 || p13 != 0.0 || p14 != 0.0 || p15 != 0.0 || p16 != 0.0 || p17 != 0.0 || p18 != 0.0 || p19 != 0.0 || p20 != 0.0 || p21 != 0.0 || p22 != 0.0 || p23 != 0.0 || p24 != 0.0 || p25 != 0.0 || p26 != 0.0 || p27 != 0.0 || p28 != 0.0 || p29 != 0.0 || p30 != 0.0 || p31 != 0.0 || p32 != 0.0 || p33 != 0.0 || p34 != 0.0 || p35 != 0.0 || p36 != 0.0) {
584         pv    = ba + 36 * diag_offset[row];
585         pj    = bj + diag_offset[row] + 1;
586         x1    = pv[0];
587         x2    = pv[1];
588         x3    = pv[2];
589         x4    = pv[3];
590         x5    = pv[4];
591         x6    = pv[5];
592         x7    = pv[6];
593         x8    = pv[7];
594         x9    = pv[8];
595         x10   = pv[9];
596         x11   = pv[10];
597         x12   = pv[11];
598         x13   = pv[12];
599         x14   = pv[13];
600         x15   = pv[14];
601         x16   = pv[15];
602         x17   = pv[16];
603         x18   = pv[17];
604         x19   = pv[18];
605         x20   = pv[19];
606         x21   = pv[20];
607         x22   = pv[21];
608         x23   = pv[22];
609         x24   = pv[23];
610         x25   = pv[24];
611         x26   = pv[25];
612         x27   = pv[26];
613         x28   = pv[27];
614         x29   = pv[28];
615         x30   = pv[29];
616         x31   = pv[30];
617         x32   = pv[31];
618         x33   = pv[32];
619         x34   = pv[33];
620         x35   = pv[34];
621         x36   = pv[35];
622         pc[0] = m1 = p1 * x1 + p7 * x2 + p13 * x3 + p19 * x4 + p25 * x5 + p31 * x6;
623         pc[1] = m2 = p2 * x1 + p8 * x2 + p14 * x3 + p20 * x4 + p26 * x5 + p32 * x6;
624         pc[2] = m3 = p3 * x1 + p9 * x2 + p15 * x3 + p21 * x4 + p27 * x5 + p33 * x6;
625         pc[3] = m4 = p4 * x1 + p10 * x2 + p16 * x3 + p22 * x4 + p28 * x5 + p34 * x6;
626         pc[4] = m5 = p5 * x1 + p11 * x2 + p17 * x3 + p23 * x4 + p29 * x5 + p35 * x6;
627         pc[5] = m6 = p6 * x1 + p12 * x2 + p18 * x3 + p24 * x4 + p30 * x5 + p36 * x6;
628 
629         pc[6] = m7 = p1 * x7 + p7 * x8 + p13 * x9 + p19 * x10 + p25 * x11 + p31 * x12;
630         pc[7] = m8 = p2 * x7 + p8 * x8 + p14 * x9 + p20 * x10 + p26 * x11 + p32 * x12;
631         pc[8] = m9 = p3 * x7 + p9 * x8 + p15 * x9 + p21 * x10 + p27 * x11 + p33 * x12;
632         pc[9] = m10 = p4 * x7 + p10 * x8 + p16 * x9 + p22 * x10 + p28 * x11 + p34 * x12;
633         pc[10] = m11 = p5 * x7 + p11 * x8 + p17 * x9 + p23 * x10 + p29 * x11 + p35 * x12;
634         pc[11] = m12 = p6 * x7 + p12 * x8 + p18 * x9 + p24 * x10 + p30 * x11 + p36 * x12;
635 
636         pc[12] = m13 = p1 * x13 + p7 * x14 + p13 * x15 + p19 * x16 + p25 * x17 + p31 * x18;
637         pc[13] = m14 = p2 * x13 + p8 * x14 + p14 * x15 + p20 * x16 + p26 * x17 + p32 * x18;
638         pc[14] = m15 = p3 * x13 + p9 * x14 + p15 * x15 + p21 * x16 + p27 * x17 + p33 * x18;
639         pc[15] = m16 = p4 * x13 + p10 * x14 + p16 * x15 + p22 * x16 + p28 * x17 + p34 * x18;
640         pc[16] = m17 = p5 * x13 + p11 * x14 + p17 * x15 + p23 * x16 + p29 * x17 + p35 * x18;
641         pc[17] = m18 = p6 * x13 + p12 * x14 + p18 * x15 + p24 * x16 + p30 * x17 + p36 * x18;
642 
643         pc[18] = m19 = p1 * x19 + p7 * x20 + p13 * x21 + p19 * x22 + p25 * x23 + p31 * x24;
644         pc[19] = m20 = p2 * x19 + p8 * x20 + p14 * x21 + p20 * x22 + p26 * x23 + p32 * x24;
645         pc[20] = m21 = p3 * x19 + p9 * x20 + p15 * x21 + p21 * x22 + p27 * x23 + p33 * x24;
646         pc[21] = m22 = p4 * x19 + p10 * x20 + p16 * x21 + p22 * x22 + p28 * x23 + p34 * x24;
647         pc[22] = m23 = p5 * x19 + p11 * x20 + p17 * x21 + p23 * x22 + p29 * x23 + p35 * x24;
648         pc[23] = m24 = p6 * x19 + p12 * x20 + p18 * x21 + p24 * x22 + p30 * x23 + p36 * x24;
649 
650         pc[24] = m25 = p1 * x25 + p7 * x26 + p13 * x27 + p19 * x28 + p25 * x29 + p31 * x30;
651         pc[25] = m26 = p2 * x25 + p8 * x26 + p14 * x27 + p20 * x28 + p26 * x29 + p32 * x30;
652         pc[26] = m27 = p3 * x25 + p9 * x26 + p15 * x27 + p21 * x28 + p27 * x29 + p33 * x30;
653         pc[27] = m28 = p4 * x25 + p10 * x26 + p16 * x27 + p22 * x28 + p28 * x29 + p34 * x30;
654         pc[28] = m29 = p5 * x25 + p11 * x26 + p17 * x27 + p23 * x28 + p29 * x29 + p35 * x30;
655         pc[29] = m30 = p6 * x25 + p12 * x26 + p18 * x27 + p24 * x28 + p30 * x29 + p36 * x30;
656 
657         pc[30] = m31 = p1 * x31 + p7 * x32 + p13 * x33 + p19 * x34 + p25 * x35 + p31 * x36;
658         pc[31] = m32 = p2 * x31 + p8 * x32 + p14 * x33 + p20 * x34 + p26 * x35 + p32 * x36;
659         pc[32] = m33 = p3 * x31 + p9 * x32 + p15 * x33 + p21 * x34 + p27 * x35 + p33 * x36;
660         pc[33] = m34 = p4 * x31 + p10 * x32 + p16 * x33 + p22 * x34 + p28 * x35 + p34 * x36;
661         pc[34] = m35 = p5 * x31 + p11 * x32 + p17 * x33 + p23 * x34 + p29 * x35 + p35 * x36;
662         pc[35] = m36 = p6 * x31 + p12 * x32 + p18 * x33 + p24 * x34 + p30 * x35 + p36 * x36;
663 
664         nz = bi[row + 1] - diag_offset[row] - 1;
665         pv += 36;
666         for (j = 0; j < nz; j++) {
667           x1  = pv[0];
668           x2  = pv[1];
669           x3  = pv[2];
670           x4  = pv[3];
671           x5  = pv[4];
672           x6  = pv[5];
673           x7  = pv[6];
674           x8  = pv[7];
675           x9  = pv[8];
676           x10 = pv[9];
677           x11 = pv[10];
678           x12 = pv[11];
679           x13 = pv[12];
680           x14 = pv[13];
681           x15 = pv[14];
682           x16 = pv[15];
683           x17 = pv[16];
684           x18 = pv[17];
685           x19 = pv[18];
686           x20 = pv[19];
687           x21 = pv[20];
688           x22 = pv[21];
689           x23 = pv[22];
690           x24 = pv[23];
691           x25 = pv[24];
692           x26 = pv[25];
693           x27 = pv[26];
694           x28 = pv[27];
695           x29 = pv[28];
696           x30 = pv[29];
697           x31 = pv[30];
698           x32 = pv[31];
699           x33 = pv[32];
700           x34 = pv[33];
701           x35 = pv[34];
702           x36 = pv[35];
703           x   = rtmp + 36 * pj[j];
704           x[0] -= m1 * x1 + m7 * x2 + m13 * x3 + m19 * x4 + m25 * x5 + m31 * x6;
705           x[1] -= m2 * x1 + m8 * x2 + m14 * x3 + m20 * x4 + m26 * x5 + m32 * x6;
706           x[2] -= m3 * x1 + m9 * x2 + m15 * x3 + m21 * x4 + m27 * x5 + m33 * x6;
707           x[3] -= m4 * x1 + m10 * x2 + m16 * x3 + m22 * x4 + m28 * x5 + m34 * x6;
708           x[4] -= m5 * x1 + m11 * x2 + m17 * x3 + m23 * x4 + m29 * x5 + m35 * x6;
709           x[5] -= m6 * x1 + m12 * x2 + m18 * x3 + m24 * x4 + m30 * x5 + m36 * x6;
710 
711           x[6] -= m1 * x7 + m7 * x8 + m13 * x9 + m19 * x10 + m25 * x11 + m31 * x12;
712           x[7] -= m2 * x7 + m8 * x8 + m14 * x9 + m20 * x10 + m26 * x11 + m32 * x12;
713           x[8] -= m3 * x7 + m9 * x8 + m15 * x9 + m21 * x10 + m27 * x11 + m33 * x12;
714           x[9] -= m4 * x7 + m10 * x8 + m16 * x9 + m22 * x10 + m28 * x11 + m34 * x12;
715           x[10] -= m5 * x7 + m11 * x8 + m17 * x9 + m23 * x10 + m29 * x11 + m35 * x12;
716           x[11] -= m6 * x7 + m12 * x8 + m18 * x9 + m24 * x10 + m30 * x11 + m36 * x12;
717 
718           x[12] -= m1 * x13 + m7 * x14 + m13 * x15 + m19 * x16 + m25 * x17 + m31 * x18;
719           x[13] -= m2 * x13 + m8 * x14 + m14 * x15 + m20 * x16 + m26 * x17 + m32 * x18;
720           x[14] -= m3 * x13 + m9 * x14 + m15 * x15 + m21 * x16 + m27 * x17 + m33 * x18;
721           x[15] -= m4 * x13 + m10 * x14 + m16 * x15 + m22 * x16 + m28 * x17 + m34 * x18;
722           x[16] -= m5 * x13 + m11 * x14 + m17 * x15 + m23 * x16 + m29 * x17 + m35 * x18;
723           x[17] -= m6 * x13 + m12 * x14 + m18 * x15 + m24 * x16 + m30 * x17 + m36 * x18;
724 
725           x[18] -= m1 * x19 + m7 * x20 + m13 * x21 + m19 * x22 + m25 * x23 + m31 * x24;
726           x[19] -= m2 * x19 + m8 * x20 + m14 * x21 + m20 * x22 + m26 * x23 + m32 * x24;
727           x[20] -= m3 * x19 + m9 * x20 + m15 * x21 + m21 * x22 + m27 * x23 + m33 * x24;
728           x[21] -= m4 * x19 + m10 * x20 + m16 * x21 + m22 * x22 + m28 * x23 + m34 * x24;
729           x[22] -= m5 * x19 + m11 * x20 + m17 * x21 + m23 * x22 + m29 * x23 + m35 * x24;
730           x[23] -= m6 * x19 + m12 * x20 + m18 * x21 + m24 * x22 + m30 * x23 + m36 * x24;
731 
732           x[24] -= m1 * x25 + m7 * x26 + m13 * x27 + m19 * x28 + m25 * x29 + m31 * x30;
733           x[25] -= m2 * x25 + m8 * x26 + m14 * x27 + m20 * x28 + m26 * x29 + m32 * x30;
734           x[26] -= m3 * x25 + m9 * x26 + m15 * x27 + m21 * x28 + m27 * x29 + m33 * x30;
735           x[27] -= m4 * x25 + m10 * x26 + m16 * x27 + m22 * x28 + m28 * x29 + m34 * x30;
736           x[28] -= m5 * x25 + m11 * x26 + m17 * x27 + m23 * x28 + m29 * x29 + m35 * x30;
737           x[29] -= m6 * x25 + m12 * x26 + m18 * x27 + m24 * x28 + m30 * x29 + m36 * x30;
738 
739           x[30] -= m1 * x31 + m7 * x32 + m13 * x33 + m19 * x34 + m25 * x35 + m31 * x36;
740           x[31] -= m2 * x31 + m8 * x32 + m14 * x33 + m20 * x34 + m26 * x35 + m32 * x36;
741           x[32] -= m3 * x31 + m9 * x32 + m15 * x33 + m21 * x34 + m27 * x35 + m33 * x36;
742           x[33] -= m4 * x31 + m10 * x32 + m16 * x33 + m22 * x34 + m28 * x35 + m34 * x36;
743           x[34] -= m5 * x31 + m11 * x32 + m17 * x33 + m23 * x34 + m29 * x35 + m35 * x36;
744           x[35] -= m6 * x31 + m12 * x32 + m18 * x33 + m24 * x34 + m30 * x35 + m36 * x36;
745 
746           pv += 36;
747         }
748         PetscCall(PetscLogFlops(432.0 * nz + 396.0));
749       }
750       row = *ajtmp++;
751     }
752     /* finished row so stick it into b->a */
753     pv = ba + 36 * bi[i];
754     pj = bj + bi[i];
755     nz = bi[i + 1] - bi[i];
756     for (j = 0; j < nz; j++) {
757       x      = rtmp + 36 * pj[j];
758       pv[0]  = x[0];
759       pv[1]  = x[1];
760       pv[2]  = x[2];
761       pv[3]  = x[3];
762       pv[4]  = x[4];
763       pv[5]  = x[5];
764       pv[6]  = x[6];
765       pv[7]  = x[7];
766       pv[8]  = x[8];
767       pv[9]  = x[9];
768       pv[10] = x[10];
769       pv[11] = x[11];
770       pv[12] = x[12];
771       pv[13] = x[13];
772       pv[14] = x[14];
773       pv[15] = x[15];
774       pv[16] = x[16];
775       pv[17] = x[17];
776       pv[18] = x[18];
777       pv[19] = x[19];
778       pv[20] = x[20];
779       pv[21] = x[21];
780       pv[22] = x[22];
781       pv[23] = x[23];
782       pv[24] = x[24];
783       pv[25] = x[25];
784       pv[26] = x[26];
785       pv[27] = x[27];
786       pv[28] = x[28];
787       pv[29] = x[29];
788       pv[30] = x[30];
789       pv[31] = x[31];
790       pv[32] = x[32];
791       pv[33] = x[33];
792       pv[34] = x[34];
793       pv[35] = x[35];
794       pv += 36;
795     }
796     /* invert diagonal block */
797     w = ba + 36 * diag_offset[i];
798     PetscCall(PetscKernel_A_gets_inverse_A_6(w, shift, allowzeropivot, &zeropivotdetected));
799     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
800   }
801 
802   PetscCall(PetscFree(rtmp));
803 
804   C->ops->solve          = MatSolve_SeqBAIJ_6_NaturalOrdering_inplace;
805   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_6_NaturalOrdering_inplace;
806   C->assembled           = PETSC_TRUE;
807 
808   PetscCall(PetscLogFlops(1.333333333333 * 6 * 6 * 6 * b->mbs)); /* from inverting diagonal blocks */
809   PetscFunctionReturn(PETSC_SUCCESS);
810 }
811 
812 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering(Mat B, Mat A, const MatFactorInfo *info)
813 {
814   Mat             C = B;
815   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data;
816   PetscInt        i, j, k, nz, nzL, row;
817   const PetscInt  n = a->mbs, *ai = a->i, *aj = a->j, *bi = b->i, *bj = b->j;
818   const PetscInt *ajtmp, *bjtmp, *bdiag = b->diag, *pj, bs2 = a->bs2;
819   MatScalar      *rtmp, *pc, *mwork, *v, *pv, *aa = a->a;
820   PetscInt        flg;
821   PetscReal       shift = info->shiftamount;
822   PetscBool       allowzeropivot, zeropivotdetected;
823 
824   PetscFunctionBegin;
825   allowzeropivot = PetscNot(A->erroriffailure);
826 
827   /* generate work space needed by the factorization */
828   PetscCall(PetscMalloc2(bs2 * n, &rtmp, bs2, &mwork));
829   PetscCall(PetscArrayzero(rtmp, bs2 * n));
830 
831   for (i = 0; i < n; i++) {
832     /* zero rtmp */
833     /* L part */
834     nz    = bi[i + 1] - bi[i];
835     bjtmp = bj + bi[i];
836     for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2));
837 
838     /* U part */
839     nz    = bdiag[i] - bdiag[i + 1];
840     bjtmp = bj + bdiag[i + 1] + 1;
841     for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2));
842 
843     /* load in initial (unfactored row) */
844     nz    = ai[i + 1] - ai[i];
845     ajtmp = aj + ai[i];
846     v     = aa + bs2 * ai[i];
847     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(rtmp + bs2 * ajtmp[j], v + bs2 * j, bs2));
848 
849     /* elimination */
850     bjtmp = bj + bi[i];
851     nzL   = bi[i + 1] - bi[i];
852     for (k = 0; k < nzL; k++) {
853       row = bjtmp[k];
854       pc  = rtmp + bs2 * row;
855       for (flg = 0, j = 0; j < bs2; j++) {
856         if (pc[j] != 0.0) {
857           flg = 1;
858           break;
859         }
860       }
861       if (flg) {
862         pv = b->a + bs2 * bdiag[row];
863         /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
864         PetscCall(PetscKernel_A_gets_A_times_B_6(pc, pv, mwork));
865 
866         pj = b->j + bdiag[row + 1] + 1; /* beginning of U(row,:) */
867         pv = b->a + bs2 * (bdiag[row + 1] + 1);
868         nz = bdiag[row] - bdiag[row + 1] - 1; /* num of entries inU(row,:), excluding diag */
869         for (j = 0; j < nz; j++) {
870           /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
871           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
872           v = rtmp + bs2 * pj[j];
873           PetscCall(PetscKernel_A_gets_A_minus_B_times_C_6(v, pc, pv));
874           pv += bs2;
875         }
876         PetscCall(PetscLogFlops(432.0 * nz + 396)); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
877       }
878     }
879 
880     /* finished row so stick it into b->a */
881     /* L part */
882     pv = b->a + bs2 * bi[i];
883     pj = b->j + bi[i];
884     nz = bi[i + 1] - bi[i];
885     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2));
886 
887     /* Mark diagonal and invert diagonal for simpler triangular solves */
888     pv = b->a + bs2 * bdiag[i];
889     pj = b->j + bdiag[i];
890     PetscCall(PetscArraycpy(pv, rtmp + bs2 * pj[0], bs2));
891     PetscCall(PetscKernel_A_gets_inverse_A_6(pv, shift, allowzeropivot, &zeropivotdetected));
892     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
893 
894     /* U part */
895     pv = b->a + bs2 * (bdiag[i + 1] + 1);
896     pj = b->j + bdiag[i + 1] + 1;
897     nz = bdiag[i] - bdiag[i + 1] - 1;
898     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2));
899   }
900   PetscCall(PetscFree2(rtmp, mwork));
901 
902   C->ops->solve          = MatSolve_SeqBAIJ_6_NaturalOrdering;
903   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_6_NaturalOrdering;
904   C->assembled           = PETSC_TRUE;
905 
906   PetscCall(PetscLogFlops(1.333333333333 * 6 * 6 * 6 * n)); /* from inverting diagonal blocks */
907   PetscFunctionReturn(PETSC_SUCCESS);
908 }
909