xref: /petsc/src/mat/impls/baij/seq/baijfact7.c (revision 503c0ea9b45bcfbcebbb1ea5341243bbc69f0bea)
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];  x[1] =  v[1];  x[2] =  v[2];  x[3] =  v[3];
56       x[4]  =  v[4];  x[5] =  v[5];  x[6] =  v[6];  x[7] =  v[7];
57       x[8]  =  v[8];  x[9] =  v[9];  x[10] = v[10]; x[11] = v[11];
58       x[12] = v[12]; x[13] = v[13]; x[14] = v[14]; x[15] = v[15];
59       x[16] = v[16]; x[17] = v[17]; x[18] = v[18]; x[19] = v[19];
60       x[20] = v[20]; x[21] = v[21]; x[22] = v[22]; x[23] = v[23];
61       x[24] = v[24]; x[25] = v[25]; x[26] = v[26]; x[27] = v[27];
62       x[28] = v[28]; x[29] = v[29]; x[30] = v[30]; x[31] = v[31];
63       x[32] = v[32]; x[33] = v[33]; x[34] = v[34]; x[35] = v[35];
64       v    += 36;
65     }
66     row = *ajtmp++;
67     while (row < i) {
68       pc  =  rtmp + 36*row;
69       p1  = pc[0];  p2  = pc[1];  p3  = pc[2];  p4  = pc[3];
70       p5  = pc[4];  p6  = pc[5];  p7  = pc[6];  p8  = pc[7];
71       p9  = pc[8];  p10 = pc[9];  p11 = pc[10]; p12 = pc[11];
72       p13 = pc[12]; p14 = pc[13]; p15 = pc[14]; p16 = pc[15];
73       p17 = pc[16]; p18 = pc[17]; p19 = pc[18]; p20 = pc[19];
74       p21 = pc[20]; p22 = pc[21]; p23 = pc[22]; p24 = pc[23];
75       p25 = pc[24]; p26 = pc[25]; p27 = pc[26]; p28 = pc[27];
76       p29 = pc[28]; p30 = pc[29]; p31 = pc[30]; p32 = pc[31];
77       p33 = pc[32]; p34 = pc[33]; p35 = pc[34]; p36 = pc[35];
78       if (p1  != 0.0 || p2  != 0.0 || p3  != 0.0 || p4  != 0.0 ||
79           p5  != 0.0 || p6  != 0.0 || p7  != 0.0 || p8  != 0.0 ||
80           p9  != 0.0 || p10 != 0.0 || p11 != 0.0 || p12 != 0.0 ||
81           p13 != 0.0 || p14 != 0.0 || p15 != 0.0 || p16 != 0.0 ||
82           p17 != 0.0 || p18 != 0.0 || p19 != 0.0 || p20 != 0.0 ||
83           p21 != 0.0 || p22 != 0.0 || p23 != 0.0 || p24 != 0.0 ||
84           p25 != 0.0 || p26 != 0.0 || p27 != 0.0 || p28 != 0.0 ||
85           p29 != 0.0 || p30 != 0.0 || p31 != 0.0 || p32 != 0.0 ||
86           p33 != 0.0 || p34 != 0.0 || p35 != 0.0 || p36 != 0.0) {
87         pv    = ba + 36*diag_offset[row];
88         pj    = bj + diag_offset[row] + 1;
89         x1    = pv[0];  x2  = pv[1];  x3  = pv[2];  x4  = pv[3];
90         x5    = pv[4];  x6  = pv[5];  x7  = pv[6];  x8  = pv[7];
91         x9    = pv[8];  x10 = pv[9];  x11 = pv[10]; x12 = pv[11];
92         x13   = pv[12]; x14 = pv[13]; x15 = pv[14]; x16 = pv[15];
93         x17   = pv[16]; x18 = pv[17]; x19 = pv[18]; x20 = pv[19];
94         x21   = pv[20]; x22 = pv[21]; x23 = pv[22]; x24 = pv[23];
95         x25   = pv[24]; x26 = pv[25]; x27 = pv[26]; x28 = pv[27];
96         x29   = pv[28]; x30 = pv[29]; x31 = pv[30]; x32 = pv[31];
97         x33   = pv[32]; x34 = pv[33]; x35 = pv[34]; x36 = pv[35];
98         pc[0] = m1  = p1*x1  + p7*x2   + p13*x3  + p19*x4  + p25*x5  + p31*x6;
99         pc[1] = m2  = p2*x1  + p8*x2   + p14*x3  + p20*x4  + p26*x5  + p32*x6;
100         pc[2] = m3  = p3*x1  + p9*x2   + p15*x3  + p21*x4  + p27*x5  + p33*x6;
101         pc[3] = m4  = p4*x1  + p10*x2  + p16*x3  + p22*x4  + p28*x5  + p34*x6;
102         pc[4] = m5  = p5*x1  + p11*x2  + p17*x3  + p23*x4  + p29*x5  + p35*x6;
103         pc[5] = m6  = p6*x1  + p12*x2  + p18*x3  + p24*x4  + p30*x5  + p36*x6;
104 
105         pc[6]  = m7  = p1*x7  + p7*x8   + p13*x9  + p19*x10 + p25*x11 + p31*x12;
106         pc[7]  = m8  = p2*x7  + p8*x8   + p14*x9  + p20*x10 + p26*x11 + p32*x12;
107         pc[8]  = m9  = p3*x7  + p9*x8   + p15*x9  + p21*x10 + p27*x11 + p33*x12;
108         pc[9]  = m10 = p4*x7  + p10*x8  + p16*x9  + p22*x10 + p28*x11 + p34*x12;
109         pc[10] = m11 = p5*x7  + p11*x8  + p17*x9  + p23*x10 + p29*x11 + p35*x12;
110         pc[11] = m12 = p6*x7  + p12*x8  + p18*x9  + p24*x10 + p30*x11 + p36*x12;
111 
112         pc[12] = m13 = p1*x13 + p7*x14  + p13*x15 + p19*x16 + p25*x17 + p31*x18;
113         pc[13] = m14 = p2*x13 + p8*x14  + p14*x15 + p20*x16 + p26*x17 + p32*x18;
114         pc[14] = m15 = p3*x13 + p9*x14  + p15*x15 + p21*x16 + p27*x17 + p33*x18;
115         pc[15] = m16 = p4*x13 + p10*x14 + p16*x15 + p22*x16 + p28*x17 + p34*x18;
116         pc[16] = m17 = p5*x13 + p11*x14 + p17*x15 + p23*x16 + p29*x17 + p35*x18;
117         pc[17] = m18 = p6*x13 + p12*x14 + p18*x15 + p24*x16 + p30*x17 + p36*x18;
118 
119         pc[18] = m19 = p1*x19 + p7*x20  + p13*x21 + p19*x22 + p25*x23 + p31*x24;
120         pc[19] = m20 = p2*x19 + p8*x20  + p14*x21 + p20*x22 + p26*x23 + p32*x24;
121         pc[20] = m21 = p3*x19 + p9*x20  + p15*x21 + p21*x22 + p27*x23 + p33*x24;
122         pc[21] = m22 = p4*x19 + p10*x20 + p16*x21 + p22*x22 + p28*x23 + p34*x24;
123         pc[22] = m23 = p5*x19 + p11*x20 + p17*x21 + p23*x22 + p29*x23 + p35*x24;
124         pc[23] = m24 = p6*x19 + p12*x20 + p18*x21 + p24*x22 + p30*x23 + p36*x24;
125 
126         pc[24] = m25 = p1*x25 + p7*x26  + p13*x27 + p19*x28 + p25*x29 + p31*x30;
127         pc[25] = m26 = p2*x25 + p8*x26  + p14*x27 + p20*x28 + p26*x29 + p32*x30;
128         pc[26] = m27 = p3*x25 + p9*x26  + p15*x27 + p21*x28 + p27*x29 + p33*x30;
129         pc[27] = m28 = p4*x25 + p10*x26 + p16*x27 + p22*x28 + p28*x29 + p34*x30;
130         pc[28] = m29 = p5*x25 + p11*x26 + p17*x27 + p23*x28 + p29*x29 + p35*x30;
131         pc[29] = m30 = p6*x25 + p12*x26 + p18*x27 + p24*x28 + p30*x29 + p36*x30;
132 
133         pc[30] = m31 = p1*x31 + p7*x32  + p13*x33 + p19*x34 + p25*x35 + p31*x36;
134         pc[31] = m32 = p2*x31 + p8*x32  + p14*x33 + p20*x34 + p26*x35 + p32*x36;
135         pc[32] = m33 = p3*x31 + p9*x32  + p15*x33 + p21*x34 + p27*x35 + p33*x36;
136         pc[33] = m34 = p4*x31 + p10*x32 + p16*x33 + p22*x34 + p28*x35 + p34*x36;
137         pc[34] = m35 = p5*x31 + p11*x32 + p17*x33 + p23*x34 + p29*x35 + p35*x36;
138         pc[35] = m36 = p6*x31 + p12*x32 + p18*x33 + p24*x34 + p30*x35 + p36*x36;
139 
140         nz  = bi[row+1] - diag_offset[row] - 1;
141         pv += 36;
142         for (j=0; j<nz; j++) {
143           x1    = pv[0];  x2  = pv[1];  x3  = pv[2];  x4  = pv[3];
144           x5    = pv[4];  x6  = pv[5];  x7  = pv[6];  x8  = pv[7];
145           x9    = pv[8];  x10 = pv[9];  x11 = pv[10]; x12 = pv[11];
146           x13   = pv[12]; x14 = pv[13]; x15 = pv[14]; x16 = pv[15];
147           x17   = pv[16]; x18 = pv[17]; x19 = pv[18]; x20 = pv[19];
148           x21   = pv[20]; x22 = pv[21]; x23 = pv[22]; x24 = pv[23];
149           x25   = pv[24]; x26 = pv[25]; x27 = pv[26]; x28 = pv[27];
150           x29   = pv[28]; x30 = pv[29]; x31 = pv[30]; x32 = pv[31];
151           x33   = pv[32]; x34 = pv[33]; x35 = pv[34]; x36 = pv[35];
152           x     = rtmp + 36*pj[j];
153           x[0] -= m1*x1  + m7*x2   + m13*x3  + m19*x4  + m25*x5  + m31*x6;
154           x[1] -= m2*x1  + m8*x2   + m14*x3  + m20*x4  + m26*x5  + m32*x6;
155           x[2] -= m3*x1  + m9*x2   + m15*x3  + m21*x4  + m27*x5  + m33*x6;
156           x[3] -= m4*x1  + m10*x2  + m16*x3  + m22*x4  + m28*x5  + m34*x6;
157           x[4] -= m5*x1  + m11*x2  + m17*x3  + m23*x4  + m29*x5  + m35*x6;
158           x[5] -= m6*x1  + m12*x2  + m18*x3  + m24*x4  + m30*x5  + m36*x6;
159 
160           x[6]  -= m1*x7  + m7*x8   + m13*x9  + m19*x10 + m25*x11 + m31*x12;
161           x[7]  -= m2*x7  + m8*x8   + m14*x9  + m20*x10 + m26*x11 + m32*x12;
162           x[8]  -= m3*x7  + m9*x8   + m15*x9  + m21*x10 + m27*x11 + m33*x12;
163           x[9]  -= m4*x7  + m10*x8  + m16*x9  + m22*x10 + m28*x11 + m34*x12;
164           x[10] -= m5*x7  + m11*x8  + m17*x9  + m23*x10 + m29*x11 + m35*x12;
165           x[11] -= m6*x7  + m12*x8  + m18*x9  + m24*x10 + m30*x11 + m36*x12;
166 
167           x[12] -= m1*x13 + m7*x14  + m13*x15 + m19*x16 + m25*x17 + m31*x18;
168           x[13] -= m2*x13 + m8*x14  + m14*x15 + m20*x16 + m26*x17 + m32*x18;
169           x[14] -= m3*x13 + m9*x14  + m15*x15 + m21*x16 + m27*x17 + m33*x18;
170           x[15] -= m4*x13 + m10*x14 + m16*x15 + m22*x16 + m28*x17 + m34*x18;
171           x[16] -= m5*x13 + m11*x14 + m17*x15 + m23*x16 + m29*x17 + m35*x18;
172           x[17] -= m6*x13 + m12*x14 + m18*x15 + m24*x16 + m30*x17 + m36*x18;
173 
174           x[18] -= m1*x19 + m7*x20  + m13*x21 + m19*x22 + m25*x23 + m31*x24;
175           x[19] -= m2*x19 + m8*x20  + m14*x21 + m20*x22 + m26*x23 + m32*x24;
176           x[20] -= m3*x19 + m9*x20  + m15*x21 + m21*x22 + m27*x23 + m33*x24;
177           x[21] -= m4*x19 + m10*x20 + m16*x21 + m22*x22 + m28*x23 + m34*x24;
178           x[22] -= m5*x19 + m11*x20 + m17*x21 + m23*x22 + m29*x23 + m35*x24;
179           x[23] -= m6*x19 + m12*x20 + m18*x21 + m24*x22 + m30*x23 + m36*x24;
180 
181           x[24] -= m1*x25 + m7*x26  + m13*x27 + m19*x28 + m25*x29 + m31*x30;
182           x[25] -= m2*x25 + m8*x26  + m14*x27 + m20*x28 + m26*x29 + m32*x30;
183           x[26] -= m3*x25 + m9*x26  + m15*x27 + m21*x28 + m27*x29 + m33*x30;
184           x[27] -= m4*x25 + m10*x26 + m16*x27 + m22*x28 + m28*x29 + m34*x30;
185           x[28] -= m5*x25 + m11*x26 + m17*x27 + m23*x28 + m29*x29 + m35*x30;
186           x[29] -= m6*x25 + m12*x26 + m18*x27 + m24*x28 + m30*x29 + m36*x30;
187 
188           x[30] -= m1*x31 + m7*x32  + m13*x33 + m19*x34 + m25*x35 + m31*x36;
189           x[31] -= m2*x31 + m8*x32  + m14*x33 + m20*x34 + m26*x35 + m32*x36;
190           x[32] -= m3*x31 + m9*x32  + m15*x33 + m21*x34 + m27*x35 + m33*x36;
191           x[33] -= m4*x31 + m10*x32 + m16*x33 + m22*x34 + m28*x35 + m34*x36;
192           x[34] -= m5*x31 + m11*x32 + m17*x33 + m23*x34 + m29*x35 + m35*x36;
193           x[35] -= m6*x31 + m12*x32 + m18*x33 + m24*x34 + m30*x35 + m36*x36;
194 
195           pv += 36;
196         }
197         PetscCall(PetscLogFlops(432.0*nz+396.0));
198       }
199       row = *ajtmp++;
200     }
201     /* finished row so stick it into b->a */
202     pv = ba + 36*bi[i];
203     pj = bj + bi[i];
204     nz = bi[i+1] - bi[i];
205     for (j=0; j<nz; j++) {
206       x      = rtmp+36*pj[j];
207       pv[0]  = x[0];  pv[1]  = x[1];  pv[2]  = x[2];  pv[3]  = x[3];
208       pv[4]  = x[4];  pv[5]  = x[5];  pv[6]  = x[6];  pv[7]  = x[7];
209       pv[8]  = x[8];  pv[9]  = x[9];  pv[10] = x[10]; pv[11] = x[11];
210       pv[12] = x[12]; pv[13] = x[13]; pv[14] = x[14]; pv[15] = x[15];
211       pv[16] = x[16]; pv[17] = x[17]; pv[18] = x[18]; pv[19] = x[19];
212       pv[20] = x[20]; pv[21] = x[21]; pv[22] = x[22]; pv[23] = x[23];
213       pv[24] = x[24]; pv[25] = x[25]; pv[26] = x[26]; pv[27] = x[27];
214       pv[28] = x[28]; pv[29] = x[29]; pv[30] = x[30]; pv[31] = x[31];
215       pv[32] = x[32]; pv[33] = x[33]; pv[34] = x[34]; pv[35] = x[35];
216       pv    += 36;
217     }
218     /* invert diagonal block */
219     w    = ba + 36*diag_offset[i];
220     PetscCall(PetscKernel_A_gets_inverse_A_6(w,shift,allowzeropivot,&zeropivotdetected));
221     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
222   }
223 
224   PetscCall(PetscFree(rtmp));
225   PetscCall(ISRestoreIndices(isicol,&ic));
226   PetscCall(ISRestoreIndices(isrow,&r));
227 
228   C->ops->solve          = MatSolve_SeqBAIJ_6_inplace;
229   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_6_inplace;
230   C->assembled           = PETSC_TRUE;
231 
232   PetscCall(PetscLogFlops(1.333333333333*6*6*6*b->mbs)); /* from inverting diagonal blocks */
233   PetscFunctionReturn(0);
234 }
235 
236 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_6(Mat B,Mat A,const MatFactorInfo *info)
237 {
238   Mat            C     = B;
239   Mat_SeqBAIJ    *a    = (Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data;
240   IS             isrow = b->row,isicol = b->icol;
241   const PetscInt *r,*ic;
242   PetscInt       i,j,k,nz,nzL,row;
243   const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
244   const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag,*pj,bs2=a->bs2;
245   MatScalar      *rtmp,*pc,*mwork,*v,*pv,*aa=a->a;
246   PetscInt       flg;
247   PetscReal      shift = info->shiftamount;
248   PetscBool      allowzeropivot,zeropivotdetected;
249 
250   PetscFunctionBegin;
251   allowzeropivot = PetscNot(A->erroriffailure);
252   PetscCall(ISGetIndices(isrow,&r));
253   PetscCall(ISGetIndices(isicol,&ic));
254 
255   /* generate work space needed by the factorization */
256   PetscCall(PetscMalloc2(bs2*n,&rtmp,bs2,&mwork));
257   PetscCall(PetscArrayzero(rtmp,bs2*n));
258 
259   for (i=0; i<n; i++) {
260     /* zero rtmp */
261     /* L part */
262     nz    = bi[i+1] - bi[i];
263     bjtmp = bj + bi[i];
264     for  (j=0; j<nz; j++) {
265       PetscCall(PetscArrayzero(rtmp+bs2*bjtmp[j],bs2));
266     }
267 
268     /* U part */
269     nz    = bdiag[i] - bdiag[i+1];
270     bjtmp = bj + bdiag[i+1]+1;
271     for  (j=0; j<nz; j++) {
272       PetscCall(PetscArrayzero(rtmp+bs2*bjtmp[j],bs2));
273     }
274 
275     /* load in initial (unfactored row) */
276     nz    = ai[r[i]+1] - ai[r[i]];
277     ajtmp = aj + ai[r[i]];
278     v     = aa + bs2*ai[r[i]];
279     for (j=0; j<nz; j++) {
280       PetscCall(PetscArraycpy(rtmp+bs2*ic[ajtmp[j]],v+bs2*j,bs2));
281     }
282 
283     /* elimination */
284     bjtmp = bj + bi[i];
285     nzL   = bi[i+1] - bi[i];
286     for (k=0; k < nzL; k++) {
287       row = bjtmp[k];
288       pc  = rtmp + bs2*row;
289       for (flg=0,j=0; j<bs2; j++) {
290         if (pc[j]!=0.0) {
291           flg = 1;
292           break;
293         }
294       }
295       if (flg) {
296         pv = b->a + bs2*bdiag[row];
297         /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
298         PetscCall(PetscKernel_A_gets_A_times_B_6(pc,pv,mwork));
299 
300         pj = b->j + bdiag[row+1]+1; /* beginning of U(row,:) */
301         pv = b->a + bs2*(bdiag[row+1]+1);
302         nz = bdiag[row] - bdiag[row+1] -  1; /* num of entries inU(row,:), excluding diag */
303         for (j=0; j<nz; j++) {
304           /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
305           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
306           v    = rtmp + bs2*pj[j];
307           PetscCall(PetscKernel_A_gets_A_minus_B_times_C_6(v,pc,pv));
308           pv  += bs2;
309         }
310         PetscCall(PetscLogFlops(432.0*nz+396)); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
311       }
312     }
313 
314     /* finished row so stick it into b->a */
315     /* L part */
316     pv = b->a + bs2*bi[i];
317     pj = b->j + bi[i];
318     nz = bi[i+1] - bi[i];
319     for (j=0; j<nz; j++) {
320       PetscCall(PetscArraycpy(pv+bs2*j,rtmp+bs2*pj[j],bs2));
321     }
322 
323     /* Mark diagonal and invert diagonal for simpler triangular solves */
324     pv   = b->a + bs2*bdiag[i];
325     pj   = b->j + bdiag[i];
326     PetscCall(PetscArraycpy(pv,rtmp+bs2*pj[0],bs2));
327     PetscCall(PetscKernel_A_gets_inverse_A_6(pv,shift,allowzeropivot,&zeropivotdetected));
328     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
329 
330     /* U part */
331     pv = b->a + bs2*(bdiag[i+1]+1);
332     pj = b->j + bdiag[i+1]+1;
333     nz = bdiag[i] - bdiag[i+1] - 1;
334     for (j=0; j<nz; j++) {
335       PetscCall(PetscArraycpy(pv+bs2*j,rtmp+bs2*pj[j],bs2));
336     }
337   }
338 
339   PetscCall(PetscFree2(rtmp,mwork));
340   PetscCall(ISRestoreIndices(isicol,&ic));
341   PetscCall(ISRestoreIndices(isrow,&r));
342 
343   C->ops->solve          = MatSolve_SeqBAIJ_6;
344   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_6;
345   C->assembled           = PETSC_TRUE;
346 
347   PetscCall(PetscLogFlops(1.333333333333*6*6*6*n)); /* from inverting diagonal blocks */
348   PetscFunctionReturn(0);
349 }
350 
351 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering_inplace(Mat C,Mat A,const MatFactorInfo *info)
352 {
353   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)C->data;
354   PetscInt       i,j,n = a->mbs,*bi = b->i,*bj = b->j;
355   PetscInt       *ajtmpold,*ajtmp,nz,row;
356   PetscInt       *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj;
357   MatScalar      *pv,*v,*rtmp,*pc,*w,*x;
358   MatScalar      x1,x2,x3,x4,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15;
359   MatScalar      x16,x17,x18,x19,x20,x21,x22,x23,x24,x25;
360   MatScalar      p1,p2,p3,p4,p5,p6,p7,p8,p9,p10,p11,p12,p13,p14,p15;
361   MatScalar      p16,p17,p18,p19,p20,p21,p22,p23,p24,p25;
362   MatScalar      m1,m2,m3,m4,m5,m6,m7,m8,m9,m10,m11,m12,m13,m14,m15;
363   MatScalar      m16,m17,m18,m19,m20,m21,m22,m23,m24,m25;
364   MatScalar      p26,p27,p28,p29,p30,p31,p32,p33,p34,p35,p36;
365   MatScalar      x26,x27,x28,x29,x30,x31,x32,x33,x34,x35,x36;
366   MatScalar      m26,m27,m28,m29,m30,m31,m32,m33,m34,m35,m36;
367   MatScalar      *ba   = b->a,*aa = a->a;
368   PetscReal      shift = info->shiftamount;
369   PetscBool      allowzeropivot,zeropivotdetected;
370 
371   PetscFunctionBegin;
372   allowzeropivot = PetscNot(A->erroriffailure);
373   PetscCall(PetscMalloc1(36*(n+1),&rtmp));
374   for (i=0; i<n; i++) {
375     nz    = bi[i+1] - bi[i];
376     ajtmp = bj + bi[i];
377     for  (j=0; j<nz; j++) {
378       x     = rtmp+36*ajtmp[j];
379       x[0]  = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = x[9] = 0.0;
380       x[10] = x[11] = x[12] = x[13] = x[14] = x[15] = x[16] = x[17] = 0.0;
381       x[18] = x[19] = x[20] = x[21] = x[22] = x[23] = x[24] = x[25] = 0.0;
382       x[26] = x[27] = x[28] = x[29] = x[30] = x[31] = x[32] = x[33] = 0.0;
383       x[34] = x[35] = 0.0;
384     }
385     /* load in initial (unfactored row) */
386     nz       = ai[i+1] - ai[i];
387     ajtmpold = aj + ai[i];
388     v        = aa + 36*ai[i];
389     for (j=0; j<nz; j++) {
390       x     = rtmp+36*ajtmpold[j];
391       x[0]  =  v[0];  x[1] =  v[1];  x[2] =  v[2];  x[3] =  v[3];
392       x[4]  =  v[4];  x[5] =  v[5];  x[6] =  v[6];  x[7] =  v[7];
393       x[8]  =  v[8];  x[9] =  v[9];  x[10] = v[10]; x[11] = v[11];
394       x[12] = v[12]; x[13] = v[13]; x[14] = v[14]; x[15] = v[15];
395       x[16] = v[16]; x[17] = v[17]; x[18] = v[18]; x[19] = v[19];
396       x[20] = v[20]; x[21] = v[21]; x[22] = v[22]; x[23] = v[23];
397       x[24] = v[24]; x[25] = v[25]; x[26] = v[26]; x[27] = v[27];
398       x[28] = v[28]; x[29] = v[29]; x[30] = v[30]; x[31] = v[31];
399       x[32] = v[32]; x[33] = v[33]; x[34] = v[34]; x[35] = v[35];
400       v    += 36;
401     }
402     row = *ajtmp++;
403     while (row < i) {
404       pc  = rtmp + 36*row;
405       p1  = pc[0];  p2  = pc[1];  p3  = pc[2];  p4  = pc[3];
406       p5  = pc[4];  p6  = pc[5];  p7  = pc[6];  p8  = pc[7];
407       p9  = pc[8];  p10 = pc[9];  p11 = pc[10]; p12 = pc[11];
408       p13 = pc[12]; p14 = pc[13]; p15 = pc[14]; p16 = pc[15];
409       p17 = pc[16]; p18 = pc[17]; p19 = pc[18]; p20 = pc[19];
410       p21 = pc[20]; p22 = pc[21]; p23 = pc[22]; p24 = pc[23];
411       p25 = pc[24]; p26 = pc[25]; p27 = pc[26]; p28 = pc[27];
412       p29 = pc[28]; p30 = pc[29]; p31 = pc[30]; p32 = pc[31];
413       p33 = pc[32]; p34 = pc[33]; p35 = pc[34]; p36 = pc[35];
414       if (p1  != 0.0 || p2  != 0.0 || p3  != 0.0 || p4  != 0.0 ||
415           p5  != 0.0 || p6  != 0.0 || p7  != 0.0 || p8  != 0.0 ||
416           p9  != 0.0 || p10 != 0.0 || p11 != 0.0 || p12 != 0.0 ||
417           p13 != 0.0 || p14 != 0.0 || p15 != 0.0 || p16 != 0.0 ||
418           p17 != 0.0 || p18 != 0.0 || p19 != 0.0 || p20 != 0.0 ||
419           p21 != 0.0 || p22 != 0.0 || p23 != 0.0 || p24 != 0.0 ||
420           p25 != 0.0 || p26 != 0.0 || p27 != 0.0 || p28 != 0.0 ||
421           p29 != 0.0 || p30 != 0.0 || p31 != 0.0 || p32 != 0.0 ||
422           p33 != 0.0 || p34 != 0.0 || p35 != 0.0 || p36 != 0.0) {
423         pv    = ba + 36*diag_offset[row];
424         pj    = bj + diag_offset[row] + 1;
425         x1    = pv[0];  x2  = pv[1];  x3  = pv[2];  x4  = pv[3];
426         x5    = pv[4];  x6  = pv[5];  x7  = pv[6];  x8  = pv[7];
427         x9    = pv[8];  x10 = pv[9];  x11 = pv[10]; x12 = pv[11];
428         x13   = pv[12]; x14 = pv[13]; x15 = pv[14]; x16 = pv[15];
429         x17   = pv[16]; x18 = pv[17]; x19 = pv[18]; x20 = pv[19];
430         x21   = pv[20]; x22 = pv[21]; x23 = pv[22]; x24 = pv[23];
431         x25   = pv[24]; x26 = pv[25]; x27 = pv[26]; x28 = pv[27];
432         x29   = pv[28]; x30 = pv[29]; x31 = pv[30]; x32 = pv[31];
433         x33   = pv[32]; x34 = pv[33]; x35 = pv[34]; x36 = pv[35];
434         pc[0] = m1  = p1*x1  + p7*x2   + p13*x3  + p19*x4  + p25*x5  + p31*x6;
435         pc[1] = m2  = p2*x1  + p8*x2   + p14*x3  + p20*x4  + p26*x5  + p32*x6;
436         pc[2] = m3  = p3*x1  + p9*x2   + p15*x3  + p21*x4  + p27*x5  + p33*x6;
437         pc[3] = m4  = p4*x1  + p10*x2  + p16*x3  + p22*x4  + p28*x5  + p34*x6;
438         pc[4] = m5  = p5*x1  + p11*x2  + p17*x3  + p23*x4  + p29*x5  + p35*x6;
439         pc[5] = m6  = p6*x1  + p12*x2  + p18*x3  + p24*x4  + p30*x5  + p36*x6;
440 
441         pc[6]  = m7  = p1*x7  + p7*x8   + p13*x9  + p19*x10 + p25*x11 + p31*x12;
442         pc[7]  = m8  = p2*x7  + p8*x8   + p14*x9  + p20*x10 + p26*x11 + p32*x12;
443         pc[8]  = m9  = p3*x7  + p9*x8   + p15*x9  + p21*x10 + p27*x11 + p33*x12;
444         pc[9]  = m10 = p4*x7  + p10*x8  + p16*x9  + p22*x10 + p28*x11 + p34*x12;
445         pc[10] = m11 = p5*x7  + p11*x8  + p17*x9  + p23*x10 + p29*x11 + p35*x12;
446         pc[11] = m12 = p6*x7  + p12*x8  + p18*x9  + p24*x10 + p30*x11 + p36*x12;
447 
448         pc[12] = m13 = p1*x13 + p7*x14  + p13*x15 + p19*x16 + p25*x17 + p31*x18;
449         pc[13] = m14 = p2*x13 + p8*x14  + p14*x15 + p20*x16 + p26*x17 + p32*x18;
450         pc[14] = m15 = p3*x13 + p9*x14  + p15*x15 + p21*x16 + p27*x17 + p33*x18;
451         pc[15] = m16 = p4*x13 + p10*x14 + p16*x15 + p22*x16 + p28*x17 + p34*x18;
452         pc[16] = m17 = p5*x13 + p11*x14 + p17*x15 + p23*x16 + p29*x17 + p35*x18;
453         pc[17] = m18 = p6*x13 + p12*x14 + p18*x15 + p24*x16 + p30*x17 + p36*x18;
454 
455         pc[18] = m19 = p1*x19 + p7*x20  + p13*x21 + p19*x22 + p25*x23 + p31*x24;
456         pc[19] = m20 = p2*x19 + p8*x20  + p14*x21 + p20*x22 + p26*x23 + p32*x24;
457         pc[20] = m21 = p3*x19 + p9*x20  + p15*x21 + p21*x22 + p27*x23 + p33*x24;
458         pc[21] = m22 = p4*x19 + p10*x20 + p16*x21 + p22*x22 + p28*x23 + p34*x24;
459         pc[22] = m23 = p5*x19 + p11*x20 + p17*x21 + p23*x22 + p29*x23 + p35*x24;
460         pc[23] = m24 = p6*x19 + p12*x20 + p18*x21 + p24*x22 + p30*x23 + p36*x24;
461 
462         pc[24] = m25 = p1*x25 + p7*x26  + p13*x27 + p19*x28 + p25*x29 + p31*x30;
463         pc[25] = m26 = p2*x25 + p8*x26  + p14*x27 + p20*x28 + p26*x29 + p32*x30;
464         pc[26] = m27 = p3*x25 + p9*x26  + p15*x27 + p21*x28 + p27*x29 + p33*x30;
465         pc[27] = m28 = p4*x25 + p10*x26 + p16*x27 + p22*x28 + p28*x29 + p34*x30;
466         pc[28] = m29 = p5*x25 + p11*x26 + p17*x27 + p23*x28 + p29*x29 + p35*x30;
467         pc[29] = m30 = p6*x25 + p12*x26 + p18*x27 + p24*x28 + p30*x29 + p36*x30;
468 
469         pc[30] = m31 = p1*x31 + p7*x32  + p13*x33 + p19*x34 + p25*x35 + p31*x36;
470         pc[31] = m32 = p2*x31 + p8*x32  + p14*x33 + p20*x34 + p26*x35 + p32*x36;
471         pc[32] = m33 = p3*x31 + p9*x32  + p15*x33 + p21*x34 + p27*x35 + p33*x36;
472         pc[33] = m34 = p4*x31 + p10*x32 + p16*x33 + p22*x34 + p28*x35 + p34*x36;
473         pc[34] = m35 = p5*x31 + p11*x32 + p17*x33 + p23*x34 + p29*x35 + p35*x36;
474         pc[35] = m36 = p6*x31 + p12*x32 + p18*x33 + p24*x34 + p30*x35 + p36*x36;
475 
476         nz  = bi[row+1] - diag_offset[row] - 1;
477         pv += 36;
478         for (j=0; j<nz; j++) {
479           x1    = pv[0];  x2  = pv[1];  x3  = pv[2];  x4  = pv[3];
480           x5    = pv[4];  x6  = pv[5];  x7  = pv[6];  x8  = pv[7];
481           x9    = pv[8];  x10 = pv[9];  x11 = pv[10]; x12 = pv[11];
482           x13   = pv[12]; x14 = pv[13]; x15 = pv[14]; x16 = pv[15];
483           x17   = pv[16]; x18 = pv[17]; x19 = pv[18]; x20 = pv[19];
484           x21   = pv[20]; x22 = pv[21]; x23 = pv[22]; x24 = pv[23];
485           x25   = pv[24]; x26 = pv[25]; x27 = pv[26]; x28 = pv[27];
486           x29   = pv[28]; x30 = pv[29]; x31 = pv[30]; x32 = pv[31];
487           x33   = pv[32]; x34 = pv[33]; x35 = pv[34]; x36 = pv[35];
488           x     = rtmp + 36*pj[j];
489           x[0] -= m1*x1  + m7*x2   + m13*x3  + m19*x4  + m25*x5  + m31*x6;
490           x[1] -= m2*x1  + m8*x2   + m14*x3  + m20*x4  + m26*x5  + m32*x6;
491           x[2] -= m3*x1  + m9*x2   + m15*x3  + m21*x4  + m27*x5  + m33*x6;
492           x[3] -= m4*x1  + m10*x2  + m16*x3  + m22*x4  + m28*x5  + m34*x6;
493           x[4] -= m5*x1  + m11*x2  + m17*x3  + m23*x4  + m29*x5  + m35*x6;
494           x[5] -= m6*x1  + m12*x2  + m18*x3  + m24*x4  + m30*x5  + m36*x6;
495 
496           x[6]  -= m1*x7  + m7*x8   + m13*x9  + m19*x10 + m25*x11 + m31*x12;
497           x[7]  -= m2*x7  + m8*x8   + m14*x9  + m20*x10 + m26*x11 + m32*x12;
498           x[8]  -= m3*x7  + m9*x8   + m15*x9  + m21*x10 + m27*x11 + m33*x12;
499           x[9]  -= m4*x7  + m10*x8  + m16*x9  + m22*x10 + m28*x11 + m34*x12;
500           x[10] -= m5*x7  + m11*x8  + m17*x9  + m23*x10 + m29*x11 + m35*x12;
501           x[11] -= m6*x7  + m12*x8  + m18*x9  + m24*x10 + m30*x11 + m36*x12;
502 
503           x[12] -= m1*x13 + m7*x14  + m13*x15 + m19*x16 + m25*x17 + m31*x18;
504           x[13] -= m2*x13 + m8*x14  + m14*x15 + m20*x16 + m26*x17 + m32*x18;
505           x[14] -= m3*x13 + m9*x14  + m15*x15 + m21*x16 + m27*x17 + m33*x18;
506           x[15] -= m4*x13 + m10*x14 + m16*x15 + m22*x16 + m28*x17 + m34*x18;
507           x[16] -= m5*x13 + m11*x14 + m17*x15 + m23*x16 + m29*x17 + m35*x18;
508           x[17] -= m6*x13 + m12*x14 + m18*x15 + m24*x16 + m30*x17 + m36*x18;
509 
510           x[18] -= m1*x19 + m7*x20  + m13*x21 + m19*x22 + m25*x23 + m31*x24;
511           x[19] -= m2*x19 + m8*x20  + m14*x21 + m20*x22 + m26*x23 + m32*x24;
512           x[20] -= m3*x19 + m9*x20  + m15*x21 + m21*x22 + m27*x23 + m33*x24;
513           x[21] -= m4*x19 + m10*x20 + m16*x21 + m22*x22 + m28*x23 + m34*x24;
514           x[22] -= m5*x19 + m11*x20 + m17*x21 + m23*x22 + m29*x23 + m35*x24;
515           x[23] -= m6*x19 + m12*x20 + m18*x21 + m24*x22 + m30*x23 + m36*x24;
516 
517           x[24] -= m1*x25 + m7*x26  + m13*x27 + m19*x28 + m25*x29 + m31*x30;
518           x[25] -= m2*x25 + m8*x26  + m14*x27 + m20*x28 + m26*x29 + m32*x30;
519           x[26] -= m3*x25 + m9*x26  + m15*x27 + m21*x28 + m27*x29 + m33*x30;
520           x[27] -= m4*x25 + m10*x26 + m16*x27 + m22*x28 + m28*x29 + m34*x30;
521           x[28] -= m5*x25 + m11*x26 + m17*x27 + m23*x28 + m29*x29 + m35*x30;
522           x[29] -= m6*x25 + m12*x26 + m18*x27 + m24*x28 + m30*x29 + m36*x30;
523 
524           x[30] -= m1*x31 + m7*x32  + m13*x33 + m19*x34 + m25*x35 + m31*x36;
525           x[31] -= m2*x31 + m8*x32  + m14*x33 + m20*x34 + m26*x35 + m32*x36;
526           x[32] -= m3*x31 + m9*x32  + m15*x33 + m21*x34 + m27*x35 + m33*x36;
527           x[33] -= m4*x31 + m10*x32 + m16*x33 + m22*x34 + m28*x35 + m34*x36;
528           x[34] -= m5*x31 + m11*x32 + m17*x33 + m23*x34 + m29*x35 + m35*x36;
529           x[35] -= m6*x31 + m12*x32 + m18*x33 + m24*x34 + m30*x35 + m36*x36;
530 
531           pv += 36;
532         }
533         PetscCall(PetscLogFlops(432.0*nz+396.0));
534       }
535       row = *ajtmp++;
536     }
537     /* finished row so stick it into b->a */
538     pv = ba + 36*bi[i];
539     pj = bj + bi[i];
540     nz = bi[i+1] - bi[i];
541     for (j=0; j<nz; j++) {
542       x      = rtmp+36*pj[j];
543       pv[0]  = x[0];  pv[1]  = x[1];  pv[2]  = x[2];  pv[3]  = x[3];
544       pv[4]  = x[4];  pv[5]  = x[5];  pv[6]  = x[6];  pv[7]  = x[7];
545       pv[8]  = x[8];  pv[9]  = x[9];  pv[10] = x[10]; pv[11] = x[11];
546       pv[12] = x[12]; pv[13] = x[13]; pv[14] = x[14]; pv[15] = x[15];
547       pv[16] = x[16]; pv[17] = x[17]; pv[18] = x[18]; pv[19] = x[19];
548       pv[20] = x[20]; pv[21] = x[21]; pv[22] = x[22]; pv[23] = x[23];
549       pv[24] = x[24]; pv[25] = x[25]; pv[26] = x[26]; pv[27] = x[27];
550       pv[28] = x[28]; pv[29] = x[29]; pv[30] = x[30]; pv[31] = x[31];
551       pv[32] = x[32]; pv[33] = x[33]; pv[34] = x[34]; pv[35] = x[35];
552       pv    += 36;
553     }
554     /* invert diagonal block */
555     w    = ba + 36*diag_offset[i];
556     PetscCall(PetscKernel_A_gets_inverse_A_6(w,shift,allowzeropivot,&zeropivotdetected));
557     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
558   }
559 
560   PetscCall(PetscFree(rtmp));
561 
562   C->ops->solve          = MatSolve_SeqBAIJ_6_NaturalOrdering_inplace;
563   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_6_NaturalOrdering_inplace;
564   C->assembled           = PETSC_TRUE;
565 
566   PetscCall(PetscLogFlops(1.333333333333*6*6*6*b->mbs)); /* from inverting diagonal blocks */
567   PetscFunctionReturn(0);
568 }
569 
570 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info)
571 {
572   Mat            C =B;
573   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data;
574   PetscInt       i,j,k,nz,nzL,row;
575   const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
576   const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag,*pj,bs2=a->bs2;
577   MatScalar      *rtmp,*pc,*mwork,*v,*pv,*aa=a->a;
578   PetscInt       flg;
579   PetscReal      shift = info->shiftamount;
580   PetscBool      allowzeropivot,zeropivotdetected;
581 
582   PetscFunctionBegin;
583   allowzeropivot = PetscNot(A->erroriffailure);
584 
585   /* generate work space needed by the factorization */
586   PetscCall(PetscMalloc2(bs2*n,&rtmp,bs2,&mwork));
587   PetscCall(PetscArrayzero(rtmp,bs2*n));
588 
589   for (i=0; i<n; i++) {
590     /* zero rtmp */
591     /* L part */
592     nz    = bi[i+1] - bi[i];
593     bjtmp = bj + bi[i];
594     for  (j=0; j<nz; j++) {
595       PetscCall(PetscArrayzero(rtmp+bs2*bjtmp[j],bs2));
596     }
597 
598     /* U part */
599     nz    = bdiag[i] - bdiag[i+1];
600     bjtmp = bj + bdiag[i+1]+1;
601     for  (j=0; j<nz; j++) {
602       PetscCall(PetscArrayzero(rtmp+bs2*bjtmp[j],bs2));
603     }
604 
605     /* load in initial (unfactored row) */
606     nz    = ai[i+1] - ai[i];
607     ajtmp = aj + ai[i];
608     v     = aa + bs2*ai[i];
609     for (j=0; j<nz; j++) {
610       PetscCall(PetscArraycpy(rtmp+bs2*ajtmp[j],v+bs2*j,bs2));
611     }
612 
613     /* elimination */
614     bjtmp = bj + bi[i];
615     nzL   = bi[i+1] - bi[i];
616     for (k=0; k < nzL; k++) {
617       row = bjtmp[k];
618       pc  = rtmp + bs2*row;
619       for (flg=0,j=0; j<bs2; j++) {
620         if (pc[j]!=0.0) {
621           flg = 1;
622           break;
623         }
624       }
625       if (flg) {
626         pv = b->a + bs2*bdiag[row];
627         /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
628         PetscCall(PetscKernel_A_gets_A_times_B_6(pc,pv,mwork));
629 
630         pj = b->j + bdiag[row+1]+1; /* beginning of U(row,:) */
631         pv = b->a + bs2*(bdiag[row+1]+1);
632         nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries inU(row,:), excluding diag */
633         for (j=0; j<nz; j++) {
634           /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
635           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
636           v    = rtmp + bs2*pj[j];
637           PetscCall(PetscKernel_A_gets_A_minus_B_times_C_6(v,pc,pv));
638           pv  += bs2;
639         }
640         PetscCall(PetscLogFlops(432.0*nz+396)); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
641       }
642     }
643 
644     /* finished row so stick it into b->a */
645     /* L part */
646     pv = b->a + bs2*bi[i];
647     pj = b->j + bi[i];
648     nz = bi[i+1] - bi[i];
649     for (j=0; j<nz; j++) {
650       PetscCall(PetscArraycpy(pv+bs2*j,rtmp+bs2*pj[j],bs2));
651     }
652 
653     /* Mark diagonal and invert diagonal for simpler triangular solves */
654     pv   = b->a + bs2*bdiag[i];
655     pj   = b->j + bdiag[i];
656     PetscCall(PetscArraycpy(pv,rtmp+bs2*pj[0],bs2));
657     PetscCall(PetscKernel_A_gets_inverse_A_6(pv,shift,allowzeropivot,&zeropivotdetected));
658     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
659 
660     /* U part */
661     pv = b->a + bs2*(bdiag[i+1]+1);
662     pj = b->j + bdiag[i+1]+1;
663     nz = bdiag[i] - bdiag[i+1] - 1;
664     for (j=0; j<nz; j++) {
665       PetscCall(PetscArraycpy(pv+bs2*j,rtmp+bs2*pj[j],bs2));
666     }
667   }
668   PetscCall(PetscFree2(rtmp,mwork));
669 
670   C->ops->solve          = MatSolve_SeqBAIJ_6_NaturalOrdering;
671   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_6_NaturalOrdering;
672   C->assembled           = PETSC_TRUE;
673 
674   PetscCall(PetscLogFlops(1.333333333333*6*6*6*n)); /* from inverting diagonal blocks */
675   PetscFunctionReturn(0);
676 }
677