xref: /petsc/src/mat/impls/baij/seq/baijfact13.c (revision cd545bacbb74b28ee04a97aa2a6bca79b7d5aa2e)
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 
28   PetscFunctionBegin;
29   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
30   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
31   ierr = PetscMalloc1(9*(n+1),&rtmp);CHKERRQ(ierr);
32 
33   for (i=0; i<n; i++) {
34     nz    = bi[i+1] - bi[i];
35     ajtmp = bj + bi[i];
36     for  (j=0; j<nz; j++) {
37       x    = rtmp + 9*ajtmp[j];
38       x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = 0.0;
39     }
40     /* load in initial (unfactored row) */
41     idx      = r[i];
42     nz       = ai[idx+1] - ai[idx];
43     ajtmpold = aj + ai[idx];
44     v        = aa + 9*ai[idx];
45     for (j=0; j<nz; j++) {
46       x    = rtmp + 9*ic[ajtmpold[j]];
47       x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3];
48       x[4] = v[4]; x[5] = v[5]; x[6] = v[6]; x[7] = v[7]; x[8] = v[8];
49       v   += 9;
50     }
51     row = *ajtmp++;
52     while (row < i) {
53       pc = rtmp + 9*row;
54       p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3];
55       p5 = pc[4]; p6 = pc[5]; p7 = pc[6]; p8 = pc[7]; p9 = pc[8];
56       if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 || p5 != 0.0 ||
57           p6 != 0.0 || p7 != 0.0 || p8 != 0.0 || p9 != 0.0) {
58         pv    = ba + 9*diag_offset[row];
59         pj    = bj + diag_offset[row] + 1;
60         x1    = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
61         x5    = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8];
62         pc[0] = m1 = p1*x1 + p4*x2 + p7*x3;
63         pc[1] = m2 = p2*x1 + p5*x2 + p8*x3;
64         pc[2] = m3 = p3*x1 + p6*x2 + p9*x3;
65 
66         pc[3] = m4 = p1*x4 + p4*x5 + p7*x6;
67         pc[4] = m5 = p2*x4 + p5*x5 + p8*x6;
68         pc[5] = m6 = p3*x4 + p6*x5 + p9*x6;
69 
70         pc[6] = m7 = p1*x7 + p4*x8 + p7*x9;
71         pc[7] = m8 = p2*x7 + p5*x8 + p8*x9;
72         pc[8] = m9 = p3*x7 + p6*x8 + p9*x9;
73         nz    = bi[row+1] - diag_offset[row] - 1;
74         pv   += 9;
75         for (j=0; j<nz; j++) {
76           x1    = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
77           x5    = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8];
78           x     = rtmp + 9*pj[j];
79           x[0] -= m1*x1 + m4*x2 + m7*x3;
80           x[1] -= m2*x1 + m5*x2 + m8*x3;
81           x[2] -= m3*x1 + m6*x2 + m9*x3;
82 
83           x[3] -= m1*x4 + m4*x5 + m7*x6;
84           x[4] -= m2*x4 + m5*x5 + m8*x6;
85           x[5] -= m3*x4 + m6*x5 + m9*x6;
86 
87           x[6] -= m1*x7 + m4*x8 + m7*x9;
88           x[7] -= m2*x7 + m5*x8 + m8*x9;
89           x[8] -= m3*x7 + m6*x8 + m9*x9;
90           pv   += 9;
91         }
92         ierr = PetscLogFlops(54.0*nz+36.0);CHKERRQ(ierr);
93       }
94       row = *ajtmp++;
95     }
96     /* finished row so stick it into b->a */
97     pv = ba + 9*bi[i];
98     pj = bj + bi[i];
99     nz = bi[i+1] - bi[i];
100     for (j=0; j<nz; j++) {
101       x     = rtmp + 9*pj[j];
102       pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3];
103       pv[4] = x[4]; pv[5] = x[5]; pv[6] = x[6]; pv[7] = x[7]; pv[8] = x[8];
104       pv   += 9;
105     }
106     /* invert diagonal block */
107     w    = ba + 9*diag_offset[i];
108     PetscBool wouldcrash;
109     //C
110     ierr = PetscKernel_A_gets_inverse_A_3(w,shift,A->erroriffailure,&wouldcrash);CHKERRQ(ierr);
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 
147   PetscFunctionBegin;
148   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
149   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
150 
151   /* generate work space needed by the factorization */
152   ierr = PetscMalloc2(bs2*n,&rtmp,bs2,&mwork);CHKERRQ(ierr);
153   ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr);
154 
155   for (i=0; i<n; i++) {
156     /* zero rtmp */
157     /* L part */
158     nz    = bi[i+1] - bi[i];
159     bjtmp = bj + bi[i];
160     for  (j=0; j<nz; j++) {
161       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
162     }
163 
164     /* U part */
165     nz    = bdiag[i] - bdiag[i+1];
166     bjtmp = bj + bdiag[i+1]+1;
167     for  (j=0; j<nz; j++) {
168       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
169     }
170 
171     /* load in initial (unfactored row) */
172     nz    = ai[r[i]+1] - ai[r[i]];
173     ajtmp = aj + ai[r[i]];
174     v     = aa + bs2*ai[r[i]];
175     for (j=0; j<nz; j++) {
176       ierr = PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr);
177     }
178 
179     /* elimination */
180     bjtmp = bj + bi[i];
181     nzL   = bi[i+1] - bi[i];
182     for (k = 0; k < nzL; k++) {
183       row = bjtmp[k];
184       pc  = rtmp + bs2*row;
185       for (flg=0,j=0; j<bs2; j++) {
186         if (pc[j]!=0.0) {
187           flg = 1;
188           break;
189         }
190       }
191       if (flg) {
192         pv = b->a + bs2*bdiag[row];
193         /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
194         ierr = PetscKernel_A_gets_A_times_B_3(pc,pv,mwork);CHKERRQ(ierr);
195 
196         pj = b->j + bdiag[row+1] + 1; /* beginning of U(row,:) */
197         pv = b->a + bs2*(bdiag[row+1]+1);
198         nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:) excluding diag */
199         for (j=0; j<nz; j++) {
200           /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
201           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
202           v    = rtmp + bs2*pj[j];
203           ierr = PetscKernel_A_gets_A_minus_B_times_C_3(v,pc,pv);CHKERRQ(ierr);
204           pv  += bs2;
205         }
206         ierr = PetscLogFlops(54*nz+45);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
207       }
208     }
209 
210     /* finished row so stick it into b->a */
211     /* L part */
212     pv = b->a + bs2*bi[i];
213     pj = b->j + bi[i];
214     nz = bi[i+1] - bi[i];
215     for (j=0; j<nz; j++) {
216       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
217     }
218 
219     /* Mark diagonal and invert diagonal for simplier triangular solves */
220     pv   = b->a + bs2*bdiag[i];
221     pj   = b->j + bdiag[i];
222     ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr);
223     /* ierr = PetscKernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); */
224     PetscBool wouldcrash;
225     //B
226     ierr = PetscKernel_A_gets_inverse_A_3(pv,shift,A->erroriffailure,&wouldcrash);CHKERRQ(ierr);
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 
264   PetscFunctionBegin;
265   ierr = PetscMalloc1(9*(n+1),&rtmp);CHKERRQ(ierr);
266 
267   for (i=0; i<n; i++) {
268     nz    = bi[i+1] - bi[i];
269     ajtmp = bj + bi[i];
270     for  (j=0; j<nz; j++) {
271       x    = rtmp+9*ajtmp[j];
272       x[0] = x[1]  = x[2]  = x[3]  = x[4]  = x[5]  = x[6] = x[7] = x[8] = 0.0;
273     }
274     /* load in initial (unfactored row) */
275     nz       = ai[i+1] - ai[i];
276     ajtmpold = aj + ai[i];
277     v        = aa + 9*ai[i];
278     for (j=0; j<nz; j++) {
279       x    = rtmp+9*ajtmpold[j];
280       x[0] = v[0];  x[1]  = v[1];  x[2]  = v[2];  x[3]  = v[3];
281       x[4] = v[4];  x[5]  = v[5];  x[6]  = v[6];  x[7]  = v[7];  x[8]  = v[8];
282       v   += 9;
283     }
284     row = *ajtmp++;
285     while (row < i) {
286       pc = rtmp + 9*row;
287       p1 = pc[0];  p2  = pc[1];  p3  = pc[2];  p4  = pc[3];
288       p5 = pc[4];  p6  = pc[5];  p7  = pc[6];  p8  = pc[7];  p9  = pc[8];
289       if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 || p5 != 0.0 ||
290           p6 != 0.0 || p7 != 0.0 || p8 != 0.0 || p9 != 0.0) {
291         pv    = ba + 9*diag_offset[row];
292         pj    = bj + diag_offset[row] + 1;
293         x1    = pv[0];  x2  = pv[1];  x3  = pv[2];  x4  = pv[3];
294         x5    = pv[4];  x6  = pv[5];  x7  = pv[6];  x8  = pv[7];  x9  = pv[8];
295         pc[0] = m1 = p1*x1 + p4*x2 + p7*x3;
296         pc[1] = m2 = p2*x1 + p5*x2 + p8*x3;
297         pc[2] = m3 = p3*x1 + p6*x2 + p9*x3;
298 
299         pc[3] = m4 = p1*x4 + p4*x5 + p7*x6;
300         pc[4] = m5 = p2*x4 + p5*x5 + p8*x6;
301         pc[5] = m6 = p3*x4 + p6*x5 + p9*x6;
302 
303         pc[6] = m7 = p1*x7 + p4*x8 + p7*x9;
304         pc[7] = m8 = p2*x7 + p5*x8 + p8*x9;
305         pc[8] = m9 = p3*x7 + p6*x8 + p9*x9;
306 
307         nz  = bi[row+1] - diag_offset[row] - 1;
308         pv += 9;
309         for (j=0; j<nz; j++) {
310           x1    = pv[0];  x2  = pv[1];   x3 = pv[2];  x4  = pv[3];
311           x5    = pv[4];  x6  = pv[5];   x7 = pv[6];  x8  = pv[7]; x9 = pv[8];
312           x     = rtmp + 9*pj[j];
313           x[0] -= m1*x1 + m4*x2 + m7*x3;
314           x[1] -= m2*x1 + m5*x2 + m8*x3;
315           x[2] -= m3*x1 + m6*x2 + m9*x3;
316 
317           x[3] -= m1*x4 + m4*x5 + m7*x6;
318           x[4] -= m2*x4 + m5*x5 + m8*x6;
319           x[5] -= m3*x4 + m6*x5 + m9*x6;
320 
321           x[6] -= m1*x7 + m4*x8 + m7*x9;
322           x[7] -= m2*x7 + m5*x8 + m8*x9;
323           x[8] -= m3*x7 + m6*x8 + m9*x9;
324           pv   += 9;
325         }
326         ierr = PetscLogFlops(54.0*nz+36.0);CHKERRQ(ierr);
327       }
328       row = *ajtmp++;
329     }
330     /* finished row so stick it into b->a */
331     pv = ba + 9*bi[i];
332     pj = bj + bi[i];
333     nz = bi[i+1] - bi[i];
334     for (j=0; j<nz; j++) {
335       x     = rtmp+9*pj[j];
336       pv[0] = x[0];  pv[1]  = x[1];  pv[2]  = x[2];  pv[3]  = x[3];
337       pv[4] = x[4];  pv[5]  = x[5];  pv[6]  = x[6];  pv[7]  = x[7]; pv[8] = x[8];
338       pv   += 9;
339     }
340     /* invert diagonal block */
341     w    = ba + 9*diag_offset[i];
342     PetscBool wouldcrash;
343     //C
344     ierr = PetscKernel_A_gets_inverse_A_3(w,shift,A->erroriffailure,&wouldcrash);CHKERRQ(ierr);
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 
375   PetscFunctionBegin;
376   /* generate work space needed by the factorization */
377   ierr = PetscMalloc2(bs2*n,&rtmp,bs2,&mwork);CHKERRQ(ierr);
378   ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr);
379 
380   for (i=0; i<n; i++) {
381     /* zero rtmp */
382     /* L part */
383     nz    = bi[i+1] - bi[i];
384     bjtmp = bj + bi[i];
385     for  (j=0; j<nz; j++) {
386       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
387     }
388 
389     /* U part */
390     nz    = bdiag[i] - bdiag[i+1];
391     bjtmp = bj + bdiag[i+1] + 1;
392     for  (j=0; j<nz; j++) {
393       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
394     }
395 
396     /* load in initial (unfactored row) */
397     nz    = ai[i+1] - ai[i];
398     ajtmp = aj + ai[i];
399     v     = aa + bs2*ai[i];
400     for (j=0; j<nz; j++) {
401       ierr = PetscMemcpy(rtmp+bs2*ajtmp[j],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr);
402     }
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         ierr = PetscKernel_A_gets_A_times_B_3(pc,pv,mwork);CHKERRQ(ierr);
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 in U(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           ierr = PetscKernel_A_gets_A_minus_B_times_C_3(v,pc,pv);CHKERRQ(ierr);
429           pv  += bs2;
430         }
431         ierr = PetscLogFlops(54*nz+45);CHKERRQ(ierr); /* 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++) {
441       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
442     }
443 
444     /* Mark diagonal and invert diagonal for simplier triangular solves */
445     pv   = b->a + bs2*bdiag[i];
446     pj   = b->j + bdiag[i];
447     ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr);
448     /* ierr = PetscKernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); */
449     PetscBool wouldcrash;
450     //B
451     ierr = PetscKernel_A_gets_inverse_A_3(pv,shift,A->erroriffailure,&wouldcrash);CHKERRQ(ierr);
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