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