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