xref: /petsc/src/mat/impls/baij/seq/baijfact13.c (revision 2205254efee3a00a594e5e2a3a70f74dcb40bc03)
1 
2 /*
3     Factorization code for BAIJ format.
4 */
5 #include <../src/mat/impls/baij/seq/baij.h>
6 #include <../src/mat/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 = PetscMalloc(9*(n+1)*sizeof(MatScalar),&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     ierr = PetscKernel_A_gets_inverse_A_3(w,shift);CHKERRQ(ierr);
109   }
110 
111   ierr = PetscFree(rtmp);CHKERRQ(ierr);
112   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
113   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
114   C->ops->solve          = MatSolve_SeqBAIJ_3_inplace;
115   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_3_inplace;
116   C->assembled = PETSC_TRUE;
117   ierr = PetscLogFlops(1.333333333333*3*3*3*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
118   PetscFunctionReturn(0);
119 }
120 
121 /* MatLUFactorNumeric_SeqBAIJ_3 -
122      copied from MatLUFactorNumeric_SeqBAIJ_N_inplace() and manually re-implemented
123        PetscKernel_A_gets_A_times_B()
124        PetscKernel_A_gets_A_minus_B_times_C()
125        PetscKernel_A_gets_inverse_A()
126 */
127 #undef __FUNCT__
128 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_3"
129 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3(Mat B,Mat A,const MatFactorInfo *info)
130 {
131   Mat            C=B;
132   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ *)C->data;
133   IS             isrow = b->row,isicol = b->icol;
134   PetscErrorCode ierr;
135   const PetscInt *r,*ic;
136   PetscInt       i,j,k,nz,nzL,row;
137   const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
138   const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag,*pj,bs2=a->bs2;
139   MatScalar      *rtmp,*pc,*mwork,*v,*pv,*aa=a->a;
140   PetscInt       flg;
141   PetscReal      shift = info->shiftamount;
142 
143   PetscFunctionBegin;
144   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
145   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
146 
147   /* generate work space needed by the factorization */
148   ierr = PetscMalloc2(bs2*n,MatScalar,&rtmp,bs2,MatScalar,&mwork);CHKERRQ(ierr);
149   ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr);
150 
151   for (i=0; i<n; i++) {
152     /* zero rtmp */
153     /* L part */
154     nz    = bi[i+1] - bi[i];
155     bjtmp = bj + bi[i];
156     for  (j=0; j<nz; j++) {
157       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
158     }
159 
160     /* U part */
161     nz = bdiag[i] - bdiag[i+1];
162     bjtmp = bj + bdiag[i+1]+1;
163     for  (j=0; j<nz; j++) {
164       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
165     }
166 
167     /* load in initial (unfactored row) */
168     nz    = ai[r[i]+1] - ai[r[i]];
169     ajtmp = aj + ai[r[i]];
170     v     = aa + bs2*ai[r[i]];
171     for (j=0; j<nz; j++) {
172       ierr = PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr);
173     }
174 
175     /* elimination */
176     bjtmp = bj + bi[i];
177     nzL   = bi[i+1] - bi[i];
178     for (k = 0;k < nzL;k++) {
179       row = bjtmp[k];
180       pc = rtmp + bs2*row;
181       for (flg=0,j=0; j<bs2; j++) {
182         if (pc[j]!=0.0) {
183           flg = 1;
184           break;
185         }
186       }
187       if (flg) {
188         pv = b->a + bs2*bdiag[row];
189         /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
190         ierr = PetscKernel_A_gets_A_times_B_3(pc,pv,mwork);CHKERRQ(ierr);
191 
192         pj = b->j + bdiag[row+1] + 1; /* beginning of U(row,:) */
193         pv = b->a + bs2*(bdiag[row+1]+1);
194         nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:) excluding diag */
195         for (j=0; j<nz; j++) {
196           /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
197           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
198           v    = rtmp + bs2*pj[j];
199           ierr = PetscKernel_A_gets_A_minus_B_times_C_3(v,pc,pv);CHKERRQ(ierr);
200           pv  += bs2;
201         }
202         ierr = PetscLogFlops(54*nz+45);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
203       }
204     }
205 
206     /* finished row so stick it into b->a */
207     /* L part */
208     pv   = b->a + bs2*bi[i] ;
209     pj   = b->j + bi[i] ;
210     nz   = bi[i+1] - bi[i];
211     for (j=0; j<nz; j++) {
212       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
213     }
214 
215     /* Mark diagonal and invert diagonal for simplier triangular solves */
216     pv   = b->a + bs2*bdiag[i];
217     pj   = b->j + bdiag[i];
218     ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr);
219     /* ierr = PetscKernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); */
220     ierr = PetscKernel_A_gets_inverse_A_3(pv,shift);CHKERRQ(ierr);
221 
222     /* U part */
223     pj = b->j + bdiag[i+1] + 1;
224     pv = b->a + bs2*(bdiag[i+1]+1);
225     nz = bdiag[i] - bdiag[i+1] - 1;
226     for (j=0; j<nz; j++) {
227       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
228     }
229   }
230 
231   ierr = PetscFree2(rtmp,mwork);CHKERRQ(ierr);
232   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
233   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
234   C->ops->solve = MatSolve_SeqBAIJ_3;
235   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_3;
236 
237   C->assembled = PETSC_TRUE;
238   ierr = PetscLogFlops(1.333333333333*3*3*3*n);CHKERRQ(ierr); /* from inverting diagonal blocks */
239   PetscFunctionReturn(0);
240 }
241 
242 #undef __FUNCT__
243 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering_inplace"
244 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering_inplace(Mat C,Mat A,const MatFactorInfo *info)
245 {
246   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data;
247   PetscErrorCode ierr;
248   PetscInt       i,j,n = a->mbs,*bi = b->i,*bj = b->j;
249   PetscInt       *ajtmpold,*ajtmp,nz,row;
250   PetscInt       *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj;
251   MatScalar      *pv,*v,*rtmp,*pc,*w,*x;
252   MatScalar      p1,p2,p3,p4,m1,m2,m3,m4,m5,m6,m7,m8,m9,x1,x2,x3,x4;
253   MatScalar      p5,p6,p7,p8,p9,x5,x6,x7,x8,x9;
254   MatScalar      *ba = b->a,*aa = a->a;
255   PetscReal      shift = info->shiftamount;
256 
257   PetscFunctionBegin;
258   ierr = PetscMalloc(9*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
259 
260   for (i=0; i<n; i++) {
261     nz    = bi[i+1] - bi[i];
262     ajtmp = bj + bi[i];
263     for  (j=0; j<nz; j++) {
264       x = rtmp+9*ajtmp[j];
265       x[0]  = x[1]  = x[2]  = x[3]  = x[4]  = x[5]  = x[6] = x[7] = x[8] = 0.0;
266     }
267     /* load in initial (unfactored row) */
268     nz       = ai[i+1] - ai[i];
269     ajtmpold = aj + ai[i];
270     v        = aa + 9*ai[i];
271     for (j=0; j<nz; j++) {
272       x    = rtmp+9*ajtmpold[j];
273       x[0]  = v[0];  x[1]  = v[1];  x[2]  = v[2];  x[3]  = v[3];
274       x[4]  = v[4];  x[5]  = v[5];  x[6]  = v[6];  x[7]  = v[7];  x[8]  = v[8];
275       v    += 9;
276     }
277     row = *ajtmp++;
278     while (row < i) {
279       pc  = rtmp + 9*row;
280       p1  = pc[0];  p2  = pc[1];  p3  = pc[2];  p4  = pc[3];
281       p5  = pc[4];  p6  = pc[5];  p7  = pc[6];  p8  = pc[7];  p9  = pc[8];
282       if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 || p5 != 0.0 ||
283           p6 != 0.0 || p7 != 0.0 || p8 != 0.0 || p9 != 0.0) {
284         pv = ba + 9*diag_offset[row];
285         pj = bj + diag_offset[row] + 1;
286         x1  = pv[0];  x2  = pv[1];  x3  = pv[2];  x4  = pv[3];
287         x5  = pv[4];  x6  = pv[5];  x7  = pv[6];  x8  = pv[7];  x9  = pv[8];
288         pc[0] = m1 = p1*x1 + p4*x2 + p7*x3;
289         pc[1] = m2 = p2*x1 + p5*x2 + p8*x3;
290         pc[2] = m3 = p3*x1 + p6*x2 + p9*x3;
291 
292         pc[3] = m4 = p1*x4 + p4*x5 + p7*x6;
293         pc[4] = m5 = p2*x4 + p5*x5 + p8*x6;
294         pc[5] = m6 = p3*x4 + p6*x5 + p9*x6;
295 
296         pc[6] = m7 = p1*x7 + p4*x8 + p7*x9;
297         pc[7] = m8 = p2*x7 + p5*x8 + p8*x9;
298         pc[8] = m9 = p3*x7 + p6*x8 + p9*x9;
299 
300         nz = bi[row+1] - diag_offset[row] - 1;
301         pv += 9;
302         for (j=0; j<nz; j++) {
303           x1   = pv[0];  x2  = pv[1];   x3 = pv[2];  x4  = pv[3];
304           x5   = pv[4];  x6  = pv[5];   x7 = pv[6];  x8  = pv[7]; x9 = pv[8];
305           x    = rtmp + 9*pj[j];
306           x[0] -= m1*x1 + m4*x2 + m7*x3;
307           x[1] -= m2*x1 + m5*x2 + m8*x3;
308           x[2] -= m3*x1 + m6*x2 + m9*x3;
309 
310           x[3] -= m1*x4 + m4*x5 + m7*x6;
311           x[4] -= m2*x4 + m5*x5 + m8*x6;
312           x[5] -= m3*x4 + m6*x5 + m9*x6;
313 
314           x[6] -= m1*x7 + m4*x8 + m7*x9;
315           x[7] -= m2*x7 + m5*x8 + m8*x9;
316           x[8] -= m3*x7 + m6*x8 + m9*x9;
317           pv   += 9;
318         }
319         ierr = PetscLogFlops(54.0*nz+36.0);CHKERRQ(ierr);
320       }
321       row = *ajtmp++;
322     }
323     /* finished row so stick it into b->a */
324     pv = ba + 9*bi[i];
325     pj = bj + bi[i];
326     nz = bi[i+1] - bi[i];
327     for (j=0; j<nz; j++) {
328       x      = rtmp+9*pj[j];
329       pv[0]  = x[0];  pv[1]  = x[1];  pv[2]  = x[2];  pv[3]  = x[3];
330       pv[4]  = x[4];  pv[5]  = x[5];  pv[6]  = x[6];  pv[7]  = x[7]; pv[8] = x[8];
331       pv   += 9;
332     }
333     /* invert diagonal block */
334     w = ba + 9*diag_offset[i];
335     ierr = PetscKernel_A_gets_inverse_A_3(w,shift);CHKERRQ(ierr);
336   }
337 
338   ierr = PetscFree(rtmp);CHKERRQ(ierr);
339   C->ops->solve          = MatSolve_SeqBAIJ_3_NaturalOrdering_inplace;
340   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_3_NaturalOrdering_inplace;
341   C->assembled = PETSC_TRUE;
342   ierr = PetscLogFlops(1.333333333333*3*3*3*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
343   PetscFunctionReturn(0);
344 }
345 
346 /*
347   MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering -
348     copied from MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_inplace()
349 */
350 #undef __FUNCT__
351 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering"
352 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info)
353 {
354   Mat            C=B;
355   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ *)C->data;
356   PetscErrorCode ierr;
357   PetscInt       i,j,k,nz,nzL,row;
358   const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
359   const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag,*pj,bs2=a->bs2;
360   MatScalar      *rtmp,*pc,*mwork,*v,*pv,*aa=a->a;
361   PetscInt       flg;
362   PetscReal      shift = info->shiftamount;
363 
364   PetscFunctionBegin;
365   /* generate work space needed by the factorization */
366   ierr = PetscMalloc2(bs2*n,MatScalar,&rtmp,bs2,MatScalar,&mwork);CHKERRQ(ierr);
367   ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr);
368 
369   for (i=0; i<n; i++) {
370     /* zero rtmp */
371     /* L part */
372     nz    = bi[i+1] - bi[i];
373     bjtmp = bj + bi[i];
374     for  (j=0; j<nz; j++) {
375       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
376     }
377 
378     /* U part */
379     nz = bdiag[i] - bdiag[i+1];
380     bjtmp = bj + bdiag[i+1] + 1;
381     for  (j=0; j<nz; j++) {
382       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
383     }
384 
385     /* load in initial (unfactored row) */
386     nz    = ai[i+1] - ai[i];
387     ajtmp = aj + ai[i];
388     v     = aa + bs2*ai[i];
389     for (j=0; j<nz; j++) {
390       ierr = PetscMemcpy(rtmp+bs2*ajtmp[j],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr);
391     }
392 
393     /* elimination */
394     bjtmp = bj + bi[i];
395     nzL   = bi[i+1] - bi[i];
396     for (k=0;k<nzL;k++) {
397       row = bjtmp[k];
398       pc = rtmp + bs2*row;
399       for (flg=0,j=0; j<bs2; j++) {
400         if (pc[j]!=0.0) {
401           flg = 1;
402           break;
403         }
404       }
405       if (flg) {
406         pv = b->a + bs2*bdiag[row];
407         /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
408         ierr = PetscKernel_A_gets_A_times_B_3(pc,pv,mwork);CHKERRQ(ierr);
409 
410         pj = b->j + bdiag[row+1]+1; /* beginning of U(row,:) */
411         pv = b->a + bs2*(bdiag[row+1]+1);
412         nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:) excluding diag */
413         for (j=0; j<nz; j++) {
414           /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
415           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
416           v    = rtmp + bs2*pj[j];
417           ierr = PetscKernel_A_gets_A_minus_B_times_C_3(v,pc,pv);CHKERRQ(ierr);
418           pv  += bs2;
419         }
420         ierr = PetscLogFlops(54*nz+45);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
421       }
422     }
423 
424     /* finished row so stick it into b->a */
425     /* L part */
426     pv   = b->a + bs2*bi[i] ;
427     pj   = b->j + bi[i] ;
428     nz   = bi[i+1] - bi[i];
429     for (j=0; j<nz; j++) {
430       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
431     }
432 
433     /* Mark diagonal and invert diagonal for simplier triangular solves */
434     pv   = b->a + bs2*bdiag[i];
435     pj   = b->j + bdiag[i];
436     ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr);
437     /* ierr = PetscKernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); */
438     ierr = PetscKernel_A_gets_inverse_A_3(pv,shift);CHKERRQ(ierr);
439 
440     /* U part */
441     pv = b->a + bs2*(bdiag[i+1]+1);
442     pj = b->j + bdiag[i+1]+1;
443     nz = bdiag[i] - bdiag[i+1] - 1;
444     for (j=0; j<nz; j++) {
445       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
446     }
447   }
448   ierr = PetscFree2(rtmp,mwork);CHKERRQ(ierr);
449   C->ops->solve          = MatSolve_SeqBAIJ_3_NaturalOrdering;
450   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_3_NaturalOrdering;
451   C->assembled = PETSC_TRUE;
452   ierr = PetscLogFlops(1.333333333333*3*3*3*n);CHKERRQ(ierr); /* from inverting diagonal blocks */
453   PetscFunctionReturn(0);
454 }
455 
456