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