xref: /petsc/src/mat/impls/baij/seq/baijfact13.c (revision 922032d7fd3c7caf4f4a8a5beb94a62d36a9c4f5)
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       Version for when blocks are 3 by 3
10 */
11 #undef __FUNCT__
12 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_3_inplace"
13 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3_inplace(Mat C,Mat A,const MatFactorInfo *info)
14 {
15   Mat_SeqBAIJ    *a    = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)C->data;
16   IS             isrow = b->row,isicol = b->icol;
17   PetscErrorCode ierr;
18   const PetscInt *r,*ic;
19   PetscInt       i,j,n = a->mbs,*bi = b->i,*bj = b->j;
20   PetscInt       *ajtmpold,*ajtmp,nz,row,*ai=a->i,*aj=a->j;
21   PetscInt       *diag_offset = b->diag,idx,*pj;
22   MatScalar      *pv,*v,*rtmp,*pc,*w,*x;
23   MatScalar      p1,p2,p3,p4,m1,m2,m3,m4,m5,m6,m7,m8,m9,x1,x2,x3,x4;
24   MatScalar      p5,p6,p7,p8,p9,x5,x6,x7,x8,x9;
25   MatScalar      *ba   = b->a,*aa = a->a;
26   PetscReal      shift = info->shiftamount;
27   PetscBool      zeropivotdetected;
28 
29   PetscFunctionBegin;
30   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
31   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
32   ierr = PetscMalloc1(9*(n+1),&rtmp);CHKERRQ(ierr);
33 
34   for (i=0; i<n; i++) {
35     nz    = bi[i+1] - bi[i];
36     ajtmp = bj + bi[i];
37     for  (j=0; j<nz; j++) {
38       x    = rtmp + 9*ajtmp[j];
39       x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = 0.0;
40     }
41     /* load in initial (unfactored row) */
42     idx      = r[i];
43     nz       = ai[idx+1] - ai[idx];
44     ajtmpold = aj + ai[idx];
45     v        = aa + 9*ai[idx];
46     for (j=0; j<nz; j++) {
47       x    = rtmp + 9*ic[ajtmpold[j]];
48       x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3];
49       x[4] = v[4]; x[5] = v[5]; x[6] = v[6]; x[7] = v[7]; x[8] = v[8];
50       v   += 9;
51     }
52     row = *ajtmp++;
53     while (row < i) {
54       pc = rtmp + 9*row;
55       p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3];
56       p5 = pc[4]; p6 = pc[5]; p7 = pc[6]; p8 = pc[7]; p9 = pc[8];
57       if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 || p5 != 0.0 ||
58           p6 != 0.0 || p7 != 0.0 || p8 != 0.0 || p9 != 0.0) {
59         pv    = ba + 9*diag_offset[row];
60         pj    = bj + diag_offset[row] + 1;
61         x1    = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
62         x5    = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8];
63         pc[0] = m1 = p1*x1 + p4*x2 + p7*x3;
64         pc[1] = m2 = p2*x1 + p5*x2 + p8*x3;
65         pc[2] = m3 = p3*x1 + p6*x2 + p9*x3;
66 
67         pc[3] = m4 = p1*x4 + p4*x5 + p7*x6;
68         pc[4] = m5 = p2*x4 + p5*x5 + p8*x6;
69         pc[5] = m6 = p3*x4 + p6*x5 + p9*x6;
70 
71         pc[6] = m7 = p1*x7 + p4*x8 + p7*x9;
72         pc[7] = m8 = p2*x7 + p5*x8 + p8*x9;
73         pc[8] = m9 = p3*x7 + p6*x8 + p9*x9;
74         nz    = bi[row+1] - diag_offset[row] - 1;
75         pv   += 9;
76         for (j=0; j<nz; j++) {
77           x1    = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
78           x5    = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8];
79           x     = rtmp + 9*pj[j];
80           x[0] -= m1*x1 + m4*x2 + m7*x3;
81           x[1] -= m2*x1 + m5*x2 + m8*x3;
82           x[2] -= m3*x1 + m6*x2 + m9*x3;
83 
84           x[3] -= m1*x4 + m4*x5 + m7*x6;
85           x[4] -= m2*x4 + m5*x5 + m8*x6;
86           x[5] -= m3*x4 + m6*x5 + m9*x6;
87 
88           x[6] -= m1*x7 + m4*x8 + m7*x9;
89           x[7] -= m2*x7 + m5*x8 + m8*x9;
90           x[8] -= m3*x7 + m6*x8 + m9*x9;
91           pv   += 9;
92         }
93         ierr = PetscLogFlops(54.0*nz+36.0);CHKERRQ(ierr);
94       }
95       row = *ajtmp++;
96     }
97     /* finished row so stick it into b->a */
98     pv = ba + 9*bi[i];
99     pj = bj + bi[i];
100     nz = bi[i+1] - bi[i];
101     for (j=0; j<nz; j++) {
102       x     = rtmp + 9*pj[j];
103       pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3];
104       pv[4] = x[4]; pv[5] = x[5]; pv[6] = x[6]; pv[7] = x[7]; pv[8] = x[8];
105       pv   += 9;
106     }
107     /* invert diagonal block */
108     w    = ba + 9*diag_offset[i];
109     ierr = PetscKernel_A_gets_inverse_A_3(w,shift,!A->erroriffailure,&zeropivotdetected);CHKERRQ(ierr);
110     if (zeropivotdetected) C->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
111   }
112 
113   ierr = PetscFree(rtmp);CHKERRQ(ierr);
114   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
115   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
116 
117   C->ops->solve          = MatSolve_SeqBAIJ_3_inplace;
118   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_3_inplace;
119   C->assembled           = PETSC_TRUE;
120 
121   ierr = PetscLogFlops(1.333333333333*3*3*3*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
122   PetscFunctionReturn(0);
123 }
124 
125 /* MatLUFactorNumeric_SeqBAIJ_3 -
126      copied from MatLUFactorNumeric_SeqBAIJ_N_inplace() and manually re-implemented
127        PetscKernel_A_gets_A_times_B()
128        PetscKernel_A_gets_A_minus_B_times_C()
129        PetscKernel_A_gets_inverse_A()
130 */
131 #undef __FUNCT__
132 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_3"
133 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3(Mat B,Mat A,const MatFactorInfo *info)
134 {
135   Mat            C     =B;
136   Mat_SeqBAIJ    *a    =(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data;
137   IS             isrow = b->row,isicol = b->icol;
138   PetscErrorCode ierr;
139   const PetscInt *r,*ic;
140   PetscInt       i,j,k,nz,nzL,row;
141   const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
142   const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag,*pj,bs2=a->bs2;
143   MatScalar      *rtmp,*pc,*mwork,*v,*pv,*aa=a->a;
144   PetscInt       flg;
145   PetscReal      shift = info->shiftamount;
146   PetscBool      zeropivotdetected;
147 
148   PetscFunctionBegin;
149   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
150   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
151 
152   /* generate work space needed by the factorization */
153   ierr = PetscMalloc2(bs2*n,&rtmp,bs2,&mwork);CHKERRQ(ierr);
154   ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr);
155 
156   for (i=0; i<n; i++) {
157     /* zero rtmp */
158     /* L part */
159     nz    = bi[i+1] - bi[i];
160     bjtmp = bj + bi[i];
161     for  (j=0; j<nz; j++) {
162       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
163     }
164 
165     /* U part */
166     nz    = bdiag[i] - bdiag[i+1];
167     bjtmp = bj + bdiag[i+1]+1;
168     for  (j=0; j<nz; j++) {
169       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
170     }
171 
172     /* load in initial (unfactored row) */
173     nz    = ai[r[i]+1] - ai[r[i]];
174     ajtmp = aj + ai[r[i]];
175     v     = aa + bs2*ai[r[i]];
176     for (j=0; j<nz; j++) {
177       ierr = PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr);
178     }
179 
180     /* elimination */
181     bjtmp = bj + bi[i];
182     nzL   = bi[i+1] - bi[i];
183     for (k = 0; k < nzL; k++) {
184       row = bjtmp[k];
185       pc  = rtmp + bs2*row;
186       for (flg=0,j=0; j<bs2; j++) {
187         if (pc[j]!=0.0) {
188           flg = 1;
189           break;
190         }
191       }
192       if (flg) {
193         pv = b->a + bs2*bdiag[row];
194         /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
195         ierr = PetscKernel_A_gets_A_times_B_3(pc,pv,mwork);CHKERRQ(ierr);
196 
197         pj = b->j + bdiag[row+1] + 1; /* beginning of U(row,:) */
198         pv = b->a + bs2*(bdiag[row+1]+1);
199         nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:) excluding diag */
200         for (j=0; j<nz; j++) {
201           /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
202           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
203           v    = rtmp + bs2*pj[j];
204           ierr = PetscKernel_A_gets_A_minus_B_times_C_3(v,pc,pv);CHKERRQ(ierr);
205           pv  += bs2;
206         }
207         ierr = PetscLogFlops(54*nz+45);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
208       }
209     }
210 
211     /* finished row so stick it into b->a */
212     /* L part */
213     pv = b->a + bs2*bi[i];
214     pj = b->j + bi[i];
215     nz = bi[i+1] - bi[i];
216     for (j=0; j<nz; j++) {
217       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
218     }
219 
220     /* Mark diagonal and invert diagonal for simplier triangular solves */
221     pv   = b->a + bs2*bdiag[i];
222     pj   = b->j + bdiag[i];
223     ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr);
224     /* ierr = PetscKernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); */
225     ierr = PetscKernel_A_gets_inverse_A_3(pv,shift,!A->erroriffailure,&zeropivotdetected);CHKERRQ(ierr);
226     if (zeropivotdetected) B->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
227 
228     /* U part */
229     pj = b->j + bdiag[i+1] + 1;
230     pv = b->a + bs2*(bdiag[i+1]+1);
231     nz = bdiag[i] - bdiag[i+1] - 1;
232     for (j=0; j<nz; j++) {
233       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
234     }
235   }
236 
237   ierr = PetscFree2(rtmp,mwork);CHKERRQ(ierr);
238   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
239   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
240 
241   C->ops->solve          = MatSolve_SeqBAIJ_3;
242   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_3;
243   C->assembled           = PETSC_TRUE;
244 
245   ierr = PetscLogFlops(1.333333333333*3*3*3*n);CHKERRQ(ierr); /* from inverting diagonal blocks */
246   PetscFunctionReturn(0);
247 }
248 
249 #undef __FUNCT__
250 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering_inplace"
251 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering_inplace(Mat C,Mat A,const MatFactorInfo *info)
252 {
253   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)C->data;
254   PetscErrorCode ierr;
255   PetscInt       i,j,n = a->mbs,*bi = b->i,*bj = b->j;
256   PetscInt       *ajtmpold,*ajtmp,nz,row;
257   PetscInt       *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj;
258   MatScalar      *pv,*v,*rtmp,*pc,*w,*x;
259   MatScalar      p1,p2,p3,p4,m1,m2,m3,m4,m5,m6,m7,m8,m9,x1,x2,x3,x4;
260   MatScalar      p5,p6,p7,p8,p9,x5,x6,x7,x8,x9;
261   MatScalar      *ba   = b->a,*aa = a->a;
262   PetscReal      shift = info->shiftamount;
263   PetscBool      zeropivotdetected;
264 
265   PetscFunctionBegin;
266   ierr = PetscMalloc1(9*(n+1),&rtmp);CHKERRQ(ierr);
267 
268   for (i=0; i<n; i++) {
269     nz    = bi[i+1] - bi[i];
270     ajtmp = bj + bi[i];
271     for  (j=0; j<nz; j++) {
272       x    = rtmp+9*ajtmp[j];
273       x[0] = x[1]  = x[2]  = x[3]  = x[4]  = x[5]  = x[6] = x[7] = x[8] = 0.0;
274     }
275     /* load in initial (unfactored row) */
276     nz       = ai[i+1] - ai[i];
277     ajtmpold = aj + ai[i];
278     v        = aa + 9*ai[i];
279     for (j=0; j<nz; j++) {
280       x    = rtmp+9*ajtmpold[j];
281       x[0] = v[0];  x[1]  = v[1];  x[2]  = v[2];  x[3]  = v[3];
282       x[4] = v[4];  x[5]  = v[5];  x[6]  = v[6];  x[7]  = v[7];  x[8]  = v[8];
283       v   += 9;
284     }
285     row = *ajtmp++;
286     while (row < i) {
287       pc = rtmp + 9*row;
288       p1 = pc[0];  p2  = pc[1];  p3  = pc[2];  p4  = pc[3];
289       p5 = pc[4];  p6  = pc[5];  p7  = pc[6];  p8  = pc[7];  p9  = pc[8];
290       if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 || p5 != 0.0 ||
291           p6 != 0.0 || p7 != 0.0 || p8 != 0.0 || p9 != 0.0) {
292         pv    = ba + 9*diag_offset[row];
293         pj    = bj + diag_offset[row] + 1;
294         x1    = pv[0];  x2  = pv[1];  x3  = pv[2];  x4  = pv[3];
295         x5    = pv[4];  x6  = pv[5];  x7  = pv[6];  x8  = pv[7];  x9  = pv[8];
296         pc[0] = m1 = p1*x1 + p4*x2 + p7*x3;
297         pc[1] = m2 = p2*x1 + p5*x2 + p8*x3;
298         pc[2] = m3 = p3*x1 + p6*x2 + p9*x3;
299 
300         pc[3] = m4 = p1*x4 + p4*x5 + p7*x6;
301         pc[4] = m5 = p2*x4 + p5*x5 + p8*x6;
302         pc[5] = m6 = p3*x4 + p6*x5 + p9*x6;
303 
304         pc[6] = m7 = p1*x7 + p4*x8 + p7*x9;
305         pc[7] = m8 = p2*x7 + p5*x8 + p8*x9;
306         pc[8] = m9 = p3*x7 + p6*x8 + p9*x9;
307 
308         nz  = bi[row+1] - diag_offset[row] - 1;
309         pv += 9;
310         for (j=0; j<nz; j++) {
311           x1    = pv[0];  x2  = pv[1];   x3 = pv[2];  x4  = pv[3];
312           x5    = pv[4];  x6  = pv[5];   x7 = pv[6];  x8  = pv[7]; x9 = pv[8];
313           x     = rtmp + 9*pj[j];
314           x[0] -= m1*x1 + m4*x2 + m7*x3;
315           x[1] -= m2*x1 + m5*x2 + m8*x3;
316           x[2] -= m3*x1 + m6*x2 + m9*x3;
317 
318           x[3] -= m1*x4 + m4*x5 + m7*x6;
319           x[4] -= m2*x4 + m5*x5 + m8*x6;
320           x[5] -= m3*x4 + m6*x5 + m9*x6;
321 
322           x[6] -= m1*x7 + m4*x8 + m7*x9;
323           x[7] -= m2*x7 + m5*x8 + m8*x9;
324           x[8] -= m3*x7 + m6*x8 + m9*x9;
325           pv   += 9;
326         }
327         ierr = PetscLogFlops(54.0*nz+36.0);CHKERRQ(ierr);
328       }
329       row = *ajtmp++;
330     }
331     /* finished row so stick it into b->a */
332     pv = ba + 9*bi[i];
333     pj = bj + bi[i];
334     nz = bi[i+1] - bi[i];
335     for (j=0; j<nz; j++) {
336       x     = rtmp+9*pj[j];
337       pv[0] = x[0];  pv[1]  = x[1];  pv[2]  = x[2];  pv[3]  = x[3];
338       pv[4] = x[4];  pv[5]  = x[5];  pv[6]  = x[6];  pv[7]  = x[7]; pv[8] = x[8];
339       pv   += 9;
340     }
341     /* invert diagonal block */
342     w    = ba + 9*diag_offset[i];
343     ierr = PetscKernel_A_gets_inverse_A_3(w,shift,!A->erroriffailure,&zeropivotdetected);CHKERRQ(ierr);
344     if (zeropivotdetected) C->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
345   }
346 
347   ierr = PetscFree(rtmp);CHKERRQ(ierr);
348 
349   C->ops->solve          = MatSolve_SeqBAIJ_3_NaturalOrdering_inplace;
350   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_3_NaturalOrdering_inplace;
351   C->assembled           = PETSC_TRUE;
352 
353   ierr = PetscLogFlops(1.333333333333*3*3*3*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
354   PetscFunctionReturn(0);
355 }
356 
357 /*
358   MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering -
359     copied from MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_inplace()
360 */
361 #undef __FUNCT__
362 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering"
363 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering(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   PetscErrorCode ierr;
368   PetscInt       i,j,k,nz,nzL,row;
369   const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
370   const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag,*pj,bs2=a->bs2;
371   MatScalar      *rtmp,*pc,*mwork,*v,*pv,*aa=a->a;
372   PetscInt       flg;
373   PetscReal      shift = info->shiftamount;
374   PetscBool      zeropivotdetected;
375 
376   PetscFunctionBegin;
377   /* generate work space needed by the factorization */
378   ierr = PetscMalloc2(bs2*n,&rtmp,bs2,&mwork);CHKERRQ(ierr);
379   ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr);
380 
381   for (i=0; i<n; i++) {
382     /* zero rtmp */
383     /* L part */
384     nz    = bi[i+1] - bi[i];
385     bjtmp = bj + bi[i];
386     for  (j=0; j<nz; j++) {
387       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
388     }
389 
390     /* U part */
391     nz    = bdiag[i] - bdiag[i+1];
392     bjtmp = bj + bdiag[i+1] + 1;
393     for  (j=0; j<nz; j++) {
394       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
395     }
396 
397     /* load in initial (unfactored row) */
398     nz    = ai[i+1] - ai[i];
399     ajtmp = aj + ai[i];
400     v     = aa + bs2*ai[i];
401     for (j=0; j<nz; j++) {
402       ierr = PetscMemcpy(rtmp+bs2*ajtmp[j],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr);
403     }
404 
405     /* elimination */
406     bjtmp = bj + bi[i];
407     nzL   = bi[i+1] - bi[i];
408     for (k=0; k<nzL; k++) {
409       row = bjtmp[k];
410       pc  = rtmp + bs2*row;
411       for (flg=0,j=0; j<bs2; j++) {
412         if (pc[j]!=0.0) {
413           flg = 1;
414           break;
415         }
416       }
417       if (flg) {
418         pv = b->a + bs2*bdiag[row];
419         /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
420         ierr = PetscKernel_A_gets_A_times_B_3(pc,pv,mwork);CHKERRQ(ierr);
421 
422         pj = b->j + bdiag[row+1]+1; /* beginning of U(row,:) */
423         pv = b->a + bs2*(bdiag[row+1]+1);
424         nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:) excluding diag */
425         for (j=0; j<nz; j++) {
426           /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
427           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
428           v    = rtmp + bs2*pj[j];
429           ierr = PetscKernel_A_gets_A_minus_B_times_C_3(v,pc,pv);CHKERRQ(ierr);
430           pv  += bs2;
431         }
432         ierr = PetscLogFlops(54*nz+45);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
433       }
434     }
435 
436     /* finished row so stick it into b->a */
437     /* L part */
438     pv = b->a + bs2*bi[i];
439     pj = b->j + bi[i];
440     nz = bi[i+1] - bi[i];
441     for (j=0; j<nz; j++) {
442       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
443     }
444 
445     /* Mark diagonal and invert diagonal for simplier triangular solves */
446     pv   = b->a + bs2*bdiag[i];
447     pj   = b->j + bdiag[i];
448     ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr);
449     /* ierr = PetscKernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); */
450     ierr = PetscKernel_A_gets_inverse_A_3(pv,shift,!A->erroriffailure,&zeropivotdetected);CHKERRQ(ierr);
451      if (zeropivotdetected) B->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
452 
453     /* U part */
454     pv = b->a + bs2*(bdiag[i+1]+1);
455     pj = b->j + bdiag[i+1]+1;
456     nz = bdiag[i] - bdiag[i+1] - 1;
457     for (j=0; j<nz; j++) {
458       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
459     }
460   }
461   ierr = PetscFree2(rtmp,mwork);CHKERRQ(ierr);
462 
463   C->ops->solve          = MatSolve_SeqBAIJ_3_NaturalOrdering;
464   C->ops->forwardsolve   = MatForwardSolve_SeqBAIJ_3_NaturalOrdering;
465   C->ops->backwardsolve  = MatBackwardSolve_SeqBAIJ_3_NaturalOrdering;
466   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_3_NaturalOrdering;
467   C->assembled           = PETSC_TRUE;
468 
469   ierr = PetscLogFlops(1.333333333333*3*3*3*n);CHKERRQ(ierr); /* from inverting diagonal blocks */
470   PetscFunctionReturn(0);
471 }
472 
473