xref: /petsc/src/mat/impls/baij/seq/baijfact.c (revision 3a40ed3dce77c081171d005ae1a6ff4bb9d13b6f)
1 #ifdef PETSC_RCS_HEADER
2 static char vcid[] = "$Id: baijfact.c,v 1.52 1997/07/26 19:03:52 bsmith Exp bsmith $";
3 #endif
4 /*
5     Factorization code for BAIJ format.
6 */
7 
8 #include "src/mat/impls/baij/seq/baij.h"
9 #include "src/vec/vecimpl.h"
10 #include "src/inline/ilu.h"
11 
12 
13 /*
14     The symbolic factorization code is identical to that for AIJ format,
15   except for very small changes since this is now a SeqBAIJ datastructure.
16   NOT good code reuse.
17 */
18 #undef __FUNC__
19 #define __FUNC__ "MatLUFactorSymbolic_SeqBAIJ"
20 int MatLUFactorSymbolic_SeqBAIJ(Mat A,IS isrow,IS iscol,double f,Mat *B)
21 {
22   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data, *b;
23   IS          isicol;
24   int         *r,*ic, ierr, i, n = a->mbs, *ai = a->i, *aj = a->j;
25   int         *ainew,*ajnew, jmax,*fill, *ajtmp, nz, bs = a->bs, bs2=a->bs2;
26   int         *idnew, idx, row,m,fm, nnz, nzi,realloc = 0,nzbd,*im;
27 
28   PetscFunctionBegin;
29   PetscValidHeaderSpecific(isrow,IS_COOKIE);
30   PetscValidHeaderSpecific(iscol,IS_COOKIE);
31   ierr = ISInvertPermutation(iscol,&isicol); CHKERRQ(ierr);
32   ISGetIndices(isrow,&r); ISGetIndices(isicol,&ic);
33 
34   /* get new row pointers */
35   ainew = (int *) PetscMalloc( (n+1)*sizeof(int) ); CHKPTRQ(ainew);
36   ainew[0] = 0;
37   /* don't know how many column pointers are needed so estimate */
38   jmax = (int) (f*ai[n] + 1);
39   ajnew = (int *) PetscMalloc( (jmax)*sizeof(int) ); CHKPTRQ(ajnew);
40   /* fill is a linked list of nonzeros in active row */
41   fill = (int *) PetscMalloc( (2*n+1)*sizeof(int)); CHKPTRQ(fill);
42   im = fill + n + 1;
43   /* idnew is location of diagonal in factor */
44   idnew = (int *) PetscMalloc( (n+1)*sizeof(int)); CHKPTRQ(idnew);
45   idnew[0] = 0;
46 
47   for ( i=0; i<n; i++ ) {
48     /* first copy previous fill into linked list */
49     nnz     = nz    = ai[r[i]+1] - ai[r[i]];
50     if (!nz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,1,"Empty row in matrix");
51     ajtmp   = aj + ai[r[i]];
52     fill[n] = n;
53     while (nz--) {
54       fm  = n;
55       idx = ic[*ajtmp++];
56       do {
57         m  = fm;
58         fm = fill[m];
59       } while (fm < idx);
60       fill[m]   = idx;
61       fill[idx] = fm;
62     }
63     row = fill[n];
64     while ( row < i ) {
65       ajtmp = ajnew + idnew[row] + 1;
66       nzbd  = 1 + idnew[row] - ainew[row];
67       nz    = im[row] - nzbd;
68       fm    = row;
69       while (nz-- > 0) {
70         idx = *ajtmp++;
71         nzbd++;
72         if (idx == i) im[row] = nzbd;
73         do {
74           m  = fm;
75           fm = fill[m];
76         } while (fm < idx);
77         if (fm != idx) {
78           fill[m]   = idx;
79           fill[idx] = fm;
80           fm        = idx;
81           nnz++;
82         }
83       }
84       row = fill[row];
85     }
86     /* copy new filled row into permanent storage */
87     ainew[i+1] = ainew[i] + nnz;
88     if (ainew[i+1] > jmax) {
89       /* allocate a longer ajnew */
90       int maxadd;
91       maxadd = (int) ((f*(ai[n]+1)*(n-i+5))/n);
92       if (maxadd < nnz) maxadd = (n-i)*(nnz+1);
93       jmax += maxadd;
94       ajtmp = (int *) PetscMalloc( jmax*sizeof(int) );CHKPTRQ(ajtmp);
95       PetscMemcpy(ajtmp,ajnew,ainew[i]*sizeof(int));
96       PetscFree(ajnew);
97       ajnew = ajtmp;
98       realloc++; /* count how many times we realloc */
99     }
100     ajtmp = ajnew + ainew[i];
101     fm    = fill[n];
102     nzi   = 0;
103     im[i] = nnz;
104     while (nnz--) {
105       if (fm < i) nzi++;
106       *ajtmp++ = fm;
107       fm       = fill[fm];
108     }
109     idnew[i] = ainew[i] + nzi;
110   }
111 
112   if (ai[n] != 0) {
113     double af = ((double)ainew[n])/((double)ai[n]);
114     PLogInfo(A,"Info:MatLUFactorSymbolic_SeqBAIJ:Reallocs %d Fill ratio:given %g needed %g\n",
115              realloc,f,af);
116     PLogInfo(A,"Info:MatLUFactorSymbolic_SeqBAIJ:Run with -pc_lu_fill %g or use \n",af);
117     PLogInfo(A,"Info:MatLUFactorSymbolic_SeqBAIJ:PCLUSetFill(pc,%g);\n",af);
118     PLogInfo(A,"Info:MatLUFactorSymbolic_SeqBAIJ:for best performance.\n");
119   } else {
120      PLogInfo(A,"Info:MatLUFactorSymbolic_SeqBAIJ:Empty matrix.\n");
121   }
122 
123   ierr = ISRestoreIndices(isrow,&r); CHKERRQ(ierr);
124   ierr = ISRestoreIndices(isicol,&ic); CHKERRQ(ierr);
125 
126   PetscFree(fill);
127 
128   /* put together the new matrix */
129   ierr = MatCreateSeqBAIJ(A->comm,bs,bs*n,bs*n,0,PETSC_NULL,B); CHKERRQ(ierr);
130   PLogObjectParent(*B,isicol);
131   ierr = ISDestroy(isicol); CHKERRQ(ierr);
132   b = (Mat_SeqBAIJ *) (*B)->data;
133   PetscFree(b->imax);
134   b->singlemalloc = 0;
135   /* the next line frees the default space generated by the Create() */
136   PetscFree(b->a); PetscFree(b->ilen);
137   b->a          = (Scalar *) PetscMalloc((ainew[n]+1)*sizeof(Scalar)*bs2);CHKPTRQ(b->a);
138   b->j          = ajnew;
139   b->i          = ainew;
140   b->diag       = idnew;
141   b->ilen       = 0;
142   b->imax       = 0;
143   b->row        = isrow;
144   b->col        = iscol;
145   b->solve_work = (Scalar *) PetscMalloc( (bs*n+bs)*sizeof(Scalar));CHKPTRQ(b->solve_work);
146   /* In b structure:  Free imax, ilen, old a, old j.
147      Allocate idnew, solve_work, new a, new j */
148   PLogObjectMemory(*B,(ainew[n]-n)*(sizeof(int)+sizeof(Scalar)));
149   b->maxnz = b->nz = ainew[n];
150 
151   (*B)->info.factor_mallocs    = realloc;
152   (*B)->info.fill_ratio_given  = f;
153   if (ai[i] != 0) {
154     (*B)->info.fill_ratio_needed = ((double)ainew[n])/((double)ai[i]);
155   } else {
156     (*B)->info.fill_ratio_needed = 0.0;
157   }
158 
159 
160   PetscFunctionReturn(0);
161 }
162 
163 /* ----------------------------------------------------------- */
164 #undef __FUNC__
165 #define __FUNC__ "MatLUFactorNumeric_SeqBAIJ_N"
166 int MatLUFactorNumeric_SeqBAIJ_N(Mat A,Mat *B)
167 {
168   Mat             C = *B;
169   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data,*b = (Mat_SeqBAIJ *)C->data;
170   IS              iscol = b->col, isrow = b->row, isicol;
171   int             *r,*ic, ierr, i, j, n = a->mbs, *bi = b->i, *bj = b->j;
172   int             *ajtmpold, *ajtmp, nz, row, bslog,*ai=a->i,*aj=a->j,k,flg;
173   int             *diag_offset=b->diag,diag,bs=a->bs,bs2 = a->bs2,*v_pivots;
174   register int    *pj;
175   register Scalar *pv,*v,*rtmp,*multiplier,*v_work,*pc,*w;
176   Scalar          *ba = b->a,*aa = a->a;
177 
178   PetscFunctionBegin;
179   ierr  = ISInvertPermutation(iscol,&isicol); CHKERRQ(ierr);
180   PLogObjectParent(*B,isicol);
181   ierr  = ISGetIndices(isrow,&r); CHKERRQ(ierr);
182   ierr  = ISGetIndices(isicol,&ic); CHKERRQ(ierr);
183   rtmp  = (Scalar *) PetscMalloc(bs2*(n+1)*sizeof(Scalar));CHKPTRQ(rtmp);
184   PetscMemzero(rtmp,bs2*(n+1)*sizeof(Scalar));
185   /* generate work space needed by dense LU factorization */
186   v_work     = (Scalar *) PetscMalloc(bs*sizeof(int) + (bs+bs2)*sizeof(Scalar));
187                CHKPTRQ(v_work);
188   multiplier = v_work + bs;
189   v_pivots   = (int *) (multiplier + bs2);
190 
191   /* flops in while loop */
192   bslog = 2*bs*bs2;
193 
194   for ( i=0; i<n; i++ ) {
195     nz    = bi[i+1] - bi[i];
196     ajtmp = bj + bi[i];
197     for  ( j=0; j<nz; j++ ) {
198       PetscMemzero(rtmp+bs2*ajtmp[j],bs2*sizeof(Scalar));
199     }
200     /* load in initial (unfactored row) */
201     nz       = ai[r[i]+1] - ai[r[i]];
202     ajtmpold = aj + ai[r[i]];
203     v        = aa + bs2*ai[r[i]];
204     for ( j=0; j<nz; j++ ) {
205       PetscMemcpy(rtmp+bs2*ic[ajtmpold[j]],v+bs2*j,bs2*sizeof(Scalar));
206     }
207     row = *ajtmp++;
208     while (row < i) {
209       pc = rtmp + bs2*row;
210 /*      if (*pc) { */
211       for ( flg=0,k=0; k<bs2; k++ ) { if (pc[k]!=0.0) { flg =1; break; }}
212       if (flg) {
213         pv = ba + bs2*diag_offset[row];
214         pj = bj + diag_offset[row] + 1;
215         Kernel_A_gets_A_times_B(bs,pc,pv,multiplier);
216         nz = bi[row+1] - diag_offset[row] - 1;
217         pv += bs2;
218         for (j=0; j<nz; j++) {
219           Kernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j);
220         }
221         PLogFlops(bslog*(nz+1)-bs);
222       }
223         row = *ajtmp++;
224     }
225     /* finished row so stick it into b->a */
226     pv = ba + bs2*bi[i];
227     pj = bj + bi[i];
228     nz = bi[i+1] - bi[i];
229     for ( j=0; j<nz; j++ ) {
230       PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(Scalar));
231     }
232     diag = diag_offset[i] - bi[i];
233     /* invert diagonal block */
234     w = pv + bs2*diag;
235     Kernel_A_gets_inverse_A(bs,w,v_pivots,v_work);
236   }
237 
238   PetscFree(rtmp); PetscFree(v_work);
239   ierr = ISRestoreIndices(isicol,&ic); CHKERRQ(ierr);
240   ierr = ISRestoreIndices(isrow,&r); CHKERRQ(ierr);
241   ierr = ISDestroy(isicol); CHKERRQ(ierr);
242   C->factor = FACTOR_LU;
243   C->assembled = PETSC_TRUE;
244   PLogFlops(1.3333*bs*bs2*b->mbs); /* from inverting diagonal blocks */
245   PetscFunctionReturn(0);
246 }
247 /* ------------------------------------------------------------*/
248 /*
249       Version for when blocks are 5 by 5
250 */
251 #undef __FUNC__
252 #define __FUNC__ "MatLUFactorNumeric_SeqBAIJ_5"
253 int MatLUFactorNumeric_SeqBAIJ_5(Mat A,Mat *B)
254 {
255   Mat             C = *B;
256   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data,*b = (Mat_SeqBAIJ *)C->data;
257   IS              iscol = b->col, isrow = b->row, isicol;
258   int             *r,*ic, ierr, i, j, n = a->mbs, *bi = b->i, *bj = b->j;
259   int             *ajtmpold, *ajtmp, nz, row;
260   int             *diag_offset = b->diag,idx,*ai=a->i,*aj=a->j;
261   register int    *pj;
262   register Scalar *pv,*v,*rtmp,*pc,*w,*x;
263   Scalar          p1,p2,p3,p4,m1,m2,m3,m4,m5,m6,m7,m8,m9,x1,x2,x3,x4;
264   Scalar          p5,p6,p7,p8,p9,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15,x16;
265   Scalar          x17,x18,x19,x20,x21,x22,x23,x24,x25,p10,p11,p12,p13,p14;
266   Scalar          p15,p16,p17,p18,p19,p20,p21,p22,p23,p24,p25,m10,m11,m12;
267   Scalar          m13,m14,m15,m16,m17,m18,m19,m20,m21,m22,m23,m24,m25;
268   Scalar          *ba = b->a,*aa = a->a;
269 
270   PetscFunctionBegin;
271   ierr  = ISInvertPermutation(iscol,&isicol); CHKERRQ(ierr);
272   PLogObjectParent(*B,isicol);
273   ierr  = ISGetIndices(isrow,&r); CHKERRQ(ierr);
274   ierr  = ISGetIndices(isicol,&ic); CHKERRQ(ierr);
275   rtmp  = (Scalar *) PetscMalloc(25*(n+1)*sizeof(Scalar));CHKPTRQ(rtmp);
276 
277   for ( i=0; i<n; i++ ) {
278     nz    = bi[i+1] - bi[i];
279     ajtmp = bj + bi[i];
280     for  ( j=0; j<nz; j++ ) {
281       x = rtmp+25*ajtmp[j];
282       x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = x[9] = 0.0;
283       x[10] = x[11] = x[12] = x[13] = x[14] = x[15] = x[16] = x[17] = 0.0;
284       x[18] = x[19] = x[20] = x[21] = x[22] = x[23] = x[24] = 0.0;
285     }
286     /* load in initial (unfactored row) */
287     idx      = r[i];
288     nz       = ai[idx+1] - ai[idx];
289     ajtmpold = aj + ai[idx];
290     v        = aa + 25*ai[idx];
291     for ( j=0; j<nz; j++ ) {
292       x    = rtmp+25*ic[ajtmpold[j]];
293       x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3];
294       x[4] = v[4]; x[5] = v[5]; x[6] = v[6]; x[7] = v[7]; x[8] = v[8];
295       x[9] = v[9]; x[10] = v[10]; x[11] = v[11]; x[12] = v[12]; x[13] = v[13];
296       x[14] = v[14]; x[15] = v[15]; x[16] = v[16]; x[17] = v[17];
297       x[18] = v[18]; x[19] = v[19]; x[20] = v[20]; x[21] = v[21];
298       x[22] = v[22]; x[23] = v[23]; x[24] = v[24];
299       v    += 25;
300     }
301     row = *ajtmp++;
302     while (row < i) {
303       pc = rtmp + 25*row;
304       p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3];
305       p5 = pc[4]; p6 = pc[5]; p7 = pc[6]; p8 = pc[7]; p9 = pc[8];
306       p10 = pc[9]; p11 = pc[10]; p12 = pc[11]; p13 = pc[12]; p14 = pc[13];
307       p15 = pc[14]; p16 = pc[15]; p17 = pc[16]; p18 = pc[17]; p19 = pc[18];
308       p20 = pc[19]; p21 = pc[20]; p22 = pc[21]; p23 = pc[22]; p24 = pc[23];
309       p25 = pc[24];
310       if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 || p5 != 0.0 ||
311           p6 != 0.0 || p7 != 0.0 || p8 != 0.0 || p9 != 0.0 || p10 != 0.0 ||
312           p11 != 0.0 || p12 != 0.0 || p13 != 0.0 || p14 != 0.0 || p15 != 0.0
313           || p16 != 0.0 || p17 != 0.0 || p18 != 0.0 || p19 != 0.0 ||
314           p20 != 0.0 || p21 != 0.0 || p22 != 0.0 || p23 != 0.0 ||
315           p24 != 0.0 || p25 != 0.0) {
316         pv = ba + 25*diag_offset[row];
317         pj = bj + diag_offset[row] + 1;
318         x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
319         x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8];
320         x10 = pv[9]; x11 = pv[10]; x12 = pv[11]; x13 = pv[12]; x14 = pv[13];
321         x15 = pv[14]; x16 = pv[15]; x17 = pv[16]; x18 = pv[17];
322         x19 = pv[18]; x20 = pv[19]; x21 = pv[20]; x22 = pv[21];
323         x23 = pv[22]; x24 = pv[23]; x25 = pv[24];
324         pc[0] = m1 = p1*x1 + p6*x2  + p11*x3 + p16*x4 + p21*x5;
325         pc[1] = m2 = p2*x1 + p7*x2  + p12*x3 + p17*x4 + p22*x5;
326         pc[2] = m3 = p3*x1 + p8*x2  + p13*x3 + p18*x4 + p23*x5;
327         pc[3] = m4 = p4*x1 + p9*x2  + p14*x3 + p19*x4 + p24*x5;
328         pc[4] = m5 = p5*x1 + p10*x2 + p15*x3 + p20*x4 + p25*x5;
329 
330         pc[5] = m6 = p1*x6 + p6*x7  + p11*x8 + p16*x9 + p21*x10;
331         pc[6] = m7 = p2*x6 + p7*x7  + p12*x8 + p17*x9 + p22*x10;
332         pc[7] = m8 = p3*x6 + p8*x7  + p13*x8 + p18*x9 + p23*x10;
333         pc[8] = m9 = p4*x6 + p9*x7  + p14*x8 + p19*x9 + p24*x10;
334         pc[9] = m10 = p5*x6 + p10*x7 + p15*x8 + p20*x9 + p25*x10;
335 
336         pc[10] = m11 = p1*x11 + p6*x12  + p11*x13 + p16*x14 + p21*x15;
337         pc[11] = m12 = p2*x11 + p7*x12  + p12*x13 + p17*x14 + p22*x15;
338         pc[12] = m13 = p3*x11 + p8*x12  + p13*x13 + p18*x14 + p23*x15;
339         pc[13] = m14 = p4*x11 + p9*x12  + p14*x13 + p19*x14 + p24*x15;
340         pc[14] = m15 = p5*x11 + p10*x12 + p15*x13 + p20*x14 + p25*x15;
341 
342         pc[15] = m16 = p1*x16 + p6*x17  + p11*x18 + p16*x19 + p21*x20;
343         pc[16] = m17 = p2*x16 + p7*x17  + p12*x18 + p17*x19 + p22*x20;
344         pc[17] = m18 = p3*x16 + p8*x17  + p13*x18 + p18*x19 + p23*x20;
345         pc[18] = m19 = p4*x16 + p9*x17  + p14*x18 + p19*x19 + p24*x20;
346         pc[19] = m20 = p5*x16 + p10*x17 + p15*x18 + p20*x19 + p25*x20;
347 
348         pc[20] = m21 = p1*x21 + p6*x22  + p11*x23 + p16*x24 + p21*x25;
349         pc[21] = m22 = p2*x21 + p7*x22  + p12*x23 + p17*x24 + p22*x25;
350         pc[22] = m23 = p3*x21 + p8*x22  + p13*x23 + p18*x24 + p23*x25;
351         pc[23] = m24 = p4*x21 + p9*x22  + p14*x23 + p19*x24 + p24*x25;
352         pc[24] = m25 = p5*x21 + p10*x22 + p15*x23 + p20*x24 + p25*x25;
353 
354         nz = bi[row+1] - diag_offset[row] - 1;
355         pv += 25;
356         for (j=0; j<nz; j++) {
357           x1   = pv[0];  x2 = pv[1];   x3  = pv[2];  x4  = pv[3];
358           x5   = pv[4];  x6 = pv[5];   x7  = pv[6];  x8  = pv[7]; x9 = pv[8];
359           x10  = pv[9];  x11 = pv[10]; x12 = pv[11]; x13 = pv[12];
360           x14  = pv[13]; x15 = pv[14]; x16 = pv[15]; x17 = pv[16];
361           x18  = pv[17]; x19 = pv[18]; x20 = pv[19]; x21 = pv[20];
362           x22  = pv[21]; x23 = pv[22]; x24 = pv[23]; x25 = pv[24];
363           x    = rtmp + 25*pj[j];
364           x[0] -= m1*x1 + m6*x2  + m11*x3 + m16*x4 + m21*x5;
365           x[1] -= m2*x1 + m7*x2  + m12*x3 + m17*x4 + m22*x5;
366           x[2] -= m3*x1 + m8*x2  + m13*x3 + m18*x4 + m23*x5;
367           x[3] -= m4*x1 + m9*x2  + m14*x3 + m19*x4 + m24*x5;
368           x[4] -= m5*x1 + m10*x2 + m15*x3 + m20*x4 + m25*x5;
369 
370           x[5] -= m1*x6 + m6*x7  + m11*x8 + m16*x9 + m21*x10;
371           x[6] -= m2*x6 + m7*x7  + m12*x8 + m17*x9 + m22*x10;
372           x[7] -= m3*x6 + m8*x7  + m13*x8 + m18*x9 + m23*x10;
373           x[8] -= m4*x6 + m9*x7  + m14*x8 + m19*x9 + m24*x10;
374           x[9] -= m5*x6 + m10*x7 + m15*x8 + m20*x9 + m25*x10;
375 
376           x[10] -= m1*x11 + m6*x12  + m11*x13 + m16*x14 + m21*x15;
377           x[11] -= m2*x11 + m7*x12  + m12*x13 + m17*x14 + m22*x15;
378           x[12] -= m3*x11 + m8*x12  + m13*x13 + m18*x14 + m23*x15;
379           x[13] -= m4*x11 + m9*x12  + m14*x13 + m19*x14 + m24*x15;
380           x[14] -= m5*x11 + m10*x12 + m15*x13 + m20*x14 + m25*x15;
381 
382           x[15] -= m1*x16 + m6*x17  + m11*x18 + m16*x19 + m21*x20;
383           x[16] -= m2*x16 + m7*x17  + m12*x18 + m17*x19 + m22*x20;
384           x[17] -= m3*x16 + m8*x17  + m13*x18 + m18*x19 + m23*x20;
385           x[18] -= m4*x16 + m9*x17  + m14*x18 + m19*x19 + m24*x20;
386           x[19] -= m5*x16 + m10*x17 + m15*x18 + m20*x19 + m25*x20;
387 
388           x[20] -= m1*x21 + m6*x22  + m11*x23 + m16*x24 + m21*x25;
389           x[21] -= m2*x21 + m7*x22  + m12*x23 + m17*x24 + m22*x25;
390           x[22] -= m3*x21 + m8*x22  + m13*x23 + m18*x24 + m23*x25;
391           x[23] -= m4*x21 + m9*x22  + m14*x23 + m19*x24 + m24*x25;
392           x[24] -= m5*x21 + m10*x22 + m15*x23 + m20*x24 + m25*x25;
393 
394           pv   += 25;
395         }
396         PLogFlops(250*nz+225);
397       }
398       row = *ajtmp++;
399     }
400     /* finished row so stick it into b->a */
401     pv = ba + 25*bi[i];
402     pj = bj + bi[i];
403     nz = bi[i+1] - bi[i];
404     for ( j=0; j<nz; j++ ) {
405       x     = rtmp+25*pj[j];
406       pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3];
407       pv[4] = x[4]; pv[5] = x[5]; pv[6] = x[6]; pv[7] = x[7]; pv[8] = x[8];
408       pv[9] = x[9]; pv[10] = x[10]; pv[11] = x[11]; pv[12] = x[12];
409       pv[13] = x[13]; pv[14] = x[14]; pv[15] = x[15]; pv[16] = x[16];
410       pv[17] = x[17]; pv[18] = x[18]; pv[19] = x[19]; pv[20] = x[20];
411       pv[21] = x[21]; pv[22] = x[22]; pv[23] = x[23]; pv[24] = x[24];
412       pv   += 25;
413     }
414     /* invert diagonal block */
415     w = ba + 25*diag_offset[i];
416     ierr = Kernel_A_gets_inverse_A_5(w); CHKERRQ(ierr);
417   }
418 
419   PetscFree(rtmp);
420   ierr = ISRestoreIndices(isicol,&ic); CHKERRQ(ierr);
421   ierr = ISRestoreIndices(isrow,&r); CHKERRQ(ierr);
422   ierr = ISDestroy(isicol); CHKERRQ(ierr);
423   C->factor = FACTOR_LU;
424   C->assembled = PETSC_TRUE;
425   PLogFlops(1.3333*125*b->mbs); /* from inverting diagonal blocks */
426   PetscFunctionReturn(0);
427 }
428 
429 /* ------------------------------------------------------------*/
430 /*
431       Version for when blocks are 4 by 4
432 */
433 #undef __FUNC__
434 #define __FUNC__ "MatLUFactorNumeric_SeqBAIJ_4"
435 int MatLUFactorNumeric_SeqBAIJ_4(Mat A,Mat *B)
436 {
437   Mat             C = *B;
438   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data,*b = (Mat_SeqBAIJ *)C->data;
439   IS              iscol = b->col, isrow = b->row, isicol;
440   int             *r,*ic, ierr, i, j, n = a->mbs, *bi = b->i, *bj = b->j;
441   int             *ajtmpold, *ajtmp, nz, row;
442   int             *diag_offset = b->diag,idx,*ai=a->i,*aj=a->j;
443   register int    *pj;
444   register Scalar *pv,*v,*rtmp,*pc,*w,*x;
445   Scalar          p1,p2,p3,p4,m1,m2,m3,m4,m5,m6,m7,m8,m9,x1,x2,x3,x4;
446   Scalar          p5,p6,p7,p8,p9,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15,x16;
447   Scalar          p10,p11,p12,p13,p14,p15,p16,m10,m11,m12;
448   Scalar          m13,m14,m15,m16;
449   Scalar          *ba = b->a,*aa = a->a;
450 
451   PetscFunctionBegin;
452   ierr  = ISInvertPermutation(iscol,&isicol); CHKERRQ(ierr);
453   PLogObjectParent(*B,isicol);
454   ierr  = ISGetIndices(isrow,&r); CHKERRQ(ierr);
455   ierr  = ISGetIndices(isicol,&ic); CHKERRQ(ierr);
456   rtmp  = (Scalar *) PetscMalloc(16*(n+1)*sizeof(Scalar));CHKPTRQ(rtmp);
457 
458   for ( i=0; i<n; i++ ) {
459     nz    = bi[i+1] - bi[i];
460     ajtmp = bj + bi[i];
461     for  ( j=0; j<nz; j++ ) {
462       x = rtmp+16*ajtmp[j];
463       x[0]  = x[1]  = x[2]  = x[3]  = x[4]  = x[5]  = x[6] = x[7] = x[8] = x[9] = 0.0;
464       x[10] = x[11] = x[12] = x[13] = x[14] = x[15] = 0.0;
465     }
466     /* load in initial (unfactored row) */
467     idx      = r[i];
468     nz       = ai[idx+1] - ai[idx];
469     ajtmpold = aj + ai[idx];
470     v        = aa + 16*ai[idx];
471     for ( j=0; j<nz; j++ ) {
472       x    = rtmp+16*ic[ajtmpold[j]];
473       x[0]  = v[0];  x[1]  = v[1];  x[2]  = v[2];  x[3]  = v[3];
474       x[4]  = v[4];  x[5]  = v[5];  x[6]  = v[6];  x[7]  = v[7];  x[8]  = v[8];
475       x[9]  = v[9];  x[10] = v[10]; x[11] = v[11]; x[12] = v[12]; x[13] = v[13];
476       x[14] = v[14]; x[15] = v[15];
477       v    += 16;
478     }
479     row = *ajtmp++;
480     while (row < i) {
481       pc  = rtmp + 16*row;
482       p1  = pc[0];  p2  = pc[1];  p3  = pc[2];  p4  = pc[3];
483       p5  = pc[4];  p6  = pc[5];  p7  = pc[6];  p8  = pc[7];  p9  = pc[8];
484       p10 = pc[9];  p11 = pc[10]; p12 = pc[11]; p13 = pc[12]; p14 = pc[13];
485       p15 = pc[14]; p16 = pc[15];
486       if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 || p5 != 0.0 ||
487           p6 != 0.0 || p7 != 0.0 || p8 != 0.0 || p9 != 0.0 || p10 != 0.0 ||
488           p11 != 0.0 || p12 != 0.0 || p13 != 0.0 || p14 != 0.0 || p15 != 0.0
489           || p16 != 0.0) {
490         pv = ba + 16*diag_offset[row];
491         pj = bj + diag_offset[row] + 1;
492         x1  = pv[0];  x2  = pv[1];  x3  = pv[2];  x4  = pv[3];
493         x5  = pv[4];  x6  = pv[5];  x7  = pv[6];  x8  = pv[7];  x9  = pv[8];
494         x10 = pv[9];  x11 = pv[10]; x12 = pv[11]; x13 = pv[12]; x14 = pv[13];
495         x15 = pv[14]; x16 = pv[15];
496         pc[0] = m1 = p1*x1 + p5*x2  + p9*x3  + p13*x4;
497         pc[1] = m2 = p2*x1 + p6*x2  + p10*x3 + p14*x4;
498         pc[2] = m3 = p3*x1 + p7*x2  + p11*x3 + p15*x4;
499         pc[3] = m4 = p4*x1 + p8*x2  + p12*x3 + p16*x4;
500 
501         pc[4] = m5 = p1*x5 + p5*x6  + p9*x7  + p13*x8;
502         pc[5] = m6 = p2*x5 + p6*x6  + p10*x7 + p14*x8;
503         pc[6] = m7 = p3*x5 + p7*x6  + p11*x7 + p15*x8;
504         pc[7] = m8 = p4*x5 + p8*x6  + p12*x7 + p16*x8;
505 
506         pc[8]  = m9  = p1*x9 + p5*x10  + p9*x11  + p13*x12;
507         pc[9]  = m10 = p2*x9 + p6*x10  + p10*x11 + p14*x12;
508         pc[10] = m11 = p3*x9 + p7*x10  + p11*x11 + p15*x12;
509         pc[11] = m12 = p4*x9 + p8*x10  + p12*x11 + p16*x12;
510 
511         pc[12] = m13 = p1*x13 + p5*x14  + p9*x15  + p13*x16;
512         pc[13] = m14 = p2*x13 + p6*x14  + p10*x15 + p14*x16;
513         pc[14] = m15 = p3*x13 + p7*x14  + p11*x15 + p15*x16;
514         pc[15] = m16 = p4*x13 + p8*x14  + p12*x15 + p16*x16;
515 
516         nz = bi[row+1] - diag_offset[row] - 1;
517         pv += 16;
518         for (j=0; j<nz; j++) {
519           x1   = pv[0];  x2  = pv[1];   x3 = pv[2];  x4  = pv[3];
520           x5   = pv[4];  x6  = pv[5];   x7 = pv[6];  x8  = pv[7]; x9 = pv[8];
521           x10  = pv[9];  x11 = pv[10]; x12 = pv[11]; x13 = pv[12];
522           x14  = pv[13]; x15 = pv[14]; x16 = pv[15];
523           x    = rtmp + 16*pj[j];
524           x[0] -= m1*x1 + m5*x2  + m9*x3  + m13*x4;
525           x[1] -= m2*x1 + m6*x2  + m10*x3 + m14*x4;
526           x[2] -= m3*x1 + m7*x2  + m11*x3 + m15*x4;
527           x[3] -= m4*x1 + m8*x2  + m12*x3 + m16*x4;
528 
529           x[4] -= m1*x5 + m5*x6  + m9*x7  + m13*x8;
530           x[5] -= m2*x5 + m6*x6  + m10*x7 + m14*x8;
531           x[6] -= m3*x5 + m7*x6  + m11*x7 + m15*x8;
532           x[7] -= m4*x5 + m8*x6  + m12*x7 + m16*x8;
533 
534           x[8]  -= m1*x9 + m5*x10 + m9*x11  + m13*x12;
535           x[9]  -= m2*x9 + m6*x10 + m10*x11 + m14*x12;
536           x[10] -= m3*x9 + m7*x10 + m11*x11 + m15*x12;
537           x[11] -= m4*x9 + m8*x10 + m12*x11 + m16*x12;
538 
539           x[12] -= m1*x13 + m5*x14  + m9*x15  + m13*x16;
540           x[13] -= m2*x13 + m6*x14  + m10*x15 + m14*x16;
541           x[14] -= m3*x13 + m7*x14  + m11*x15 + m15*x16;
542           x[15] -= m4*x13 + m8*x14  + m12*x15 + m16*x16;
543 
544           pv   += 16;
545         }
546         PLogFlops(128*nz+112);
547       }
548       row = *ajtmp++;
549     }
550     /* finished row so stick it into b->a */
551     pv = ba + 16*bi[i];
552     pj = bj + bi[i];
553     nz = bi[i+1] - bi[i];
554     for ( j=0; j<nz; j++ ) {
555       x      = rtmp+16*pj[j];
556       pv[0]  = x[0];  pv[1]  = x[1];  pv[2]  = x[2];  pv[3]  = x[3];
557       pv[4]  = x[4];  pv[5]  = x[5];  pv[6]  = x[6];  pv[7]  = x[7]; pv[8] = x[8];
558       pv[9]  = x[9];  pv[10] = x[10]; pv[11] = x[11]; pv[12] = x[12];
559       pv[13] = x[13]; pv[14] = x[14]; pv[15] = x[15];
560       pv   += 16;
561     }
562     /* invert diagonal block */
563     w = ba + 16*diag_offset[i];
564     ierr = Kernel_A_gets_inverse_A_4(w); CHKERRQ(ierr);
565   }
566 
567   PetscFree(rtmp);
568   ierr = ISRestoreIndices(isicol,&ic); CHKERRQ(ierr);
569   ierr = ISRestoreIndices(isrow,&r); CHKERRQ(ierr);
570   ierr = ISDestroy(isicol); CHKERRQ(ierr);
571   C->factor = FACTOR_LU;
572   C->assembled = PETSC_TRUE;
573   PLogFlops(1.3333*64*b->mbs); /* from inverting diagonal blocks */
574   PetscFunctionReturn(0);
575 }
576 /*
577       Version for when blocks are 4 by 4 Using natural ordering
578 */
579 #undef __FUNC__
580 #define __FUNC__ "MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering"
581 int MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering(Mat A,Mat *B)
582 {
583   Mat             C = *B;
584   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data,*b = (Mat_SeqBAIJ *)C->data;
585   int             ierr, i, j, n = a->mbs, *bi = b->i, *bj = b->j;
586   int             *ajtmpold, *ajtmp, nz, row;
587   int             *diag_offset = b->diag,*ai=a->i,*aj=a->j;
588   register int    *pj;
589   register Scalar *pv,*v,*rtmp,*pc,*w,*x;
590   Scalar          p1,p2,p3,p4,m1,m2,m3,m4,m5,m6,m7,m8,m9,x1,x2,x3,x4;
591   Scalar          p5,p6,p7,p8,p9,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15,x16;
592   Scalar          p10,p11,p12,p13,p14,p15,p16,m10,m11,m12;
593   Scalar          m13,m14,m15,m16;
594   Scalar          *ba = b->a,*aa = a->a;
595 
596   PetscFunctionBegin;
597   rtmp  = (Scalar *) PetscMalloc(16*(n+1)*sizeof(Scalar));CHKPTRQ(rtmp);
598 
599   for ( i=0; i<n; i++ ) {
600     nz    = bi[i+1] - bi[i];
601     ajtmp = bj + bi[i];
602     for  ( j=0; j<nz; j++ ) {
603       x = rtmp+16*ajtmp[j];
604       x[0]  = x[1]  = x[2]  = x[3]  = x[4]  = x[5]  = x[6] = x[7] = x[8] = x[9] = 0.0;
605       x[10] = x[11] = x[12] = x[13] = x[14] = x[15] = 0.0;
606     }
607     /* load in initial (unfactored row) */
608     nz       = ai[i+1] - ai[i];
609     ajtmpold = aj + ai[i];
610     v        = aa + 16*ai[i];
611     for ( j=0; j<nz; j++ ) {
612       x    = rtmp+16*ajtmpold[j];
613       x[0]  = v[0];  x[1]  = v[1];  x[2]  = v[2];  x[3]  = v[3];
614       x[4]  = v[4];  x[5]  = v[5];  x[6]  = v[6];  x[7]  = v[7];  x[8]  = v[8];
615       x[9]  = v[9];  x[10] = v[10]; x[11] = v[11]; x[12] = v[12]; x[13] = v[13];
616       x[14] = v[14]; x[15] = v[15];
617       v    += 16;
618     }
619     row = *ajtmp++;
620     while (row < i) {
621       pc  = rtmp + 16*row;
622       p1  = pc[0];  p2  = pc[1];  p3  = pc[2];  p4  = pc[3];
623       p5  = pc[4];  p6  = pc[5];  p7  = pc[6];  p8  = pc[7];  p9  = pc[8];
624       p10 = pc[9];  p11 = pc[10]; p12 = pc[11]; p13 = pc[12]; p14 = pc[13];
625       p15 = pc[14]; p16 = pc[15];
626       if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 || p5 != 0.0 ||
627           p6 != 0.0 || p7 != 0.0 || p8 != 0.0 || p9 != 0.0 || p10 != 0.0 ||
628           p11 != 0.0 || p12 != 0.0 || p13 != 0.0 || p14 != 0.0 || p15 != 0.0
629           || p16 != 0.0) {
630         pv = ba + 16*diag_offset[row];
631         pj = bj + diag_offset[row] + 1;
632         x1  = pv[0];  x2  = pv[1];  x3  = pv[2];  x4  = pv[3];
633         x5  = pv[4];  x6  = pv[5];  x7  = pv[6];  x8  = pv[7];  x9  = pv[8];
634         x10 = pv[9];  x11 = pv[10]; x12 = pv[11]; x13 = pv[12]; x14 = pv[13];
635         x15 = pv[14]; x16 = pv[15];
636         pc[0] = m1 = p1*x1 + p5*x2  + p9*x3  + p13*x4;
637         pc[1] = m2 = p2*x1 + p6*x2  + p10*x3 + p14*x4;
638         pc[2] = m3 = p3*x1 + p7*x2  + p11*x3 + p15*x4;
639         pc[3] = m4 = p4*x1 + p8*x2  + p12*x3 + p16*x4;
640 
641         pc[4] = m5 = p1*x5 + p5*x6  + p9*x7  + p13*x8;
642         pc[5] = m6 = p2*x5 + p6*x6  + p10*x7 + p14*x8;
643         pc[6] = m7 = p3*x5 + p7*x6  + p11*x7 + p15*x8;
644         pc[7] = m8 = p4*x5 + p8*x6  + p12*x7 + p16*x8;
645 
646         pc[8]  = m9  = p1*x9 + p5*x10  + p9*x11  + p13*x12;
647         pc[9]  = m10 = p2*x9 + p6*x10  + p10*x11 + p14*x12;
648         pc[10] = m11 = p3*x9 + p7*x10  + p11*x11 + p15*x12;
649         pc[11] = m12 = p4*x9 + p8*x10  + p12*x11 + p16*x12;
650 
651         pc[12] = m13 = p1*x13 + p5*x14  + p9*x15  + p13*x16;
652         pc[13] = m14 = p2*x13 + p6*x14  + p10*x15 + p14*x16;
653         pc[14] = m15 = p3*x13 + p7*x14  + p11*x15 + p15*x16;
654         pc[15] = m16 = p4*x13 + p8*x14  + p12*x15 + p16*x16;
655 
656         nz = bi[row+1] - diag_offset[row] - 1;
657         pv += 16;
658         for (j=0; j<nz; j++) {
659           x1   = pv[0];  x2  = pv[1];   x3 = pv[2];  x4  = pv[3];
660           x5   = pv[4];  x6  = pv[5];   x7 = pv[6];  x8  = pv[7]; x9 = pv[8];
661           x10  = pv[9];  x11 = pv[10]; x12 = pv[11]; x13 = pv[12];
662           x14  = pv[13]; x15 = pv[14]; x16 = pv[15];
663           x    = rtmp + 16*pj[j];
664           x[0] -= m1*x1 + m5*x2  + m9*x3  + m13*x4;
665           x[1] -= m2*x1 + m6*x2  + m10*x3 + m14*x4;
666           x[2] -= m3*x1 + m7*x2  + m11*x3 + m15*x4;
667           x[3] -= m4*x1 + m8*x2  + m12*x3 + m16*x4;
668 
669           x[4] -= m1*x5 + m5*x6  + m9*x7  + m13*x8;
670           x[5] -= m2*x5 + m6*x6  + m10*x7 + m14*x8;
671           x[6] -= m3*x5 + m7*x6  + m11*x7 + m15*x8;
672           x[7] -= m4*x5 + m8*x6  + m12*x7 + m16*x8;
673 
674           x[8]  -= m1*x9 + m5*x10 + m9*x11  + m13*x12;
675           x[9]  -= m2*x9 + m6*x10 + m10*x11 + m14*x12;
676           x[10] -= m3*x9 + m7*x10 + m11*x11 + m15*x12;
677           x[11] -= m4*x9 + m8*x10 + m12*x11 + m16*x12;
678 
679           x[12] -= m1*x13 + m5*x14  + m9*x15  + m13*x16;
680           x[13] -= m2*x13 + m6*x14  + m10*x15 + m14*x16;
681           x[14] -= m3*x13 + m7*x14  + m11*x15 + m15*x16;
682           x[15] -= m4*x13 + m8*x14  + m12*x15 + m16*x16;
683 
684           pv   += 16;
685         }
686         PLogFlops(128*nz+112);
687       }
688       row = *ajtmp++;
689     }
690     /* finished row so stick it into b->a */
691     pv = ba + 16*bi[i];
692     pj = bj + bi[i];
693     nz = bi[i+1] - bi[i];
694     for ( j=0; j<nz; j++ ) {
695       x      = rtmp+16*pj[j];
696       pv[0]  = x[0];  pv[1]  = x[1];  pv[2]  = x[2];  pv[3]  = x[3];
697       pv[4]  = x[4];  pv[5]  = x[5];  pv[6]  = x[6];  pv[7]  = x[7]; pv[8] = x[8];
698       pv[9]  = x[9];  pv[10] = x[10]; pv[11] = x[11]; pv[12] = x[12];
699       pv[13] = x[13]; pv[14] = x[14]; pv[15] = x[15];
700       pv   += 16;
701     }
702     /* invert diagonal block */
703     w = ba + 16*diag_offset[i];
704     ierr = Kernel_A_gets_inverse_A_4(w); CHKERRQ(ierr);
705   }
706 
707   PetscFree(rtmp);
708   C->factor    = FACTOR_LU;
709   C->assembled = PETSC_TRUE;
710   PLogFlops(1.3333*64*b->mbs); /* from inverting diagonal blocks */
711   PetscFunctionReturn(0);
712 }
713 
714 /* ------------------------------------------------------------*/
715 /*
716       Version for when blocks are 3 by 3
717 */
718 #undef __FUNC__
719 #define __FUNC__ "MatLUFactorNumeric_SeqBAIJ_3"
720 int MatLUFactorNumeric_SeqBAIJ_3(Mat A,Mat *B)
721 {
722   Mat             C = *B;
723   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data,*b = (Mat_SeqBAIJ *)C->data;
724   IS              iscol = b->col, isrow = b->row, isicol;
725   int             *r,*ic, ierr, i, j, n = a->mbs, *bi = b->i, *bj = b->j;
726   int             *ajtmpold, *ajtmp, nz, row, *ai=a->i,*aj=a->j;
727   int             *diag_offset = b->diag,idx;
728   register int    *pj;
729   register Scalar *pv,*v,*rtmp,*pc,*w,*x;
730   Scalar          p1,p2,p3,p4,m1,m2,m3,m4,m5,m6,m7,m8,m9,x1,x2,x3,x4;
731   Scalar          p5,p6,p7,p8,p9,x5,x6,x7,x8,x9;
732   Scalar          *ba = b->a,*aa = a->a;
733 
734   PetscFunctionBegin;
735   ierr  = ISInvertPermutation(iscol,&isicol); CHKERRQ(ierr);
736   PLogObjectParent(*B,isicol);
737   ierr  = ISGetIndices(isrow,&r); CHKERRQ(ierr);
738   ierr  = ISGetIndices(isicol,&ic); CHKERRQ(ierr);
739   rtmp  = (Scalar *) PetscMalloc(9*(n+1)*sizeof(Scalar));CHKPTRQ(rtmp);
740 
741   for ( i=0; i<n; i++ ) {
742     nz    = bi[i+1] - bi[i];
743     ajtmp = bj + bi[i];
744     for  ( j=0; j<nz; j++ ) {
745       x = rtmp + 9*ajtmp[j];
746       x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = x[9] = 0.0;
747     }
748     /* load in initial (unfactored row) */
749     idx      = r[i];
750     nz       = ai[idx+1] - ai[idx];
751     ajtmpold = aj + ai[idx];
752     v        = aa + 9*ai[idx];
753     for ( j=0; j<nz; j++ ) {
754       x    = rtmp + 9*ic[ajtmpold[j]];
755       x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3];
756       x[4] = v[4]; x[5] = v[5]; x[6] = v[6]; x[7] = v[7]; x[8] = v[8];
757       v    += 9;
758     }
759     row = *ajtmp++;
760     while (row < i) {
761       pc = rtmp + 9*row;
762       p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3];
763       p5 = pc[4]; p6 = pc[5]; p7 = pc[6]; p8 = pc[7]; p9 = pc[8];
764       if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 || p5 != 0.0 ||
765           p6 != 0.0 || p7 != 0.0 || p8 != 0.0 || p9 != 0.0) {
766         pv = ba + 9*diag_offset[row];
767         pj = bj + diag_offset[row] + 1;
768         x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
769         x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8];
770         pc[0] = m1 = p1*x1 + p4*x2 + p7*x3;
771         pc[1] = m2 = p2*x1 + p5*x2 + p8*x3;
772         pc[2] = m3 = p3*x1 + p6*x2 + p9*x3;
773 
774         pc[3] = m4 = p1*x4 + p4*x5 + p7*x6;
775         pc[4] = m5 = p2*x4 + p5*x5 + p8*x6;
776         pc[5] = m6 = p3*x4 + p6*x5 + p9*x6;
777 
778         pc[6] = m7 = p1*x7 + p4*x8 + p7*x9;
779         pc[7] = m8 = p2*x7 + p5*x8 + p8*x9;
780         pc[8] = m9 = p3*x7 + p6*x8 + p9*x9;
781         nz = bi[row+1] - diag_offset[row] - 1;
782         pv += 9;
783         for (j=0; j<nz; j++) {
784           x1   = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
785           x5   = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8];
786           x    = rtmp + 9*pj[j];
787           x[0] -= m1*x1 + m4*x2 + m7*x3;
788           x[1] -= m2*x1 + m5*x2 + m8*x3;
789           x[2] -= m3*x1 + m6*x2 + m9*x3;
790 
791           x[3] -= m1*x4 + m4*x5 + m7*x6;
792           x[4] -= m2*x4 + m5*x5 + m8*x6;
793           x[5] -= m3*x4 + m6*x5 + m9*x6;
794 
795           x[6] -= m1*x7 + m4*x8 + m7*x9;
796           x[7] -= m2*x7 + m5*x8 + m8*x9;
797           x[8] -= m3*x7 + m6*x8 + m9*x9;
798           pv   += 9;
799         }
800         PLogFlops(54*nz+36);
801       }
802       row = *ajtmp++;
803     }
804     /* finished row so stick it into b->a */
805     pv = ba + 9*bi[i];
806     pj = bj + bi[i];
807     nz = bi[i+1] - bi[i];
808     for ( j=0; j<nz; j++ ) {
809       x     = rtmp + 9*pj[j];
810       pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3];
811       pv[4] = x[4]; pv[5] = x[5]; pv[6] = x[6]; pv[7] = x[7]; pv[8] = x[8];
812       pv   += 9;
813     }
814     /* invert diagonal block */
815     w = ba + 9*diag_offset[i];
816     ierr = Kernel_A_gets_inverse_A_3(w); CHKERRQ(ierr);
817   }
818 
819   PetscFree(rtmp);
820   ierr = ISRestoreIndices(isicol,&ic); CHKERRQ(ierr);
821   ierr = ISRestoreIndices(isrow,&r); CHKERRQ(ierr);
822   ierr = ISDestroy(isicol); CHKERRQ(ierr);
823   C->factor = FACTOR_LU;
824   C->assembled = PETSC_TRUE;
825   PLogFlops(1.3333*27*b->mbs); /* from inverting diagonal blocks */
826   PetscFunctionReturn(0);
827 }
828 
829 /* ------------------------------------------------------------*/
830 /*
831       Version for when blocks are 2 by 2
832 */
833 #undef __FUNC__
834 #define __FUNC__ "MatLUFactorNumeric_SeqBAIJ_2"
835 int MatLUFactorNumeric_SeqBAIJ_2(Mat A,Mat *B)
836 {
837   Mat             C = *B;
838   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data,*b = (Mat_SeqBAIJ *)C->data;
839   IS              iscol = b->col, isrow = b->row, isicol;
840   int             *r,*ic, ierr, i, j, n = a->mbs, *bi = b->i, *bj = b->j;
841   int             *ajtmpold, *ajtmp, nz, row, v_pivots[2];
842   int             *diag_offset=b->diag,bs = 2,idx,*ai=a->i,*aj=a->j;
843   register int    *pj;
844   register Scalar *pv,*v,*rtmp,m1,m2,m3,m4,*pc,*w,*x,x1,x2,x3,x4;
845   Scalar          p1,p2,p3,p4,v_work[2];
846   Scalar          *ba = b->a,*aa = a->a;
847 
848   PetscFunctionBegin;
849   ierr  = ISInvertPermutation(iscol,&isicol); CHKERRQ(ierr);
850   PLogObjectParent(*B,isicol);
851   ierr  = ISGetIndices(isrow,&r); CHKERRQ(ierr);
852   ierr  = ISGetIndices(isicol,&ic); CHKERRQ(ierr);
853   rtmp  = (Scalar *) PetscMalloc(4*(n+1)*sizeof(Scalar));CHKPTRQ(rtmp);
854 
855   for ( i=0; i<n; i++ ) {
856     nz    = bi[i+1] - bi[i];
857     ajtmp = bj + bi[i];
858     for  ( j=0; j<nz; j++ ) {
859       x = rtmp+4*ajtmp[j]; x[0] = x[1] = x[2] = x[3] = 0.0;
860     }
861     /* load in initial (unfactored row) */
862     idx      = r[i];
863     nz       = ai[idx+1] - ai[idx];
864     ajtmpold = aj + ai[idx];
865     v        = aa + 4*ai[idx];
866     for ( j=0; j<nz; j++ ) {
867       x    = rtmp+4*ic[ajtmpold[j]];
868       x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3];
869       v    += 4;
870     }
871     row = *ajtmp++;
872     while (row < i) {
873       pc = rtmp + 4*row;
874       p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3];
875       if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0) {
876         pv = ba + 4*diag_offset[row];
877         pj = bj + diag_offset[row] + 1;
878         x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
879         pc[0] = m1 = p1*x1 + p3*x2;
880         pc[1] = m2 = p2*x1 + p4*x2;
881         pc[2] = m3 = p1*x3 + p3*x4;
882         pc[3] = m4 = p2*x3 + p4*x4;
883         nz = bi[row+1] - diag_offset[row] - 1;
884         pv += 4;
885         for (j=0; j<nz; j++) {
886           x1   = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
887           x    = rtmp + 4*pj[j];
888           x[0] -= m1*x1 + m3*x2;
889           x[1] -= m2*x1 + m4*x2;
890           x[2] -= m1*x3 + m3*x4;
891           x[3] -= m2*x3 + m4*x4;
892           pv   += 4;
893         }
894         PLogFlops(16*nz+12);
895       }
896       row = *ajtmp++;
897     }
898     /* finished row so stick it into b->a */
899     pv = ba + 4*bi[i];
900     pj = bj + bi[i];
901     nz = bi[i+1] - bi[i];
902     for ( j=0; j<nz; j++ ) {
903       x     = rtmp+4*pj[j];
904       pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3];
905       pv   += 4;
906     }
907     /* invert diagonal block */
908     w = ba + 4*diag_offset[i];
909     Kernel_A_gets_inverse_A(bs,w,v_pivots,v_work);
910   }
911 
912   PetscFree(rtmp);
913   ierr = ISRestoreIndices(isicol,&ic); CHKERRQ(ierr);
914   ierr = ISRestoreIndices(isrow,&r); CHKERRQ(ierr);
915   ierr = ISDestroy(isicol); CHKERRQ(ierr);
916   C->factor = FACTOR_LU;
917   C->assembled = PETSC_TRUE;
918   PLogFlops(1.3333*8*b->mbs); /* from inverting diagonal blocks */
919   PetscFunctionReturn(0);
920 }
921 
922 /* ----------------------------------------------------------- */
923 /*
924      Version for when blocks are 1 by 1.
925 */
926 #undef __FUNC__
927 #define __FUNC__ "MatLUFactorNumeric_SeqBAIJ_1"
928 int MatLUFactorNumeric_SeqBAIJ_1(Mat A,Mat *B)
929 {
930   Mat             C = *B;
931   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data, *b = (Mat_SeqBAIJ *)C->data;
932   IS              iscol = b->col, isrow = b->row, isicol;
933   int             *r,*ic, ierr, i, j, n = a->mbs, *bi = b->i, *bj = b->j;
934   int             *ajtmpold, *ajtmp, nz, row,*ai = a->i,*aj = a->j;
935   int             *diag_offset = b->diag,diag;
936   register int    *pj;
937   register Scalar *pv,*v,*rtmp,multiplier,*pc;
938   Scalar          *ba = b->a,*aa = a->a;
939 
940   PetscFunctionBegin;
941   ierr  = ISInvertPermutation(iscol,&isicol); CHKERRQ(ierr);
942   PLogObjectParent(*B,isicol);
943   ierr  = ISGetIndices(isrow,&r); CHKERRQ(ierr);
944   ierr  = ISGetIndices(isicol,&ic); CHKERRQ(ierr);
945   rtmp  = (Scalar *) PetscMalloc((n+1)*sizeof(Scalar));CHKPTRQ(rtmp);
946 
947   for ( i=0; i<n; i++ ) {
948     nz    = bi[i+1] - bi[i];
949     ajtmp = bj + bi[i];
950     for  ( j=0; j<nz; j++ ) rtmp[ajtmp[j]] = 0.0;
951 
952     /* load in initial (unfactored row) */
953     nz       = ai[r[i]+1] - ai[r[i]];
954     ajtmpold = aj + ai[r[i]];
955     v        = aa + ai[r[i]];
956     for ( j=0; j<nz; j++ ) rtmp[ic[ajtmpold[j]]] =  v[j];
957 
958     row = *ajtmp++;
959     while (row < i) {
960       pc = rtmp + row;
961       if (*pc != 0.0) {
962         pv         = ba + diag_offset[row];
963         pj         = bj + diag_offset[row] + 1;
964         multiplier = *pc * *pv++;
965         *pc        = multiplier;
966         nz         = bi[row+1] - diag_offset[row] - 1;
967         for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j];
968         PLogFlops(1+2*nz);
969       }
970       row = *ajtmp++;
971     }
972     /* finished row so stick it into b->a */
973     pv = ba + bi[i];
974     pj = bj + bi[i];
975     nz = bi[i+1] - bi[i];
976     for ( j=0; j<nz; j++ ) {pv[j] = rtmp[pj[j]];}
977     diag = diag_offset[i] - bi[i];
978     /* check pivot entry for current row */
979     if (pv[diag] == 0.0) {
980       SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,0,"Zero pivot");
981     }
982     pv[diag] = 1.0/pv[diag];
983   }
984 
985   PetscFree(rtmp);
986   ierr = ISRestoreIndices(isicol,&ic); CHKERRQ(ierr);
987   ierr = ISRestoreIndices(isrow,&r); CHKERRQ(ierr);
988   ierr = ISDestroy(isicol); CHKERRQ(ierr);
989   C->factor    = FACTOR_LU;
990   C->assembled = PETSC_TRUE;
991   PLogFlops(b->n);
992   PetscFunctionReturn(0);
993 }
994 
995 /* ----------------------------------------------------------- */
996 #undef __FUNC__
997 #define __FUNC__ "MatLUFactor_SeqBAIJ"
998 int MatLUFactor_SeqBAIJ(Mat A,IS row,IS col,double f)
999 {
1000   Mat_SeqBAIJ *mat = (Mat_SeqBAIJ *) A->data;
1001   int         ierr;
1002   Mat         C;
1003 
1004   PetscFunctionBegin;
1005   ierr = MatLUFactorSymbolic(A,row,col,f,&C); CHKERRQ(ierr);
1006   ierr = MatLUFactorNumeric(A,&C); CHKERRQ(ierr);
1007 
1008   /* free all the data structures from mat */
1009   PetscFree(mat->a);
1010   if (!mat->singlemalloc) {PetscFree(mat->i); PetscFree(mat->j);}
1011   if (mat->diag) PetscFree(mat->diag);
1012   if (mat->ilen) PetscFree(mat->ilen);
1013   if (mat->imax) PetscFree(mat->imax);
1014   if (mat->solve_work) PetscFree(mat->solve_work);
1015   if (mat->mult_work) PetscFree(mat->mult_work);
1016   PetscFree(mat);
1017 
1018   PetscMemcpy(A,C,sizeof(struct _p_Mat));
1019   PetscHeaderDestroy(C);
1020   PetscFunctionReturn(0);
1021 }
1022 /* ----------------------------------------------------------- */
1023 #undef __FUNC__
1024 #define __FUNC__ "MatSolve_SeqBAIJ_N"
1025 int MatSolve_SeqBAIJ_N(Mat A,Vec bb,Vec xx)
1026 {
1027   Mat_SeqBAIJ     *a=(Mat_SeqBAIJ *)A->data;
1028   IS              iscol=a->col,isrow=a->row;
1029   int             *r,*c,ierr,i,n=a->mbs,*vi,*ai=a->i,*aj=a->j;
1030   int             nz,bs=a->bs,bs2=a->bs2,*rout,*cout;
1031   Scalar          *aa=a->a,*sum;
1032   register Scalar *x,*b,*lsum,*tmp,*v;
1033 
1034   PetscFunctionBegin;
1035   VecGetArray_Fast(bb,b);
1036   VecGetArray_Fast(xx,x);
1037   tmp  = a->solve_work;
1038 
1039   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
1040   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
1041 
1042   /* forward solve the lower triangular */
1043   PetscMemcpy(tmp,b + bs*(*r++), bs*sizeof(Scalar));
1044   for ( i=1; i<n; i++ ) {
1045     v   = aa + bs2*ai[i];
1046     vi  = aj + ai[i];
1047     nz  = a->diag[i] - ai[i];
1048     sum = tmp + bs*i;
1049     PetscMemcpy(sum,b+bs*(*r++),bs*sizeof(Scalar));
1050     while (nz--) {
1051       Kernel_v_gets_v_minus_A_times_w(bs,sum,v,tmp+bs*(*vi++));
1052       v += bs2;
1053     }
1054   }
1055   /* backward solve the upper triangular */
1056   lsum = a->solve_work + a->n;
1057   for ( i=n-1; i>=0; i-- ){
1058     v   = aa + bs2*(a->diag[i] + 1);
1059     vi  = aj + a->diag[i] + 1;
1060     nz  = ai[i+1] - a->diag[i] - 1;
1061     PetscMemcpy(lsum,tmp+i*bs,bs*sizeof(Scalar));
1062     while (nz--) {
1063       Kernel_v_gets_v_minus_A_times_w(bs,lsum,v,tmp+bs*(*vi++));
1064       v += bs2;
1065     }
1066     Kernel_w_gets_A_times_v(bs,lsum,aa+bs2*a->diag[i],tmp+i*bs);
1067     PetscMemcpy(x + bs*(*c--),tmp+i*bs,bs*sizeof(Scalar));
1068   }
1069 
1070   ierr = ISRestoreIndices(isrow,&rout); CHKERRQ(ierr);
1071   ierr = ISRestoreIndices(iscol,&cout); CHKERRQ(ierr);
1072   VecRestoreArray_Fast(bb,b);
1073   VecRestoreArray_Fast(xx,x);
1074   PLogFlops(2*(a->bs2)*(a->nz) - a->n);
1075   PetscFunctionReturn(0);
1076 }
1077 
1078 #undef __FUNC__
1079 #define __FUNC__ "MatSolve_SeqBAIJ_7"
1080 int MatSolve_SeqBAIJ_7(Mat A,Vec bb,Vec xx)
1081 {
1082   Mat_SeqBAIJ     *a=(Mat_SeqBAIJ *)A->data;
1083   IS              iscol=a->col,isrow=a->row;
1084   int             *r,*c,ierr,i,n=a->mbs,*vi,*ai=a->i,*aj=a->j,nz,idx,idt,idc,*rout,*cout;
1085   int             *diag = a->diag;
1086   Scalar          *aa=a->a,sum1,sum2,sum3,sum4,sum5,sum6,sum7,x1,x2,x3,x4,x5,x6,x7;
1087   register Scalar *x,*b,*tmp,*v;
1088 
1089   PetscFunctionBegin;
1090   VecGetArray_Fast(bb,b);
1091   VecGetArray_Fast(xx,x);
1092   tmp  = a->solve_work;
1093 
1094   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
1095   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
1096 
1097   /* forward solve the lower triangular */
1098   idx    = 7*(*r++);
1099   tmp[0] = b[idx];   tmp[1] = b[1+idx];
1100   tmp[2] = b[2+idx]; tmp[3] = b[3+idx]; tmp[4] = b[4+idx];
1101   tmp[5] = b[5+idx]; tmp[6] = b[6+idx];
1102 
1103   for ( i=1; i<n; i++ ) {
1104     v     = aa + 49*ai[i];
1105     vi    = aj + ai[i];
1106     nz    = diag[i] - ai[i];
1107     idx   = 7*(*r++);
1108     sum1  = b[idx];sum2 = b[1+idx];sum3 = b[2+idx];sum4 = b[3+idx];
1109     sum5  = b[4+idx];sum6 = b[5+idx];sum7 = b[6+idx];
1110     while (nz--) {
1111       idx   = 7*(*vi++);
1112       x1    = tmp[idx];  x2 = tmp[1+idx];x3 = tmp[2+idx];
1113       x4    = tmp[3+idx];x5 = tmp[4+idx];
1114       x6    = tmp[5+idx];x7 = tmp[6+idx];
1115       sum1 -= v[0]*x1 + v[7]*x2  + v[14]*x3 + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7;
1116       sum2 -= v[1]*x1 + v[8]*x2  + v[15]*x3 + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7;
1117       sum3 -= v[2]*x1 + v[9]*x2  + v[16]*x3 + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7;
1118       sum4 -= v[3]*x1 + v[10]*x2 + v[17]*x3 + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7;
1119       sum5 -= v[4]*x1 + v[11]*x2 + v[18]*x3 + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7;
1120       sum6 -= v[5]*x1 + v[12]*x2 + v[19]*x3 + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7;
1121       sum7 -= v[6]*x1 + v[13]*x2 + v[20]*x3 + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7;
1122       v += 49;
1123     }
1124     idx = 7*i;
1125     tmp[idx]   = sum1;tmp[1+idx] = sum2;
1126     tmp[2+idx] = sum3;tmp[3+idx] = sum4; tmp[4+idx] = sum5;
1127     tmp[5+idx] = sum6;tmp[6+idx] = sum7;
1128   }
1129   /* backward solve the upper triangular */
1130   for ( i=n-1; i>=0; i-- ){
1131     v    = aa + 49*diag[i] + 49;
1132     vi   = aj + diag[i] + 1;
1133     nz   = ai[i+1] - diag[i] - 1;
1134     idt  = 7*i;
1135     sum1 = tmp[idt];  sum2 = tmp[1+idt];
1136     sum3 = tmp[2+idt];sum4 = tmp[3+idt]; sum5 = tmp[4+idt];
1137     sum6 = tmp[5+idt];sum7 = tmp[6+idt];
1138     while (nz--) {
1139       idx   = 7*(*vi++);
1140       x1    = tmp[idx];   x2 = tmp[1+idx];
1141       x3    = tmp[2+idx]; x4 = tmp[3+idx]; x5 = tmp[4+idx];
1142       x6    = tmp[5+idx]; x7 = tmp[6+idx];
1143       sum1 -= v[0]*x1 + v[7]*x2  + v[14]*x3 + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7;
1144       sum2 -= v[1]*x1 + v[8]*x2  + v[15]*x3 + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7;
1145       sum3 -= v[2]*x1 + v[9]*x2  + v[16]*x3 + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7;
1146       sum4 -= v[3]*x1 + v[10]*x2 + v[17]*x3 + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7;
1147       sum5 -= v[4]*x1 + v[11]*x2 + v[18]*x3 + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7;
1148       sum6 -= v[5]*x1 + v[12]*x2 + v[19]*x3 + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7;
1149       sum7 -= v[6]*x1 + v[13]*x2 + v[20]*x3 + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7;
1150       v += 49;
1151     }
1152     idc = 7*(*c--);
1153     v   = aa + 49*diag[i];
1154     x[idc]   = tmp[idt]   = v[0]*sum1+v[7]*sum2+v[14]*sum3+
1155                                  v[21]*sum4+v[28]*sum5+v[35]*sum6+v[42]*sum7;
1156     x[1+idc] = tmp[1+idt] = v[1]*sum1+v[8]*sum2+v[15]*sum3+
1157                                  v[22]*sum4+v[29]*sum5+v[36]*sum6+v[43]*sum7;
1158     x[2+idc] = tmp[2+idt] = v[2]*sum1+v[9]*sum2+v[16]*sum3+
1159                                  v[23]*sum4+v[30]*sum5+v[37]*sum6+v[44]*sum7;
1160     x[3+idc] = tmp[3+idt] = v[3]*sum1+v[10]*sum2+v[17]*sum3+
1161                                  v[24]*sum4+v[31]*sum5+v[38]*sum6+v[45]*sum7;
1162     x[4+idc] = tmp[4+idt] = v[4]*sum1+v[11]*sum2+v[18]*sum3+
1163                                  v[25]*sum4+v[32]*sum5+v[39]*sum6+v[46]*sum7;
1164     x[5+idc] = tmp[5+idt] = v[5]*sum1+v[12]*sum2+v[19]*sum3+
1165                                  v[26]*sum4+v[33]*sum5+v[40]*sum6+v[47]*sum7;
1166     x[6+idc] = tmp[6+idt] = v[6]*sum1+v[13]*sum2+v[20]*sum3+
1167                                  v[27]*sum4+v[34]*sum5+v[41]*sum6+v[48]*sum7;
1168   }
1169 
1170   ierr = ISRestoreIndices(isrow,&rout); CHKERRQ(ierr);
1171   ierr = ISRestoreIndices(iscol,&cout); CHKERRQ(ierr);
1172   VecRestoreArray_Fast(bb,b);
1173   VecRestoreArray_Fast(xx,x);
1174   PLogFlops(2*49*(a->nz) - a->n);
1175   PetscFunctionReturn(0);
1176 }
1177 
1178 #undef __FUNC__
1179 #define __FUNC__ "MatSolve_SeqBAIJ_5"
1180 int MatSolve_SeqBAIJ_5(Mat A,Vec bb,Vec xx)
1181 {
1182   Mat_SeqBAIJ     *a=(Mat_SeqBAIJ *)A->data;
1183   IS              iscol=a->col,isrow=a->row;
1184   int             *r,*c,ierr,i,n=a->mbs,*vi,*ai=a->i,*aj=a->j,nz,idx,idt,idc,*rout,*cout;
1185   int             *diag = a->diag;
1186   Scalar          *aa=a->a,sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5;
1187   register Scalar *x,*b,*tmp,*v;
1188 
1189   PetscFunctionBegin;
1190   VecGetArray_Fast(bb,b);
1191   VecGetArray_Fast(xx,x);
1192   tmp  = a->solve_work;
1193 
1194   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
1195   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
1196 
1197   /* forward solve the lower triangular */
1198   idx    = 5*(*r++);
1199   tmp[0] = b[idx];   tmp[1] = b[1+idx];
1200   tmp[2] = b[2+idx]; tmp[3] = b[3+idx]; tmp[4] = b[4+idx];
1201   for ( i=1; i<n; i++ ) {
1202     v     = aa + 25*ai[i];
1203     vi    = aj + ai[i];
1204     nz    = diag[i] - ai[i];
1205     idx   = 5*(*r++);
1206     sum1  = b[idx];sum2 = b[1+idx];sum3 = b[2+idx];sum4 = b[3+idx];
1207     sum5  = b[4+idx];
1208     while (nz--) {
1209       idx   = 5*(*vi++);
1210       x1    = tmp[idx];  x2 = tmp[1+idx];x3 = tmp[2+idx];
1211       x4    = tmp[3+idx];x5 = tmp[4+idx];
1212       sum1 -= v[0]*x1 + v[5]*x2 + v[10]*x3 + v[15]*x4 + v[20]*x5;
1213       sum2 -= v[1]*x1 + v[6]*x2 + v[11]*x3 + v[16]*x4 + v[21]*x5;
1214       sum3 -= v[2]*x1 + v[7]*x2 + v[12]*x3 + v[17]*x4 + v[22]*x5;
1215       sum4 -= v[3]*x1 + v[8]*x2 + v[13]*x3 + v[18]*x4 + v[23]*x5;
1216       sum5 -= v[4]*x1 + v[9]*x2 + v[14]*x3 + v[19]*x4 + v[24]*x5;
1217       v += 25;
1218     }
1219     idx = 5*i;
1220     tmp[idx]   = sum1;tmp[1+idx] = sum2;
1221     tmp[2+idx] = sum3;tmp[3+idx] = sum4; tmp[4+idx] = sum5;
1222   }
1223   /* backward solve the upper triangular */
1224   for ( i=n-1; i>=0; i-- ){
1225     v    = aa + 25*diag[i] + 25;
1226     vi   = aj + diag[i] + 1;
1227     nz   = ai[i+1] - diag[i] - 1;
1228     idt  = 5*i;
1229     sum1 = tmp[idt];  sum2 = tmp[1+idt];
1230     sum3 = tmp[2+idt];sum4 = tmp[3+idt]; sum5 = tmp[4+idt];
1231     while (nz--) {
1232       idx   = 5*(*vi++);
1233       x1    = tmp[idx];   x2 = tmp[1+idx];
1234       x3    = tmp[2+idx]; x4 = tmp[3+idx]; x5 = tmp[4+idx];
1235       sum1 -= v[0]*x1 + v[5]*x2 + v[10]*x3 + v[15]*x4 + v[20]*x5;
1236       sum2 -= v[1]*x1 + v[6]*x2 + v[11]*x3 + v[16]*x4 + v[21]*x5;
1237       sum3 -= v[2]*x1 + v[7]*x2 + v[12]*x3 + v[17]*x4 + v[22]*x5;
1238       sum4 -= v[3]*x1 + v[8]*x2 + v[13]*x3 + v[18]*x4 + v[23]*x5;
1239       sum5 -= v[4]*x1 + v[9]*x2 + v[14]*x3 + v[19]*x4 + v[24]*x5;
1240       v += 25;
1241     }
1242     idc = 5*(*c--);
1243     v   = aa + 25*diag[i];
1244     x[idc]   = tmp[idt]   = v[0]*sum1+v[5]*sum2+v[10]*sum3+
1245                                  v[15]*sum4+v[20]*sum5;
1246     x[1+idc] = tmp[1+idt] = v[1]*sum1+v[6]*sum2+v[11]*sum3+
1247                                  v[16]*sum4+v[21]*sum5;
1248     x[2+idc] = tmp[2+idt] = v[2]*sum1+v[7]*sum2+v[12]*sum3+
1249                                  v[17]*sum4+v[22]*sum5;
1250     x[3+idc] = tmp[3+idt] = v[3]*sum1+v[8]*sum2+v[13]*sum3+
1251                                  v[18]*sum4+v[23]*sum5;
1252     x[4+idc] = tmp[4+idt] = v[4]*sum1+v[9]*sum2+v[14]*sum3+
1253                                  v[19]*sum4+v[24]*sum5;
1254   }
1255 
1256   ierr = ISRestoreIndices(isrow,&rout); CHKERRQ(ierr);
1257   ierr = ISRestoreIndices(iscol,&cout); CHKERRQ(ierr);
1258   VecRestoreArray_Fast(bb,b);
1259   VecRestoreArray_Fast(xx,x);
1260   PLogFlops(2*25*(a->nz) - a->n);
1261   PetscFunctionReturn(0);
1262 }
1263 
1264 #undef __FUNC__
1265 #define __FUNC__ "MatSolve_SeqBAIJ_4"
1266 int MatSolve_SeqBAIJ_4(Mat A,Vec bb,Vec xx)
1267 {
1268   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *)A->data;
1269   IS              iscol=a->col,isrow=a->row;
1270   int             *r,*c,ierr,i,n=a->mbs,*vi,*ai=a->i,*aj=a->j,nz,idx,idt,idc,*rout,*cout;
1271   int             *diag = a->diag;
1272   Scalar          *aa=a->a,sum1,sum2,sum3,sum4,x1,x2,x3,x4;
1273   register Scalar *x,*b,*tmp,*v;
1274 
1275   PetscFunctionBegin;
1276   VecGetArray_Fast(bb,b);
1277   VecGetArray_Fast(xx,x);
1278   tmp  = a->solve_work;
1279 
1280   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
1281   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
1282 
1283   /* forward solve the lower triangular */
1284   idx    = 4*(*r++);
1285   tmp[0] = b[idx];   tmp[1] = b[1+idx];
1286   tmp[2] = b[2+idx]; tmp[3] = b[3+idx];
1287   for ( i=1; i<n; i++ ) {
1288     v     = aa + 16*ai[i];
1289     vi    = aj + ai[i];
1290     nz    = diag[i] - ai[i];
1291     idx   = 4*(*r++);
1292     sum1  = b[idx];sum2 = b[1+idx];sum3 = b[2+idx];sum4 = b[3+idx];
1293     while (nz--) {
1294       idx   = 4*(*vi++);
1295       x1    = tmp[idx];x2 = tmp[1+idx];x3 = tmp[2+idx];x4 = tmp[3+idx];
1296       sum1 -= v[0]*x1 + v[4]*x2 + v[8]*x3  + v[12]*x4;
1297       sum2 -= v[1]*x1 + v[5]*x2 + v[9]*x3  + v[13]*x4;
1298       sum3 -= v[2]*x1 + v[6]*x2 + v[10]*x3 + v[14]*x4;
1299       sum4 -= v[3]*x1 + v[7]*x2 + v[11]*x3 + v[15]*x4;
1300       v    += 16;
1301     }
1302     idx        = 4*i;
1303     tmp[idx]   = sum1;tmp[1+idx] = sum2;
1304     tmp[2+idx] = sum3;tmp[3+idx] = sum4;
1305   }
1306   /* backward solve the upper triangular */
1307   for ( i=n-1; i>=0; i-- ){
1308     v    = aa + 16*diag[i] + 16;
1309     vi   = aj + diag[i] + 1;
1310     nz   = ai[i+1] - diag[i] - 1;
1311     idt  = 4*i;
1312     sum1 = tmp[idt];  sum2 = tmp[1+idt];
1313     sum3 = tmp[2+idt];sum4 = tmp[3+idt];
1314     while (nz--) {
1315       idx   = 4*(*vi++);
1316       x1    = tmp[idx];   x2 = tmp[1+idx];
1317       x3    = tmp[2+idx]; x4 = tmp[3+idx];
1318       sum1 -= v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
1319       sum2 -= v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
1320       sum3 -= v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
1321       sum4 -= v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
1322       v += 16;
1323     }
1324     idc      = 4*(*c--);
1325     v        = aa + 16*diag[i];
1326     x[idc]   = tmp[idt]   = v[0]*sum1+v[4]*sum2+v[8]*sum3+v[12]*sum4;
1327     x[1+idc] = tmp[1+idt] = v[1]*sum1+v[5]*sum2+v[9]*sum3+v[13]*sum4;
1328     x[2+idc] = tmp[2+idt] = v[2]*sum1+v[6]*sum2+v[10]*sum3+v[14]*sum4;
1329     x[3+idc] = tmp[3+idt] = v[3]*sum1+v[7]*sum2+v[11]*sum3+v[15]*sum4;
1330   }
1331 
1332   ierr = ISRestoreIndices(isrow,&rout); CHKERRQ(ierr);
1333   ierr = ISRestoreIndices(iscol,&cout); CHKERRQ(ierr);
1334   VecRestoreArray_Fast(bb,b);
1335   VecRestoreArray_Fast(xx,x);
1336   PLogFlops(2*16*(a->nz) - a->n);
1337   PetscFunctionReturn(0);
1338 }
1339 
1340 /*
1341       Special case where the matrix was ILU(0) factored in the natural
1342    ordering. This eliminates the need for the column and row permutation.
1343 */
1344 #undef __FUNC__
1345 #define __FUNC__ "MatSolve_SeqBAIJ_4_NaturalOrdering"
1346 int MatSolve_SeqBAIJ_4_NaturalOrdering(Mat A,Vec bb,Vec xx)
1347 {
1348   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *)A->data;
1349   int             i,n=a->mbs,*vi,*ai=a->i,*aj=a->j,nz,idx,idt;
1350   int             *diag = a->diag,jdx;
1351   Scalar          *aa=a->a,sum1,sum2,sum3,sum4,x1,x2,x3,x4;
1352   register Scalar *x,*b,*v;
1353 
1354   PetscFunctionBegin;
1355   VecGetArray_Fast(bb,b);
1356   VecGetArray_Fast(xx,x);
1357 
1358   /* forward solve the lower triangular */
1359   idx    = 0;
1360   x[0] = b[idx]; x[1] = b[1+idx]; x[2] = b[2+idx]; x[3] = b[3+idx];
1361   for ( i=1; i<n; i++ ) {
1362     v     =  aa + 16*ai[i];
1363     vi    =  aj + ai[i];
1364     nz    =  diag[i] - ai[i];
1365     idx   =  4*i;
1366     sum1  =  b[idx];sum2 = b[1+idx];sum3 = b[2+idx];sum4 = b[3+idx];
1367     while (nz--) {
1368       jdx   = 4*(*vi++);
1369       x1    = x[jdx];x2 = x[1+jdx];x3 = x[2+jdx];x4 = x[3+jdx];
1370       sum1 -= v[0]*x1 + v[4]*x2 + v[8]*x3  + v[12]*x4;
1371       sum2 -= v[1]*x1 + v[5]*x2 + v[9]*x3  + v[13]*x4;
1372       sum3 -= v[2]*x1 + v[6]*x2 + v[10]*x3 + v[14]*x4;
1373       sum4 -= v[3]*x1 + v[7]*x2 + v[11]*x3 + v[15]*x4;
1374       v    += 16;
1375     }
1376     x[idx]   = sum1;
1377     x[1+idx] = sum2;
1378     x[2+idx] = sum3;
1379     x[3+idx] = sum4;
1380   }
1381   /* backward solve the upper triangular */
1382   for ( i=n-1; i>=0; i-- ){
1383     v    = aa + 16*diag[i] + 16;
1384     vi   = aj + diag[i] + 1;
1385     nz   = ai[i+1] - diag[i] - 1;
1386     idt  = 4*i;
1387     sum1 = x[idt];  sum2 = x[1+idt];
1388     sum3 = x[2+idt];sum4 = x[3+idt];
1389     while (nz--) {
1390       idx   = 4*(*vi++);
1391       x1    = x[idx];   x2 = x[1+idx];x3    = x[2+idx]; x4 = x[3+idx];
1392       sum1 -= v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
1393       sum2 -= v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
1394       sum3 -= v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
1395       sum4 -= v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
1396       v += 16;
1397     }
1398     v        = aa + 16*diag[i];
1399     x[idt]   = v[0]*sum1 + v[4]*sum2 + v[8]*sum3  + v[12]*sum4;
1400     x[1+idt] = v[1]*sum1 + v[5]*sum2 + v[9]*sum3  + v[13]*sum4;
1401     x[2+idt] = v[2]*sum1 + v[6]*sum2 + v[10]*sum3 + v[14]*sum4;
1402     x[3+idt] = v[3]*sum1 + v[7]*sum2 + v[11]*sum3 + v[15]*sum4;
1403   }
1404 
1405   VecRestoreArray_Fast(bb,b);
1406   VecRestoreArray_Fast(xx,x);
1407   PLogFlops(2*16*(a->nz) - a->n);
1408   PetscFunctionReturn(0);
1409 }
1410 
1411 
1412 #undef __FUNC__
1413 #define __FUNC__ "MatSolve_SeqBAIJ_3"
1414 int MatSolve_SeqBAIJ_3(Mat A,Vec bb,Vec xx)
1415 {
1416   Mat_SeqBAIJ     *a=(Mat_SeqBAIJ *)A->data;
1417   IS              iscol=a->col,isrow=a->row;
1418   int             *r,*c,ierr,i,n=a->mbs,*vi,*ai=a->i,*aj=a->j,nz,idx,idt,idc,*rout,*cout;
1419   int             *diag = a->diag;
1420   Scalar          *aa=a->a,sum1,sum2,sum3,x1,x2,x3;
1421   register Scalar *x,*b,*tmp,*v;
1422 
1423   PetscFunctionBegin;
1424   VecGetArray_Fast(bb,b);
1425   VecGetArray_Fast(xx,x);
1426   tmp  = a->solve_work;
1427 
1428   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
1429   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
1430 
1431   /* forward solve the lower triangular */
1432   idx    = 3*(*r++);
1433   tmp[0] = b[idx]; tmp[1] = b[1+idx]; tmp[2] = b[2+idx];
1434   for ( i=1; i<n; i++ ) {
1435     v     = aa + 9*ai[i];
1436     vi    = aj + ai[i];
1437     nz    = diag[i] - ai[i];
1438     idx   = 3*(*r++);
1439     sum1  = b[idx]; sum2 = b[1+idx]; sum3 = b[2+idx];
1440     while (nz--) {
1441       idx   = 3*(*vi++);
1442       x1    = tmp[idx]; x2 = tmp[1+idx]; x3 = tmp[2+idx];
1443       sum1 -= v[0]*x1 + v[3]*x2 + v[6]*x3;
1444       sum2 -= v[1]*x1 + v[4]*x2 + v[7]*x3;
1445       sum3 -= v[2]*x1 + v[5]*x2 + v[8]*x3;
1446       v += 9;
1447     }
1448     idx = 3*i;
1449     tmp[idx] = sum1; tmp[1+idx] = sum2; tmp[2+idx] = sum3;
1450   }
1451   /* backward solve the upper triangular */
1452   for ( i=n-1; i>=0; i-- ){
1453     v    = aa + 9*diag[i] + 9;
1454     vi   = aj + diag[i] + 1;
1455     nz   = ai[i+1] - diag[i] - 1;
1456     idt  = 3*i;
1457     sum1 = tmp[idt]; sum2 = tmp[1+idt]; sum3 = tmp[2+idt];
1458     while (nz--) {
1459       idx   = 3*(*vi++);
1460       x1    = tmp[idx]; x2 = tmp[1+idx]; x3 = tmp[2+idx];
1461       sum1 -= v[0]*x1 + v[3]*x2 + v[6]*x3;
1462       sum2 -= v[1]*x1 + v[4]*x2 + v[7]*x3;
1463       sum3 -= v[2]*x1 + v[5]*x2 + v[8]*x3;
1464       v += 9;
1465     }
1466     idc = 3*(*c--);
1467     v   = aa + 9*diag[i];
1468     x[idc]   = tmp[idt]   = v[0]*sum1 + v[3]*sum2 + v[6]*sum3;
1469     x[1+idc] = tmp[1+idt] = v[1]*sum1 + v[4]*sum2 + v[7]*sum3;
1470     x[2+idc] = tmp[2+idt] = v[2]*sum1 + v[5]*sum2 + v[8]*sum3;
1471   }
1472   ierr = ISRestoreIndices(isrow,&rout); CHKERRQ(ierr);
1473   ierr = ISRestoreIndices(iscol,&cout); CHKERRQ(ierr);
1474   VecRestoreArray_Fast(bb,b);
1475   VecRestoreArray_Fast(xx,x);
1476   PLogFlops(2*9*(a->nz) - a->n);
1477   PetscFunctionReturn(0);
1478 }
1479 
1480 #undef __FUNC__
1481 #define __FUNC__ "MatSolve_SeqBAIJ_2"
1482 int MatSolve_SeqBAIJ_2(Mat A,Vec bb,Vec xx)
1483 {
1484   Mat_SeqBAIJ     *a=(Mat_SeqBAIJ *)A->data;
1485   IS              iscol=a->col,isrow=a->row;
1486   int             *r,*c,ierr,i,n=a->mbs,*vi,*ai=a->i,*aj=a->j,nz,idx,idt,idc,*rout,*cout;
1487   int             *diag = a->diag;
1488   Scalar          *aa=a->a,sum1,sum2,x1,x2;
1489   register Scalar *x,*b,*tmp,*v;
1490 
1491   PetscFunctionBegin;
1492   VecGetArray_Fast(bb,b);
1493   VecGetArray_Fast(xx,x);
1494   tmp  = a->solve_work;
1495 
1496   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
1497   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
1498 
1499   /* forward solve the lower triangular */
1500   idx    = 2*(*r++);
1501   tmp[0] = b[idx]; tmp[1] = b[1+idx];
1502   for ( i=1; i<n; i++ ) {
1503     v     = aa + 4*ai[i];
1504     vi    = aj + ai[i];
1505     nz    = diag[i] - ai[i];
1506     idx   = 2*(*r++);
1507     sum1  = b[idx]; sum2 = b[1+idx];
1508     while (nz--) {
1509       idx   = 2*(*vi++);
1510       x1    = tmp[idx]; x2 = tmp[1+idx];
1511       sum1 -= v[0]*x1 + v[2]*x2;
1512       sum2 -= v[1]*x1 + v[3]*x2;
1513       v += 4;
1514     }
1515     idx = 2*i;
1516     tmp[idx] = sum1; tmp[1+idx] = sum2;
1517   }
1518   /* backward solve the upper triangular */
1519   for ( i=n-1; i>=0; i-- ){
1520     v    = aa + 4*diag[i] + 4;
1521     vi   = aj + diag[i] + 1;
1522     nz   = ai[i+1] - diag[i] - 1;
1523     idt  = 2*i;
1524     sum1 = tmp[idt]; sum2 = tmp[1+idt];
1525     while (nz--) {
1526       idx   = 2*(*vi++);
1527       x1    = tmp[idx]; x2 = tmp[1+idx];
1528       sum1 -= v[0]*x1 + v[2]*x2;
1529       sum2 -= v[1]*x1 + v[3]*x2;
1530       v += 4;
1531     }
1532     idc = 2*(*c--);
1533     v   = aa + 4*diag[i];
1534     x[idc]   = tmp[idt]   = v[0]*sum1 + v[2]*sum2;
1535     x[1+idc] = tmp[1+idt] = v[1]*sum1 + v[3]*sum2;
1536   }
1537   ierr = ISRestoreIndices(isrow,&rout); CHKERRQ(ierr);
1538   ierr = ISRestoreIndices(iscol,&cout); CHKERRQ(ierr);
1539   VecRestoreArray_Fast(bb,b);
1540   VecRestoreArray_Fast(xx,x);
1541   PLogFlops(2*4*(a->nz) - a->n);
1542   PetscFunctionReturn(0);
1543 }
1544 
1545 
1546 #undef __FUNC__
1547 #define __FUNC__ "MatSolve_SeqBAIJ_1"
1548 int MatSolve_SeqBAIJ_1(Mat A,Vec bb,Vec xx)
1549 {
1550   Mat_SeqBAIJ     *a=(Mat_SeqBAIJ *)A->data;
1551   IS              iscol=a->col,isrow=a->row;
1552   int             *r,*c,ierr,i,n=a->mbs,*vi,*ai=a->i,*aj=a->j,nz,*rout,*cout;
1553   int             *diag = a->diag;
1554   Scalar          *aa=a->a,sum1;
1555   register Scalar *x,*b,*tmp,*v;
1556 
1557   PetscFunctionBegin;
1558   if (!n) PetscFunctionReturn(0);
1559 
1560   VecGetArray_Fast(bb,b);
1561   VecGetArray_Fast(xx,x);
1562   tmp  = a->solve_work;
1563 
1564   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
1565   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1);
1566 
1567   /* forward solve the lower triangular */
1568   tmp[0] = b[*r++];
1569   for ( i=1; i<n; i++ ) {
1570     v     = aa + ai[i];
1571     vi    = aj + ai[i];
1572     nz    = diag[i] - ai[i];
1573     sum1  = b[*r++];
1574     while (nz--) {
1575       sum1 -= (*v++)*tmp[*vi++];
1576     }
1577     tmp[i] = sum1;
1578   }
1579   /* backward solve the upper triangular */
1580   for ( i=n-1; i>=0; i-- ){
1581     v    = aa + diag[i] + 1;
1582     vi   = aj + diag[i] + 1;
1583     nz   = ai[i+1] - diag[i] - 1;
1584     sum1 = tmp[i];
1585     while (nz--) {
1586       sum1 -= (*v++)*tmp[*vi++];
1587     }
1588     x[*c--] = tmp[i] = aa[diag[i]]*sum1;
1589   }
1590 
1591   ierr = ISRestoreIndices(isrow,&rout); CHKERRQ(ierr);
1592   ierr = ISRestoreIndices(iscol,&cout); CHKERRQ(ierr);
1593   VecRestoreArray_Fast(bb,b);
1594   VecRestoreArray_Fast(xx,x);
1595   PLogFlops(2*1*(a->nz) - a->n);
1596   PetscFunctionReturn(0);
1597 }
1598 
1599 extern int MatSolve_SeqBAIJ_4_NaturalOrdering(Mat,Vec,Vec);
1600 /* ----------------------------------------------------------------*/
1601 /*
1602      This code is virtually identical to MatILUFactorSymbolic_SeqAIJ
1603    except that the data structure of Mat_SeqAIJ is slightly different.
1604    Not a good example of code reuse.
1605 */
1606 #undef __FUNC__
1607 #define __FUNC__ "MatILUFactorSymbolic_SeqBAIJ"
1608 int MatILUFactorSymbolic_SeqBAIJ(Mat A,IS isrow,IS iscol,double f,int levels,
1609                                  Mat *fact)
1610 {
1611   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data, *b;
1612   IS          isicol;
1613   int         *r,*ic, ierr, prow, n = a->mbs, *ai = a->i, *aj = a->j;
1614   int         *ainew,*ajnew, jmax,*fill, *xi, nz, *im,*ajfill,*flev;
1615   int         *dloc, idx, row,m,fm, nzf, nzi,len,  realloc = 0;
1616   int         incrlev,nnz,i,bs = a->bs,bs2 = a->bs2;
1617   PetscTruth  col_identity, row_identity;
1618 
1619   PetscFunctionBegin;
1620   /* special case that simply copies fill pattern */
1621   PetscValidHeaderSpecific(isrow,IS_COOKIE);
1622   PetscValidHeaderSpecific(iscol,IS_COOKIE);
1623   ISIdentity(isrow,&row_identity); ISIdentity(iscol,&col_identity);
1624   if (levels == 0 && row_identity && col_identity) {
1625     ierr = MatConvertSameType_SeqBAIJ(A,fact,DO_NOT_COPY_VALUES); CHKERRQ(ierr);
1626     (*fact)->factor = FACTOR_LU;
1627     b               = (Mat_SeqBAIJ *) (*fact)->data;
1628     if (!b->diag) {
1629       ierr = MatMarkDiag_SeqBAIJ(*fact); CHKERRQ(ierr);
1630     }
1631     b->row        = isrow;
1632     b->col        = iscol;
1633     b->solve_work = (Scalar *) PetscMalloc((b->m+1+b->bs)*sizeof(Scalar));CHKPTRQ(b->solve_work);
1634     /*
1635         Blocksize 4 has a special faster solver for ILU(0) factorization
1636         with natural ordering
1637     */
1638     if (b->bs == 4) {
1639       (*fact)->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering;
1640       (*fact)->ops.solve           = MatSolve_SeqBAIJ_4_NaturalOrdering;
1641     }
1642     PetscFunctionReturn(0);
1643   }
1644 
1645   ierr = ISInvertPermutation(iscol,&isicol); CHKERRQ(ierr);
1646   ierr = ISGetIndices(isrow,&r); CHKERRQ(ierr);
1647   ierr = ISGetIndices(isicol,&ic); CHKERRQ(ierr);
1648 
1649   /* get new row pointers */
1650   ainew = (int *) PetscMalloc( (n+1)*sizeof(int) ); CHKPTRQ(ainew);
1651   ainew[0] = 0;
1652   /* don't know how many column pointers are needed so estimate */
1653   jmax = (int) (f*ai[n] + 1);
1654   ajnew = (int *) PetscMalloc( (jmax)*sizeof(int) ); CHKPTRQ(ajnew);
1655   /* ajfill is level of fill for each fill entry */
1656   ajfill = (int *) PetscMalloc( (jmax)*sizeof(int) ); CHKPTRQ(ajfill);
1657   /* fill is a linked list of nonzeros in active row */
1658   fill = (int *) PetscMalloc( (n+1)*sizeof(int)); CHKPTRQ(fill);
1659   /* im is level for each filled value */
1660   im = (int *) PetscMalloc( (n+1)*sizeof(int)); CHKPTRQ(im);
1661   /* dloc is location of diagonal in factor */
1662   dloc = (int *) PetscMalloc( (n+1)*sizeof(int)); CHKPTRQ(dloc);
1663   dloc[0]  = 0;
1664   for ( prow=0; prow<n; prow++ ) {
1665     /* first copy previous fill into linked list */
1666     nzf     = nz  = ai[r[prow]+1] - ai[r[prow]];
1667     if (!nz) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,1,"Empty row in matrix");
1668     xi      = aj + ai[r[prow]];
1669     fill[n] = n;
1670     while (nz--) {
1671       fm  = n;
1672       idx = ic[*xi++];
1673       do {
1674         m  = fm;
1675         fm = fill[m];
1676       } while (fm < idx);
1677       fill[m]   = idx;
1678       fill[idx] = fm;
1679       im[idx]   = 0;
1680     }
1681     nzi = 0;
1682     row = fill[n];
1683     while ( row < prow ) {
1684       incrlev = im[row] + 1;
1685       nz      = dloc[row];
1686       xi      = ajnew  + ainew[row] + nz;
1687       flev    = ajfill + ainew[row] + nz + 1;
1688       nnz     = ainew[row+1] - ainew[row] - nz - 1;
1689       if (*xi++ != row) {
1690         SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,0,"Zero pivot: try running with -pc_ilu_nonzeros_along_diagonal");
1691       }
1692       fm      = row;
1693       while (nnz-- > 0) {
1694         idx = *xi++;
1695         if (*flev + incrlev > levels) {
1696           flev++;
1697           continue;
1698         }
1699         do {
1700           m  = fm;
1701           fm = fill[m];
1702         } while (fm < idx);
1703         if (fm != idx) {
1704           im[idx]   = *flev + incrlev;
1705           fill[m]   = idx;
1706           fill[idx] = fm;
1707           fm        = idx;
1708           nzf++;
1709         }
1710         else {
1711           if (im[idx] > *flev + incrlev) im[idx] = *flev+incrlev;
1712         }
1713         flev++;
1714       }
1715       row = fill[row];
1716       nzi++;
1717     }
1718     /* copy new filled row into permanent storage */
1719     ainew[prow+1] = ainew[prow] + nzf;
1720     if (ainew[prow+1] > jmax) {
1721       /* allocate a longer ajnew */
1722       int maxadd;
1723       maxadd = (int) (((f*ai[n]+1)*(n-prow+5))/n);
1724       if (maxadd < nzf) maxadd = (n-prow)*(nzf+1);
1725       jmax += maxadd;
1726       xi = (int *) PetscMalloc( jmax*sizeof(int) );CHKPTRQ(xi);
1727       PetscMemcpy(xi,ajnew,ainew[prow]*sizeof(int));
1728       PetscFree(ajnew);
1729       ajnew = xi;
1730       /* allocate a longer ajfill */
1731       xi = (int *) PetscMalloc( jmax*sizeof(int) );CHKPTRQ(xi);
1732       PetscMemcpy(xi,ajfill,ainew[prow]*sizeof(int));
1733       PetscFree(ajfill);
1734       ajfill = xi;
1735       realloc++;
1736     }
1737     xi          = ajnew + ainew[prow];
1738     flev        = ajfill + ainew[prow];
1739     dloc[prow]  = nzi;
1740     fm          = fill[n];
1741     while (nzf--) {
1742       *xi++   = fm;
1743       *flev++ = im[fm];
1744       fm      = fill[fm];
1745     }
1746   }
1747   PetscFree(ajfill);
1748   ierr = ISRestoreIndices(isrow,&r); CHKERRQ(ierr);
1749   ierr = ISRestoreIndices(isicol,&ic); CHKERRQ(ierr);
1750   ierr = ISDestroy(isicol); CHKERRQ(ierr);
1751   PetscFree(fill); PetscFree(im);
1752 
1753   {
1754     double af = ((double)ainew[n])/((double)ai[n]);
1755     PLogInfo(A,"Info:MatILUFactorSymbolic_SeqBAIJ:Reallocs %d Fill ratio:given %g needed %g\n",
1756              realloc,f,af);
1757     PLogInfo(A,"Info:MatILUFactorSymbolic_SeqBAIJ:Run with -pc_ilu_fill %g or use \n",af);
1758     PLogInfo(A,"Info:MatILUFactorSymbolic_SeqBAIJ:PCILUSetFill(pc,%g);\n",af);
1759     PLogInfo(A,"Info:MatILUFactorSymbolic_SeqBAIJ:for best performance.\n");
1760   }
1761 
1762   /* put together the new matrix */
1763   ierr = MatCreateSeqBAIJ(A->comm,bs,bs*n,bs*n,0,PETSC_NULL,fact);CHKERRQ(ierr);
1764   b = (Mat_SeqBAIJ *) (*fact)->data;
1765   PetscFree(b->imax);
1766   b->singlemalloc = 0;
1767   len = bs2*ainew[n]*sizeof(Scalar);
1768   /* the next line frees the default space generated by the Create() */
1769   PetscFree(b->a); PetscFree(b->ilen);
1770   b->a          = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a);
1771   b->j          = ajnew;
1772   b->i          = ainew;
1773   for ( i=0; i<n; i++ ) dloc[i] += ainew[i];
1774   b->diag       = dloc;
1775   b->ilen       = 0;
1776   b->imax       = 0;
1777   b->row        = isrow;
1778   b->col        = iscol;
1779   b->solve_work = (Scalar *) PetscMalloc( (bs*n+bs)*sizeof(Scalar));CHKPTRQ(b->solve_work);
1780   /* In b structure:  Free imax, ilen, old a, old j.
1781      Allocate dloc, solve_work, new a, new j */
1782   PLogObjectMemory(*fact,(ainew[n]-n)*(sizeof(int))+bs2*ainew[n]*sizeof(Scalar));
1783   b->maxnz          = b->nz = ainew[n];
1784   (*fact)->factor   = FACTOR_LU;
1785 
1786   (*fact)->info.factor_mallocs    = realloc;
1787   (*fact)->info.fill_ratio_given  = f;
1788   (*fact)->info.fill_ratio_needed = ((double)ainew[n])/((double)ai[prow]);
1789 
1790   PetscFunctionReturn(0);
1791 }
1792 
1793 
1794 
1795 
1796