xref: /petsc/src/mat/impls/baij/seq/baijfact.c (revision c05d87d6c2fa1d0952fe676de20339d09daa79bc)
1 #define PETSCMAT_DLL
2 
3 /*
4     Factorization code for BAIJ format.
5 */
6 #include "../src/mat/impls/baij/seq/baij.h"
7 #include "../src/mat/blockinvert.h"
8 
9 #undef __FUNCT__
10 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2_newdatastruct"
11 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_newdatastruct(Mat B,Mat A,const MatFactorInfo *info)
12 {
13   Mat            C=B;
14   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ *)C->data;
15   IS             isrow = b->row,isicol = b->icol;
16   PetscErrorCode ierr;
17   const PetscInt *r,*ic,*ics;
18   PetscInt       i,j,k,n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
19   PetscInt       *ajtmp,*bjtmp,nz,nzL,row,*bdiag=b->diag,*pj;
20   MatScalar      *rtmp,*pc,*mwork,*v,*pv,*aa=a->a;
21   PetscInt       bs2 = a->bs2,flg;
22   PetscReal      shift = info->shiftinblocks;
23 
24   PetscFunctionBegin;
25   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
26   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
27 
28   /* generate work space needed by the factorization */
29   ierr = PetscMalloc2(bs2*n,MatScalar,&rtmp,bs2,MatScalar,&mwork);CHKERRQ(ierr);
30   ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr);
31   ics  = ic;
32 
33   for (i=0; i<n; i++){
34     /* zero rtmp */
35     /* L part */
36     nz    = bi[i+1] - bi[i];
37     bjtmp = bj + bi[i];
38     for  (j=0; j<nz; j++){
39       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
40     }
41 
42     /* U part */
43     nz = bdiag[i] - bdiag[i+1];
44     bjtmp = bj + bdiag[i+1]+1;
45     for  (j=0; j<nz; j++){
46       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
47     }
48 
49     /* load in initial (unfactored row) */
50     nz    = ai[r[i]+1] - ai[r[i]];
51     ajtmp = aj + ai[r[i]];
52     v     = aa + bs2*ai[r[i]];
53     for (j=0; j<nz; j++) {
54       ierr = PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr);
55     }
56 
57     /* elimination */
58     bjtmp = bj + bi[i];
59     nzL   = bi[i+1] - bi[i];
60     for(k=0;k < nzL;k++) {
61       row = bjtmp[k];
62       pc = rtmp + bs2*row;
63       for (flg=0,j=0; j<bs2; j++) { if (pc[j]!=0.0) { flg = 1; break; }}
64       if (flg) {
65         pv = b->a + bs2*bdiag[row];
66         /* Kernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
67         ierr = Kernel_A_gets_A_times_B_2(pc,pv,mwork);CHKERRQ(ierr);
68 
69         pj = b->j + bdiag[row+1]+1; /* begining of U(row,:) */
70         pv = b->a + bs2*(bdiag[row+1]+1);
71         nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries inU(row,:), excluding diag */
72         for (j=0; j<nz; j++) {
73           /* Kernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
74           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
75           v    = rtmp + 4*pj[j];
76           ierr = Kernel_A_gets_A_minus_B_times_C_2(v,pc,pv);CHKERRQ(ierr);
77           pv  += 4;
78         }
79         ierr = PetscLogFlops(16*nz+12);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
80       }
81     }
82 
83     /* finished row so stick it into b->a */
84     /* L part */
85     pv   = b->a + bs2*bi[i] ;
86     pj   = b->j + bi[i] ;
87     nz   = bi[i+1] - bi[i];
88     for (j=0; j<nz; j++) {
89       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
90     }
91 
92     /* Mark diagonal and invert diagonal for simplier triangular solves */
93     pv   = b->a + bs2*bdiag[i];
94     pj   = b->j + bdiag[i];
95     ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr);
96     /* ierr = Kernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); */
97     ierr = Kernel_A_gets_inverse_A_2(pv,shift);CHKERRQ(ierr);
98 
99     /* U part */
100     pv = b->a + bs2*(bdiag[i+1]+1);
101     pj = b->j + bdiag[i+1]+1;
102     nz = bdiag[i] - bdiag[i+1] - 1;
103     for (j=0; j<nz; j++){
104       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
105     }
106   }
107 
108   ierr = PetscFree2(rtmp,mwork);CHKERRQ(ierr);
109   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
110   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
111 
112   C->assembled = PETSC_TRUE;
113   ierr = PetscLogFlops(1.3333*bs2*n);CHKERRQ(ierr); /* from inverting diagonal blocks */
114   PetscFunctionReturn(0);
115 }
116 
117 #undef __FUNCT__
118 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_newdatastruct"
119 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_newdatastruct(Mat B,Mat A,const MatFactorInfo *info)
120 {
121   Mat            C=B;
122   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ *)C->data;
123   PetscErrorCode ierr;
124   PetscInt       i,j,k,n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
125   PetscInt       *ajtmp,*bjtmp,nz,nzL,row,*bdiag=b->diag,*pj;
126   MatScalar      *rtmp,*pc,*mwork,*v,*pv,*aa=a->a;
127   PetscInt       bs2 = a->bs2,flg;
128   PetscReal      shift = info->shiftinblocks;
129 
130   PetscFunctionBegin;
131   /* generate work space needed by the factorization */
132   ierr = PetscMalloc2(bs2*n,MatScalar,&rtmp,bs2,MatScalar,&mwork);CHKERRQ(ierr);
133   ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr);
134 
135   for (i=0; i<n; i++){
136     /* zero rtmp */
137     /* L part */
138     nz    = bi[i+1] - bi[i];
139     bjtmp = bj + bi[i];
140     for  (j=0; j<nz; j++){
141       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
142     }
143 
144     /* U part */
145     nz = bdiag[i] - bdiag[i+1];
146     bjtmp = bj + bdiag[i+1]+1;
147     for  (j=0; j<nz; j++){
148       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
149     }
150 
151     /* load in initial (unfactored row) */
152     nz    = ai[i+1] - ai[i];
153     ajtmp = aj + ai[i];
154     v     = aa + bs2*ai[i];
155     for (j=0; j<nz; j++) {
156       ierr = PetscMemcpy(rtmp+bs2*ajtmp[j],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr);
157     }
158 
159     /* elimination */
160     bjtmp = bj + bi[i];
161     nzL   = bi[i+1] - bi[i];
162     for(k=0;k < nzL;k++) {
163       row = bjtmp[k];
164       pc = rtmp + bs2*row;
165       for (flg=0,j=0; j<bs2; j++) { if (pc[j]!=0.0) { flg = 1; break; }}
166       if (flg) {
167         pv = b->a + bs2*bdiag[row];
168         /* Kernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
169         ierr = Kernel_A_gets_A_times_B_2(pc,pv,mwork);CHKERRQ(ierr);
170 
171         pj = b->j + bdiag[row+1]+1; /* beginning of U(row,:) */
172 	pv = b->a + bs2*(bdiag[row+1]+1);
173 	nz = bdiag[row]-bdiag[row+1] - 1; /* num of entries in U(row,:) excluding diag */
174         for (j=0; j<nz; j++) {
175           /* Kernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
176           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
177           v    = rtmp + 4*pj[j];
178           ierr = Kernel_A_gets_A_minus_B_times_C_2(v,pc,pv);CHKERRQ(ierr);
179           pv  += 4;
180         }
181         ierr = PetscLogFlops(16*nz+12);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
182       }
183     }
184 
185     /* finished row so stick it into b->a */
186     /* L part */
187     pv   = b->a + bs2*bi[i] ;
188     pj   = b->j + bi[i] ;
189     nz   = bi[i+1] - bi[i];
190     for (j=0; j<nz; j++) {
191       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
192     }
193 
194     /* Mark diagonal and invert diagonal for simplier triangular solves */
195     pv   = b->a + bs2*bdiag[i];
196     pj   = b->j + bdiag[i];
197     ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr);
198     /* ierr = Kernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); */
199     ierr = Kernel_A_gets_inverse_A_2(pv,shift);CHKERRQ(ierr);
200 
201     /* U part */
202     /*
203     pv = b->a + bs2*bi[2*n-i];
204     pj = b->j + bi[2*n-i];
205     nz = bi[2*n-i+1] - bi[2*n-i] - 1;
206     */
207     pv = b->a + bs2*(bdiag[i+1]+1);
208     pj = b->j + bdiag[i+1]+1;
209     nz = bdiag[i] - bdiag[i+1] - 1;
210     for (j=0; j<nz; j++){
211       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
212     }
213   }
214 
215   ierr = PetscFree2(rtmp,mwork);CHKERRQ(ierr);
216   C->assembled = PETSC_TRUE;
217   ierr = PetscLogFlops(1.3333*bs2*n);CHKERRQ(ierr); /* from inverting diagonal blocks */
218   PetscFunctionReturn(0);
219 }
220 
221 #undef __FUNCT__
222 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2"
223 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2(Mat B,Mat A,const MatFactorInfo *info)
224 {
225   Mat            C = B;
226   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data;
227   IS             isrow = b->row,isicol = b->icol;
228   PetscErrorCode ierr;
229   const PetscInt *r,*ic;
230   PetscInt       i,j,n = a->mbs,*bi = b->i,*bj = b->j;
231   PetscInt       *ajtmpold,*ajtmp,nz,row;
232   PetscInt       *diag_offset=b->diag,idx,*ai=a->i,*aj=a->j,*pj;
233   MatScalar      *pv,*v,*rtmp,m1,m2,m3,m4,*pc,*w,*x,x1,x2,x3,x4;
234   MatScalar      p1,p2,p3,p4;
235   MatScalar      *ba = b->a,*aa = a->a;
236   PetscReal      shift = info->shiftinblocks;
237 
238   PetscFunctionBegin;
239   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
240   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
241   ierr  = PetscMalloc(4*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
242 
243   for (i=0; i<n; i++) {
244     nz    = bi[i+1] - bi[i];
245     ajtmp = bj + bi[i];
246     for  (j=0; j<nz; j++) {
247       x = rtmp+4*ajtmp[j]; x[0] = x[1] = x[2] = x[3] = 0.0;
248     }
249     /* load in initial (unfactored row) */
250     idx      = r[i];
251     nz       = ai[idx+1] - ai[idx];
252     ajtmpold = aj + ai[idx];
253     v        = aa + 4*ai[idx];
254     for (j=0; j<nz; j++) {
255       x    = rtmp+4*ic[ajtmpold[j]];
256       x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3];
257       v    += 4;
258     }
259     row = *ajtmp++;
260     while (row < i) {
261       pc = rtmp + 4*row;
262       p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3];
263       if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0) {
264         pv = ba + 4*diag_offset[row];
265         pj = bj + diag_offset[row] + 1;
266         x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
267         pc[0] = m1 = p1*x1 + p3*x2;
268         pc[1] = m2 = p2*x1 + p4*x2;
269         pc[2] = m3 = p1*x3 + p3*x4;
270         pc[3] = m4 = p2*x3 + p4*x4;
271         nz = bi[row+1] - diag_offset[row] - 1;
272         pv += 4;
273         for (j=0; j<nz; j++) {
274           x1   = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
275           x    = rtmp + 4*pj[j];
276           x[0] -= m1*x1 + m3*x2;
277           x[1] -= m2*x1 + m4*x2;
278           x[2] -= m1*x3 + m3*x4;
279           x[3] -= m2*x3 + m4*x4;
280           pv   += 4;
281         }
282         ierr = PetscLogFlops(16.0*nz+12.0);CHKERRQ(ierr);
283       }
284       row = *ajtmp++;
285     }
286     /* finished row so stick it into b->a */
287     pv = ba + 4*bi[i];
288     pj = bj + bi[i];
289     nz = bi[i+1] - bi[i];
290     for (j=0; j<nz; j++) {
291       x     = rtmp+4*pj[j];
292       pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3];
293       pv   += 4;
294     }
295     /* invert diagonal block */
296     w = ba + 4*diag_offset[i];
297     ierr = Kernel_A_gets_inverse_A_2(w,shift);CHKERRQ(ierr);
298   }
299 
300   ierr = PetscFree(rtmp);CHKERRQ(ierr);
301   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
302   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
303   C->ops->solve          = MatSolve_SeqBAIJ_2;
304   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2;
305   C->assembled = PETSC_TRUE;
306   ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
307   PetscFunctionReturn(0);
308 }
309 /*
310       Version for when blocks are 2 by 2 Using natural ordering
311 */
312 #undef __FUNCT__
313 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering"
314 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info)
315 {
316   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data;
317   PetscErrorCode ierr;
318   PetscInt       i,j,n = a->mbs,*bi = b->i,*bj = b->j;
319   PetscInt       *ajtmpold,*ajtmp,nz,row;
320   PetscInt       *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj;
321   MatScalar      *pv,*v,*rtmp,*pc,*w,*x;
322   MatScalar      p1,p2,p3,p4,m1,m2,m3,m4,x1,x2,x3,x4;
323   MatScalar      *ba = b->a,*aa = a->a;
324   PetscReal      shift = info->shiftinblocks;
325 
326   PetscFunctionBegin;
327   ierr = PetscMalloc(4*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
328   for (i=0; i<n; i++) {
329     nz    = bi[i+1] - bi[i];
330     ajtmp = bj + bi[i];
331     for  (j=0; j<nz; j++) {
332       x = rtmp+4*ajtmp[j];
333       x[0]  = x[1]  = x[2]  = x[3]  = 0.0;
334     }
335     /* load in initial (unfactored row) */
336     nz       = ai[i+1] - ai[i];
337     ajtmpold = aj + ai[i];
338     v        = aa + 4*ai[i];
339     for (j=0; j<nz; j++) {
340       x    = rtmp+4*ajtmpold[j];
341       x[0]  = v[0];  x[1]  = v[1];  x[2]  = v[2];  x[3]  = v[3];
342       v    += 4;
343     }
344     row = *ajtmp++;
345     while (row < i) {
346       pc  = rtmp + 4*row;
347       p1  = pc[0];  p2  = pc[1];  p3  = pc[2];  p4  = pc[3];
348       if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0) {
349         pv = ba + 4*diag_offset[row];
350         pj = bj + diag_offset[row] + 1;
351         x1  = pv[0];  x2  = pv[1];  x3  = pv[2];  x4  = pv[3];
352         pc[0] = m1 = p1*x1 + p3*x2;
353         pc[1] = m2 = p2*x1 + p4*x2;
354         pc[2] = m3 = p1*x3 + p3*x4;
355         pc[3] = m4 = p2*x3 + p4*x4;
356         nz = bi[row+1] - diag_offset[row] - 1;
357         pv += 4;
358         for (j=0; j<nz; j++) {
359           x1   = pv[0];  x2  = pv[1];   x3 = pv[2];  x4  = pv[3];
360           x    = rtmp + 4*pj[j];
361           x[0] -= m1*x1 + m3*x2;
362           x[1] -= m2*x1 + m4*x2;
363           x[2] -= m1*x3 + m3*x4;
364           x[3] -= m2*x3 + m4*x4;
365           pv   += 4;
366         }
367         ierr = PetscLogFlops(16.0*nz+12.0);CHKERRQ(ierr);
368       }
369       row = *ajtmp++;
370     }
371     /* finished row so stick it into b->a */
372     pv = ba + 4*bi[i];
373     pj = bj + bi[i];
374     nz = bi[i+1] - bi[i];
375     for (j=0; j<nz; j++) {
376       x      = rtmp+4*pj[j];
377       pv[0]  = x[0];  pv[1]  = x[1];  pv[2]  = x[2];  pv[3]  = x[3];
378       /*
379       printf(" col %d:",pj[j]);
380       PetscInt j1;
381       for (j1=0; j1<4; j1++) printf(" %g,",*(pv+j1));
382       printf("\n");
383       */
384       pv   += 4;
385     }
386     /* invert diagonal block */
387     w = ba + 4*diag_offset[i];
388     /*
389     printf(" \n%d -th: diag: ",i);
390     for (j=0; j<4; j++){
391       printf(" %g,",w[j]);
392     }
393     printf("\n----------------------------\n");
394     */
395     ierr = Kernel_A_gets_inverse_A_2(w,shift);CHKERRQ(ierr);
396   }
397 
398   ierr = PetscFree(rtmp);CHKERRQ(ierr);
399   C->ops->solve          = MatSolve_SeqBAIJ_2_NaturalOrdering;
400   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_NaturalOrdering;
401   C->assembled = PETSC_TRUE;
402   ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
403   PetscFunctionReturn(0);
404 }
405 
406 /* ----------------------------------------------------------- */
407 /*
408      Version for when blocks are 1 by 1.
409 */
410 #undef __FUNCT__
411 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_1"
412 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_1(Mat C,Mat A,const MatFactorInfo *info)
413 {
414   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data;
415   IS             isrow = b->row,isicol = b->icol;
416   PetscErrorCode ierr;
417   const PetscInt *r,*ic;
418   PetscInt       i,j,n = a->mbs,*bi = b->i,*bj = b->j;
419   PetscInt       *ajtmpold,*ajtmp,nz,row,*ai = a->i,*aj = a->j;
420   PetscInt       *diag_offset = b->diag,diag,*pj;
421   MatScalar      *pv,*v,*rtmp,multiplier,*pc;
422   MatScalar      *ba = b->a,*aa = a->a;
423   PetscTruth     row_identity, col_identity;
424 
425   PetscFunctionBegin;
426   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
427   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
428   ierr  = PetscMalloc((n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
429 
430   for (i=0; i<n; i++) {
431     nz    = bi[i+1] - bi[i];
432     ajtmp = bj + bi[i];
433     for  (j=0; j<nz; j++) rtmp[ajtmp[j]] = 0.0;
434 
435     /* load in initial (unfactored row) */
436     nz       = ai[r[i]+1] - ai[r[i]];
437     ajtmpold = aj + ai[r[i]];
438     v        = aa + ai[r[i]];
439     for (j=0; j<nz; j++) rtmp[ic[ajtmpold[j]]] =  v[j];
440 
441     row = *ajtmp++;
442     while (row < i) {
443       pc = rtmp + row;
444       if (*pc != 0.0) {
445         pv         = ba + diag_offset[row];
446         pj         = bj + diag_offset[row] + 1;
447         multiplier = *pc * *pv++;
448         *pc        = multiplier;
449         nz         = bi[row+1] - diag_offset[row] - 1;
450         for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j];
451         ierr = PetscLogFlops(1.0+2.0*nz);CHKERRQ(ierr);
452       }
453       row = *ajtmp++;
454     }
455     /* finished row so stick it into b->a */
456     pv = ba + bi[i];
457     pj = bj + bi[i];
458     nz = bi[i+1] - bi[i];
459     for (j=0; j<nz; j++) {pv[j] = rtmp[pj[j]];}
460     diag = diag_offset[i] - bi[i];
461     /* check pivot entry for current row */
462     if (pv[diag] == 0.0) {
463       SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot: row in original ordering %D in permuted ordering %D",r[i],i);
464     }
465     pv[diag] = 1.0/pv[diag];
466   }
467 
468   ierr = PetscFree(rtmp);CHKERRQ(ierr);
469   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
470   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
471   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
472   ierr = ISIdentity(isicol,&col_identity);CHKERRQ(ierr);
473   if (row_identity && col_identity) {
474     C->ops->solve          = MatSolve_SeqBAIJ_1_NaturalOrdering;
475     C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1_NaturalOrdering;
476   } else {
477     C->ops->solve          = MatSolve_SeqBAIJ_1;
478     C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1;
479   }
480   C->assembled = PETSC_TRUE;
481   ierr = PetscLogFlops(C->cmap->n);CHKERRQ(ierr);
482   PetscFunctionReturn(0);
483 }
484 
485 EXTERN_C_BEGIN
486 #undef __FUNCT__
487 #define __FUNCT__ "MatGetFactor_seqbaij_petsc"
488 PetscErrorCode MatGetFactor_seqbaij_petsc(Mat A,MatFactorType ftype,Mat *B)
489 {
490   PetscInt           n = A->rmap->n;
491   PetscErrorCode     ierr;
492 
493   PetscFunctionBegin;
494   ierr = MatCreate(((PetscObject)A)->comm,B);CHKERRQ(ierr);
495   ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr);
496   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU || ftype == MAT_FACTOR_ILUDT) {
497     ierr = MatSetType(*B,MATSEQBAIJ);CHKERRQ(ierr);
498     (*B)->ops->lufactorsymbolic  = MatLUFactorSymbolic_SeqBAIJ;
499     (*B)->ops->ilufactorsymbolic = MatILUFactorSymbolic_SeqBAIJ;
500     (*B)->ops->iludtfactor       = MatILUDTFactor_SeqBAIJ;
501   } else if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) {
502     ierr = MatSetType(*B,MATSEQSBAIJ);CHKERRQ(ierr);
503     ierr = MatSeqSBAIJSetPreallocation(*B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
504     (*B)->ops->iccfactorsymbolic      = MatICCFactorSymbolic_SeqBAIJ;
505     (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqBAIJ;
506   } else SETERRQ(PETSC_ERR_SUP,"Factor type not supported");
507   (*B)->factor = ftype;
508   PetscFunctionReturn(0);
509 }
510 EXTERN_C_END
511 
512 EXTERN_C_BEGIN
513 #undef __FUNCT__
514 #define __FUNCT__ "MatGetFactorAvailable_seqaij_petsc"
515 PetscErrorCode MatGetFactorAvailable_seqbaij_petsc(Mat A,MatFactorType ftype,PetscTruth *flg)
516 {
517   PetscFunctionBegin;
518   *flg = PETSC_TRUE;
519   PetscFunctionReturn(0);
520 }
521 EXTERN_C_END
522 
523 /* ----------------------------------------------------------- */
524 #undef __FUNCT__
525 #define __FUNCT__ "MatLUFactor_SeqBAIJ"
526 PetscErrorCode MatLUFactor_SeqBAIJ(Mat A,IS row,IS col,const MatFactorInfo *info)
527 {
528   PetscErrorCode ierr;
529   Mat            C;
530 
531   PetscFunctionBegin;
532   ierr = MatGetFactor(A,MAT_SOLVER_PETSC,MAT_FACTOR_LU,&C);CHKERRQ(ierr);
533   ierr = MatLUFactorSymbolic(C,A,row,col,info);CHKERRQ(ierr);
534   ierr = MatLUFactorNumeric(C,A,info);CHKERRQ(ierr);
535   A->ops->solve            = C->ops->solve;
536   A->ops->solvetranspose   = C->ops->solvetranspose;
537   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
538   ierr = PetscLogObjectParent(A,((Mat_SeqBAIJ*)(A->data))->icol);CHKERRQ(ierr);
539   PetscFunctionReturn(0);
540 }
541 
542 #include "../src/mat/impls/sbaij/seq/sbaij.h"
543 #undef __FUNCT__
544 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N"
545 PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N(Mat C,Mat A,const MatFactorInfo *info)
546 {
547   PetscErrorCode ierr;
548   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data;
549   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
550   IS             ip=b->row;
551   const PetscInt *rip;
552   PetscInt       i,j,mbs=a->mbs,bs=A->rmap->bs,*bi=b->i,*bj=b->j,*bcol;
553   PetscInt       *ai=a->i,*aj=a->j;
554   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
555   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
556   PetscReal      zeropivot,rs,shiftnz;
557   PetscReal      shiftpd;
558   ChShift_Ctx    sctx;
559   PetscInt       newshift;
560 
561   PetscFunctionBegin;
562   if (bs > 1) {
563     if (!a->sbaijMat){
564       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
565     }
566     ierr = (a->sbaijMat)->ops->choleskyfactornumeric(C,a->sbaijMat,info);CHKERRQ(ierr);
567     ierr = MatDestroy(a->sbaijMat);CHKERRQ(ierr);
568     a->sbaijMat = PETSC_NULL;
569     PetscFunctionReturn(0);
570   }
571 
572   /* initialization */
573   shiftnz   = info->shiftnz;
574   shiftpd   = info->shiftpd;
575   zeropivot = info->zeropivot;
576 
577   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
578   ierr = PetscMalloc3(mbs,MatScalar,&rtmp,mbs,PetscInt,&il,mbs,PetscInt,&jl);CHKERRQ(ierr);
579 
580   sctx.shift_amount = 0;
581   sctx.nshift       = 0;
582   do {
583     sctx.chshift = PETSC_FALSE;
584     for (i=0; i<mbs; i++) {
585       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
586     }
587 
588     for (k = 0; k<mbs; k++){
589       bval = ba + bi[k];
590       /* initialize k-th row by the perm[k]-th row of A */
591       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
592       for (j = jmin; j < jmax; j++){
593         col = rip[aj[j]];
594         if (col >= k){ /* only take upper triangular entry */
595           rtmp[col] = aa[j];
596           *bval++  = 0.0; /* for in-place factorization */
597         }
598       }
599 
600       /* shift the diagonal of the matrix */
601       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
602 
603       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
604       dk = rtmp[k];
605       i = jl[k]; /* first row to be added to k_th row  */
606 
607       while (i < k){
608         nexti = jl[i]; /* next row to be added to k_th row */
609 
610         /* compute multiplier, update diag(k) and U(i,k) */
611         ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
612         uikdi = - ba[ili]*ba[bi[i]];  /* diagonal(k) */
613         dk += uikdi*ba[ili];
614         ba[ili] = uikdi; /* -U(i,k) */
615 
616         /* add multiple of row i to k-th row */
617         jmin = ili + 1; jmax = bi[i+1];
618         if (jmin < jmax){
619           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
620           /* update il and jl for row i */
621           il[i] = jmin;
622           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
623         }
624         i = nexti;
625       }
626 
627       /* shift the diagonals when zero pivot is detected */
628       /* compute rs=sum of abs(off-diagonal) */
629       rs   = 0.0;
630       jmin = bi[k]+1;
631       nz   = bi[k+1] - jmin;
632       if (nz){
633         bcol = bj + jmin;
634         while (nz--){
635           rs += PetscAbsScalar(rtmp[*bcol]);
636           bcol++;
637         }
638       }
639 
640       sctx.rs = rs;
641       sctx.pv = dk;
642       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
643       if (newshift == 1) break;
644 
645       /* copy data into U(k,:) */
646       ba[bi[k]] = 1.0/dk; /* U(k,k) */
647       jmin = bi[k]+1; jmax = bi[k+1];
648       if (jmin < jmax) {
649         for (j=jmin; j<jmax; j++){
650           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
651         }
652         /* add the k-th row into il and jl */
653         il[k] = jmin;
654         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
655       }
656     }
657   } while (sctx.chshift);
658   ierr = PetscFree3(rtmp,il,jl);CHKERRQ(ierr);
659 
660   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
661   C->assembled    = PETSC_TRUE;
662   C->preallocated = PETSC_TRUE;
663   ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr);
664   if (sctx.nshift){
665     if (shiftpd) {
666       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
667     } else if (shiftnz) {
668       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
669     }
670   }
671   PetscFunctionReturn(0);
672 }
673 
674 #undef __FUNCT__
675 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering"
676 PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info)
677 {
678   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data;
679   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
680   PetscErrorCode ierr;
681   PetscInt       i,j,am=a->mbs;
682   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
683   PetscInt       k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz;
684   MatScalar      *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval;
685   PetscReal      zeropivot,rs,shiftnz;
686   PetscReal      shiftpd;
687   ChShift_Ctx    sctx;
688   PetscInt       newshift;
689 
690   PetscFunctionBegin;
691   /* initialization */
692   shiftnz   = info->shiftnz;
693   shiftpd   = info->shiftpd;
694   zeropivot = info->zeropivot;
695 
696   ierr = PetscMalloc3(am,MatScalar,&rtmp,am,PetscInt,&il,am,PetscInt,&jl);CHKERRQ(ierr);
697 
698   sctx.shift_amount = 0;
699   sctx.nshift       = 0;
700   do {
701     sctx.chshift = PETSC_FALSE;
702     for (i=0; i<am; i++) {
703       rtmp[i] = 0.0; jl[i] = am; il[0] = 0;
704     }
705 
706     for (k = 0; k<am; k++){
707     /* initialize k-th row with elements nonzero in row perm(k) of A */
708       nz   = ai[k+1] - ai[k];
709       acol = aj + ai[k];
710       aval = aa + ai[k];
711       bval = ba + bi[k];
712       while (nz -- ){
713         if (*acol < k) { /* skip lower triangular entries */
714           acol++; aval++;
715         } else {
716           rtmp[*acol++] = *aval++;
717           *bval++       = 0.0; /* for in-place factorization */
718         }
719       }
720 
721       /* shift the diagonal of the matrix */
722       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
723 
724       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
725       dk = rtmp[k];
726       i  = jl[k]; /* first row to be added to k_th row  */
727 
728       while (i < k){
729         nexti = jl[i]; /* next row to be added to k_th row */
730         /* compute multiplier, update D(k) and U(i,k) */
731         ili   = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
732         uikdi = - ba[ili]*ba[bi[i]];
733         dk   += uikdi*ba[ili];
734         ba[ili] = uikdi; /* -U(i,k) */
735 
736         /* add multiple of row i to k-th row ... */
737         jmin = ili + 1;
738         nz   = bi[i+1] - jmin;
739         if (nz > 0){
740           bcol = bj + jmin;
741           bval = ba + jmin;
742           while (nz --) rtmp[*bcol++] += uikdi*(*bval++);
743           /* update il and jl for i-th row */
744           il[i] = jmin;
745           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
746         }
747         i = nexti;
748       }
749 
750       /* shift the diagonals when zero pivot is detected */
751       /* compute rs=sum of abs(off-diagonal) */
752       rs   = 0.0;
753       jmin = bi[k]+1;
754       nz   = bi[k+1] - jmin;
755       if (nz){
756         bcol = bj + jmin;
757         while (nz--){
758           rs += PetscAbsScalar(rtmp[*bcol]);
759           bcol++;
760         }
761       }
762 
763       sctx.rs = rs;
764       sctx.pv = dk;
765       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
766       if (newshift == 1) break;    /* sctx.shift_amount is updated */
767 
768       /* copy data into U(k,:) */
769       ba[bi[k]] = 1.0/dk;
770       jmin      = bi[k]+1;
771       nz        = bi[k+1] - jmin;
772       if (nz){
773         bcol = bj + jmin;
774         bval = ba + jmin;
775         while (nz--){
776           *bval++       = rtmp[*bcol];
777           rtmp[*bcol++] = 0.0;
778         }
779         /* add k-th row into il and jl */
780         il[k] = jmin;
781         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
782       }
783     }
784   } while (sctx.chshift);
785   ierr = PetscFree3(rtmp,il,jl);CHKERRQ(ierr);
786 
787   C->ops->solve                 = MatSolve_SeqSBAIJ_1_NaturalOrdering;
788   C->ops->solvetranspose        = MatSolve_SeqSBAIJ_1_NaturalOrdering;
789   C->assembled    = PETSC_TRUE;
790   C->preallocated = PETSC_TRUE;
791   ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr);
792     if (sctx.nshift){
793     if (shiftnz) {
794       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
795     } else if (shiftpd) {
796       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
797     }
798   }
799   PetscFunctionReturn(0);
800 }
801 
802 #include "petscbt.h"
803 #include "../src/mat/utils/freespace.h"
804 #undef __FUNCT__
805 #define __FUNCT__ "MatICCFactorSymbolic_SeqBAIJ"
806 PetscErrorCode MatICCFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
807 {
808   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data;
809   Mat_SeqSBAIJ       *b;
810   Mat                B;
811   PetscErrorCode     ierr;
812   PetscTruth         perm_identity;
813   PetscInt           reallocs=0,i,*ai=a->i,*aj=a->j,am=a->mbs,bs=A->rmap->bs,*ui;
814   const PetscInt     *rip;
815   PetscInt           jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow;
816   PetscInt           nlnk,*lnk,*lnk_lvl=PETSC_NULL,ncols,ncols_upper,*cols,*cols_lvl,*uj,**uj_ptr,**uj_lvl_ptr;
817   PetscReal          fill=info->fill,levels=info->levels;
818   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
819   PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
820   PetscBT            lnkbt;
821 
822   PetscFunctionBegin;
823   if (bs > 1){
824     if (!a->sbaijMat){
825       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
826     }
827     (fact)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqSBAIJ;  /* undue the change made in MatGetFactor_seqbaij_petsc */
828     ierr = MatICCFactorSymbolic(fact,a->sbaijMat,perm,info);CHKERRQ(ierr);
829     PetscFunctionReturn(0);
830   }
831 
832   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
833   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
834 
835   /* special case that simply copies fill pattern */
836   if (!levels && perm_identity) {
837     ierr = MatMarkDiagonal_SeqBAIJ(A);CHKERRQ(ierr);
838     ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
839     for (i=0; i<am; i++) {
840       ui[i] = ai[i+1] - a->diag[i]; /* ui: rowlengths - changes when !perm_identity */
841     }
842     B = fact;
843     ierr = MatSeqSBAIJSetPreallocation(B,1,0,ui);CHKERRQ(ierr);
844 
845 
846     b  = (Mat_SeqSBAIJ*)B->data;
847     uj = b->j;
848     for (i=0; i<am; i++) {
849       aj = a->j + a->diag[i];
850       for (j=0; j<ui[i]; j++){
851         *uj++ = *aj++;
852       }
853       b->ilen[i] = ui[i];
854     }
855     ierr = PetscFree(ui);CHKERRQ(ierr);
856     B->factor = MAT_FACTOR_NONE;
857     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
858     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
859     B->factor = MAT_FACTOR_ICC;
860 
861     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
862     PetscFunctionReturn(0);
863   }
864 
865   /* initialization */
866   ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
867   ui[0] = 0;
868   ierr  = PetscMalloc((2*am+1)*sizeof(PetscInt),&cols_lvl);CHKERRQ(ierr);
869 
870   /* jl: linked list for storing indices of the pivot rows
871      il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
872   ierr = PetscMalloc4(am,PetscInt*,&uj_ptr,am,PetscInt*,&uj_lvl_ptr,am,PetscInt,&il,am,PetscInt,&jl);CHKERRQ(ierr);
873   for (i=0; i<am; i++){
874     jl[i] = am; il[i] = 0;
875   }
876 
877   /* create and initialize a linked list for storing column indices of the active row k */
878   nlnk = am + 1;
879   ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
880 
881   /* initial FreeSpace size is fill*(ai[am]+1) */
882   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
883   current_space = free_space;
884   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr);
885   current_space_lvl = free_space_lvl;
886 
887   for (k=0; k<am; k++){  /* for each active row k */
888     /* initialize lnk by the column indices of row rip[k] of A */
889     nzk   = 0;
890     ncols = ai[rip[k]+1] - ai[rip[k]];
891     ncols_upper = 0;
892     cols        = cols_lvl + am;
893     for (j=0; j<ncols; j++){
894       i = rip[*(aj + ai[rip[k]] + j)];
895       if (i >= k){ /* only take upper triangular entry */
896         cols[ncols_upper] = i;
897         cols_lvl[ncols_upper] = -1;  /* initialize level for nonzero entries */
898         ncols_upper++;
899       }
900     }
901     ierr = PetscIncompleteLLAdd(ncols_upper,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
902     nzk += nlnk;
903 
904     /* update lnk by computing fill-in for each pivot row to be merged in */
905     prow = jl[k]; /* 1st pivot row */
906 
907     while (prow < k){
908       nextprow = jl[prow];
909 
910       /* merge prow into k-th row */
911       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
912       jmax = ui[prow+1];
913       ncols = jmax-jmin;
914       i     = jmin - ui[prow];
915       cols = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
916       for (j=0; j<ncols; j++) cols_lvl[j] = *(uj_lvl_ptr[prow] + i + j);
917       ierr = PetscIncompleteLLAddSorted(ncols,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
918       nzk += nlnk;
919 
920       /* update il and jl for prow */
921       if (jmin < jmax){
922         il[prow] = jmin;
923         j = *cols; jl[prow] = jl[j]; jl[j] = prow;
924       }
925       prow = nextprow;
926     }
927 
928     /* if free space is not available, make more free space */
929     if (current_space->local_remaining<nzk) {
930       i = am - k + 1; /* num of unfactored rows */
931       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
932       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
933       ierr = PetscFreeSpaceGet(i,&current_space_lvl);CHKERRQ(ierr);
934       reallocs++;
935     }
936 
937     /* copy data into free_space and free_space_lvl, then initialize lnk */
938     ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
939 
940     /* add the k-th row into il and jl */
941     if (nzk-1 > 0){
942       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
943       jl[k] = jl[i]; jl[i] = k;
944       il[k] = ui[k] + 1;
945     }
946     uj_ptr[k]     = current_space->array;
947     uj_lvl_ptr[k] = current_space_lvl->array;
948 
949     current_space->array           += nzk;
950     current_space->local_used      += nzk;
951     current_space->local_remaining -= nzk;
952 
953     current_space_lvl->array           += nzk;
954     current_space_lvl->local_used      += nzk;
955     current_space_lvl->local_remaining -= nzk;
956 
957     ui[k+1] = ui[k] + nzk;
958   }
959 
960 #if defined(PETSC_USE_INFO)
961   if (ai[am] != 0) {
962     PetscReal af = ((PetscReal)(2*ui[am]-am))/((PetscReal)ai[am]);
963     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
964     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
965     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
966   } else {
967     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
968   }
969 #endif
970 
971   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
972   ierr = PetscFree4(uj_ptr,uj_lvl_ptr,il,jl);CHKERRQ(ierr);
973   ierr = PetscFree(cols_lvl);CHKERRQ(ierr);
974 
975   /* destroy list of free space and other temporary array(s) */
976   ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
977   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
978   ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
979   ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
980 
981   /* put together the new matrix in MATSEQSBAIJ format */
982   B = fact;
983   ierr = MatSeqSBAIJSetPreallocation(B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
984 
985   b = (Mat_SeqSBAIJ*)B->data;
986   b->singlemalloc = PETSC_FALSE;
987   b->free_a       = PETSC_TRUE;
988   b->free_ij       = PETSC_TRUE;
989   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
990   b->j    = uj;
991   b->i    = ui;
992   b->diag = 0;
993   b->ilen = 0;
994   b->imax = 0;
995   b->row  = perm;
996   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
997   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
998   b->icol = perm;
999   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
1000   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1001   ierr    = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
1002   b->maxnz = b->nz = ui[am];
1003 
1004   B->info.factor_mallocs    = reallocs;
1005   B->info.fill_ratio_given  = fill;
1006   if (ai[am] != 0) {
1007     B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
1008   } else {
1009     B->info.fill_ratio_needed = 0.0;
1010   }
1011   if (perm_identity){
1012     B->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1013     B->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1014     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
1015   } else {
1016     (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
1017   }
1018   PetscFunctionReturn(0);
1019 }
1020 
1021 #undef __FUNCT__
1022 #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqBAIJ"
1023 PetscErrorCode MatCholeskyFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
1024 {
1025   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data;
1026   Mat_SeqSBAIJ       *b;
1027   Mat                B;
1028   PetscErrorCode     ierr;
1029   PetscTruth         perm_identity;
1030   PetscReal          fill = info->fill;
1031   const PetscInt     *rip;
1032   PetscInt           i,mbs=a->mbs,bs=A->rmap->bs,*ai=a->i,*aj=a->j,reallocs=0,prow;
1033   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
1034   PetscInt           nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr;
1035   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
1036   PetscBT            lnkbt;
1037 
1038   PetscFunctionBegin;
1039   if (bs > 1) { /* convert to seqsbaij */
1040     if (!a->sbaijMat){
1041       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
1042     }
1043     (fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqSBAIJ; /* undue the change made in MatGetFactor_seqbaij_petsc */
1044     ierr = MatCholeskyFactorSymbolic(fact,a->sbaijMat,perm,info);CHKERRQ(ierr);
1045     PetscFunctionReturn(0);
1046   }
1047 
1048   /* check whether perm is the identity mapping */
1049   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
1050   if (!perm_identity) SETERRQ(PETSC_ERR_SUP,"Matrix reordering is not supported");
1051   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
1052 
1053   /* initialization */
1054   ierr  = PetscMalloc((mbs+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
1055   ui[0] = 0;
1056 
1057   /* jl: linked list for storing indices of the pivot rows
1058      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
1059   ierr = PetscMalloc4(mbs,PetscInt*,&ui_ptr,mbs,PetscInt,&il,mbs,PetscInt,&jl,mbs,PetscInt,&cols);CHKERRQ(ierr);
1060   for (i=0; i<mbs; i++){
1061     jl[i] = mbs; il[i] = 0;
1062   }
1063 
1064   /* create and initialize a linked list for storing column indices of the active row k */
1065   nlnk = mbs + 1;
1066   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
1067 
1068   /* initial FreeSpace size is fill*(ai[mbs]+1) */
1069   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[mbs]+1)),&free_space);CHKERRQ(ierr);
1070   current_space = free_space;
1071 
1072   for (k=0; k<mbs; k++){  /* for each active row k */
1073     /* initialize lnk by the column indices of row rip[k] of A */
1074     nzk   = 0;
1075     ncols = ai[rip[k]+1] - ai[rip[k]];
1076     ncols_upper = 0;
1077     for (j=0; j<ncols; j++){
1078       i = rip[*(aj + ai[rip[k]] + j)];
1079       if (i >= k){ /* only take upper triangular entry */
1080         cols[ncols_upper] = i;
1081         ncols_upper++;
1082       }
1083     }
1084     ierr = PetscLLAdd(ncols_upper,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
1085     nzk += nlnk;
1086 
1087     /* update lnk by computing fill-in for each pivot row to be merged in */
1088     prow = jl[k]; /* 1st pivot row */
1089 
1090     while (prow < k){
1091       nextprow = jl[prow];
1092       /* merge prow into k-th row */
1093       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
1094       jmax = ui[prow+1];
1095       ncols = jmax-jmin;
1096       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
1097       ierr = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
1098       nzk += nlnk;
1099 
1100       /* update il and jl for prow */
1101       if (jmin < jmax){
1102         il[prow] = jmin;
1103         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
1104       }
1105       prow = nextprow;
1106     }
1107 
1108     /* if free space is not available, make more free space */
1109     if (current_space->local_remaining<nzk) {
1110       i = mbs - k + 1; /* num of unfactored rows */
1111       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
1112       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
1113       reallocs++;
1114     }
1115 
1116     /* copy data into free space, then initialize lnk */
1117     ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
1118 
1119     /* add the k-th row into il and jl */
1120     if (nzk-1 > 0){
1121       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
1122       jl[k] = jl[i]; jl[i] = k;
1123       il[k] = ui[k] + 1;
1124     }
1125     ui_ptr[k] = current_space->array;
1126     current_space->array           += nzk;
1127     current_space->local_used      += nzk;
1128     current_space->local_remaining -= nzk;
1129 
1130     ui[k+1] = ui[k] + nzk;
1131   }
1132 
1133 #if defined(PETSC_USE_INFO)
1134   if (ai[mbs] != 0) {
1135     PetscReal af = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
1136     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
1137     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
1138     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
1139   } else {
1140     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
1141   }
1142 #endif
1143 
1144   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
1145   ierr = PetscFree4(ui_ptr,il,jl,cols);CHKERRQ(ierr);
1146 
1147   /* destroy list of free space and other temporary array(s) */
1148   ierr = PetscMalloc((ui[mbs]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
1149   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
1150   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1151 
1152   /* put together the new matrix in MATSEQSBAIJ format */
1153   B    = fact;
1154   ierr = MatSeqSBAIJSetPreallocation(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
1155 
1156   b = (Mat_SeqSBAIJ*)B->data;
1157   b->singlemalloc = PETSC_FALSE;
1158   b->free_a       = PETSC_TRUE;
1159   b->free_ij      = PETSC_TRUE;
1160   ierr = PetscMalloc((ui[mbs]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
1161   b->j    = uj;
1162   b->i    = ui;
1163   b->diag = 0;
1164   b->ilen = 0;
1165   b->imax = 0;
1166   b->row  = perm;
1167   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
1168   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
1169   b->icol = perm;
1170   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
1171   ierr    = PetscMalloc((mbs+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1172   ierr    = PetscLogObjectMemory(B,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
1173   b->maxnz = b->nz = ui[mbs];
1174 
1175   B->info.factor_mallocs    = reallocs;
1176   B->info.fill_ratio_given  = fill;
1177   if (ai[mbs] != 0) {
1178     B->info.fill_ratio_needed = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
1179   } else {
1180     B->info.fill_ratio_needed = 0.0;
1181   }
1182   if (perm_identity){
1183     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
1184   } else {
1185     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
1186   }
1187   PetscFunctionReturn(0);
1188 }
1189 
1190 /* --------------------------------------------------------- */
1191 #undef __FUNCT__
1192 #define __FUNCT__ "MatSolve_SeqBAIJ_N_NaturalOrdering_newdatastruct"
1193 PetscErrorCode MatSolve_SeqBAIJ_N_NaturalOrdering_newdatastruct(Mat A,Vec bb,Vec xx)
1194 {
1195   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ *)A->data;
1196   PetscErrorCode ierr;
1197   const PetscInt *ai=a->i,*aj=a->j,*vi;
1198   PetscInt       i,k,n=a->mbs;
1199   PetscInt       nz,bs=A->rmap->bs,bs2=a->bs2,*adiag=a->diag;
1200   MatScalar      *aa=a->a,*v;
1201   PetscScalar    *x,*b,*s,*t,*ls;
1202 
1203   PetscFunctionBegin;
1204   /* printf("MatSolve_SeqBAIJ_NaturalOrdering_iludt..bs %d\n",bs); */
1205   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
1206   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1207   t  = a->solve_work;
1208 
1209   /* forward solve the lower triangular */
1210   ierr = PetscMemcpy(t,b,bs*sizeof(PetscScalar));CHKERRQ(ierr); /* copy 1st block of b to t */
1211 
1212   for (i=1; i<n; i++) {
1213     v   = aa + bs2*ai[i];
1214     vi  = aj + ai[i];
1215     nz = ai[i+1] - ai[i];
1216     s = t + bs*i;
1217     ierr = PetscMemcpy(s,b+bs*i,bs*sizeof(PetscScalar));CHKERRQ(ierr); /* copy i_th block of b to t */
1218     for(k=0;k<nz;k++){
1219       Kernel_v_gets_v_minus_A_times_w(bs,s,v,t+bs*vi[k]);
1220       v += bs2;
1221     }
1222   }
1223 
1224   /* backward solve the upper triangular */
1225   ls = a->solve_work + A->cmap->n;
1226   for (i=n-1; i>=0; i--){
1227     v  = aa + bs2*ai[2*n-i];
1228     vi = aj + ai[2*n-i];
1229     nz = ai[2*n-i +1] - ai[2*n-i]-1;
1230     ierr = PetscMemcpy(ls,t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr);
1231     for(k=0;k<nz;k++){
1232       Kernel_v_gets_v_minus_A_times_w(bs,ls,v,t+bs*vi[k]);
1233       v += bs2;
1234     }
1235     Kernel_w_gets_A_times_v(bs,ls,aa+bs2*adiag[i],t+i*bs); /* *inv(diagonal[i]) */
1236     ierr = PetscMemcpy(x+i*bs,t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr);
1237   }
1238 
1239   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
1240   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1241   ierr = PetscLogFlops(2.0*(a->bs2)*(a->nz) - A->rmap->bs*A->cmap->n);CHKERRQ(ierr);
1242   PetscFunctionReturn(0);
1243 }
1244 
1245 #undef __FUNCT__
1246 #define __FUNCT__ "MatSolve_SeqBAIJ_N_NaturalOrdering_newdatastruct_v2"
1247 PetscErrorCode MatSolve_SeqBAIJ_N_NaturalOrdering_newdatastruct_v2(Mat A,Vec bb,Vec xx)
1248 {
1249   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ *)A->data;
1250   PetscErrorCode ierr;
1251   const PetscInt *ai=a->i,*aj=a->j,*adiag=a->diag,*vi;
1252   PetscInt       i,k,n=a->mbs;
1253   PetscInt       nz,bs=A->rmap->bs,bs2=a->bs2;
1254   MatScalar      *aa=a->a,*v;
1255   PetscScalar    *x,*b,*s,*t,*ls;
1256 
1257   PetscFunctionBegin;
1258   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
1259   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1260   t  = a->solve_work;
1261 
1262   /* forward solve the lower triangular */
1263   ierr = PetscMemcpy(t,b,bs*sizeof(PetscScalar));CHKERRQ(ierr); /* copy 1st block of b to t */
1264 
1265   for (i=1; i<n; i++) {
1266     v   = aa + bs2*ai[i];
1267     vi  = aj + ai[i];
1268     nz = ai[i+1] - ai[i];
1269     s = t + bs*i;
1270     ierr = PetscMemcpy(s,b+bs*i,bs*sizeof(PetscScalar));CHKERRQ(ierr); /* copy i_th block of b to t */
1271     for(k=0;k<nz;k++){
1272       Kernel_v_gets_v_minus_A_times_w(bs,s,v,t+bs*vi[k]);
1273       v += bs2;
1274     }
1275   }
1276 
1277   /* backward solve the upper triangular */
1278   ls = a->solve_work + A->cmap->n;
1279   for (i=n-1; i>=0; i--){
1280     v  = aa + bs2*(adiag[i+1]+1);
1281     vi = aj + adiag[i+1]+1;
1282     nz = adiag[i] - adiag[i+1]-1;
1283     ierr = PetscMemcpy(ls,t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr);
1284     for(k=0;k<nz;k++){
1285       Kernel_v_gets_v_minus_A_times_w(bs,ls,v,t+bs*vi[k]);
1286       v += bs2;
1287     }
1288     Kernel_w_gets_A_times_v(bs,ls,aa+bs2*adiag[i],t+i*bs); /* *inv(diagonal[i]) */
1289     ierr = PetscMemcpy(x+i*bs,t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr);
1290   }
1291 
1292   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
1293   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1294   ierr = PetscLogFlops(2.0*(a->bs2)*(a->nz) - A->rmap->bs*A->cmap->n);CHKERRQ(ierr);
1295   PetscFunctionReturn(0);
1296 }
1297 
1298 #undef __FUNCT__
1299 #define __FUNCT__ "MatSolve_SeqBAIJ_N_newdatastruct"
1300 PetscErrorCode MatSolve_SeqBAIJ_N_newdatastruct(Mat A,Vec bb,Vec xx)
1301 {
1302   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ *)A->data;
1303   IS             iscol=a->col,isrow=a->row;
1304   PetscErrorCode ierr;
1305   const PetscInt *r,*c,*rout,*cout,*ai=a->i,*aj=a->j,*vi;
1306   PetscInt       i,m,n=a->mbs;
1307   PetscInt       nz,bs=A->rmap->bs,bs2=a->bs2,k;
1308   MatScalar      *aa=a->a,*v;
1309   PetscScalar    *x,*b,*s,*t,*ls;
1310 
1311   PetscFunctionBegin;
1312   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
1313   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1314   t  = a->solve_work;
1315 
1316   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
1317   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
1318 
1319   /* forward solve the lower triangular */
1320   ierr = PetscMemcpy(t,b+bs*r[0],bs*sizeof(PetscScalar));CHKERRQ(ierr);
1321   for (i=1; i<n; i++) {
1322     v   = aa + bs2*ai[i];
1323     vi  = aj + ai[i];
1324     nz = ai[i+1] - ai[i];
1325     s = t + bs*i;
1326     ierr = PetscMemcpy(s,b+bs*r[i],bs*sizeof(PetscScalar));CHKERRQ(ierr);
1327     for(m=0;m<nz;m++){
1328       Kernel_v_gets_v_minus_A_times_w(bs,s,v,t+bs*vi[m]);
1329       v += bs2;
1330     }
1331   }
1332 
1333   /* backward solve the upper triangular */
1334   ls = a->solve_work + A->cmap->n;
1335   for (i=n-1; i>=0; i--){
1336     k  = 2*n-i;
1337     v  = aa + bs2*ai[k];
1338     vi = aj + ai[k];
1339     nz = ai[k+1] - ai[k] - 1;
1340     ierr = PetscMemcpy(ls,t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr);
1341     for(m=0;m<nz;m++){
1342       Kernel_v_gets_v_minus_A_times_w(bs,ls,v,t+bs*vi[m]);
1343       v += bs2;
1344     }
1345     Kernel_w_gets_A_times_v(bs,ls,v,t+i*bs); /* *inv(diagonal[i]) */
1346     ierr = PetscMemcpy(x + bs*c[i],t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr);
1347   }
1348   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
1349   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
1350   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
1351   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1352   ierr = PetscLogFlops(2.0*(a->bs2)*(a->nz) - A->rmap->bs*A->cmap->n);CHKERRQ(ierr);
1353   PetscFunctionReturn(0);
1354 }
1355 
1356 #undef __FUNCT__
1357 #define __FUNCT__ "MatSolve_SeqBAIJ_N_newdatastruct_v2"
1358 PetscErrorCode MatSolve_SeqBAIJ_N_newdatastruct_v2(Mat A,Vec bb,Vec xx)
1359 {
1360   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ *)A->data;
1361   IS             iscol=a->col,isrow=a->row;
1362   PetscErrorCode ierr;
1363   const PetscInt *r,*c,*rout,*cout,*ai=a->i,*aj=a->j,*adiag=a->diag,*vi;
1364   PetscInt       i,m,n=a->mbs;
1365   PetscInt       nz,bs=A->rmap->bs,bs2=a->bs2;
1366   MatScalar      *aa=a->a,*v;
1367   PetscScalar    *x,*b,*s,*t,*ls;
1368 
1369   PetscFunctionBegin;
1370   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
1371   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1372   t  = a->solve_work;
1373 
1374   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
1375   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
1376 
1377   /* forward solve the lower triangular */
1378   ierr = PetscMemcpy(t,b+bs*r[0],bs*sizeof(PetscScalar));CHKERRQ(ierr);
1379   for (i=1; i<n; i++) {
1380     v   = aa + bs2*ai[i];
1381     vi  = aj + ai[i];
1382     nz = ai[i+1] - ai[i];
1383     s = t + bs*i;
1384     ierr = PetscMemcpy(s,b+bs*r[i],bs*sizeof(PetscScalar));CHKERRQ(ierr);
1385     for(m=0;m<nz;m++){
1386       Kernel_v_gets_v_minus_A_times_w(bs,s,v,t+bs*vi[m]);
1387       v += bs2;
1388     }
1389   }
1390 
1391   /* backward solve the upper triangular */
1392   ls = a->solve_work + A->cmap->n;
1393   for (i=n-1; i>=0; i--){
1394     v  = aa + bs2*(adiag[i+1]+1);
1395     vi = aj + adiag[i+1]+1;
1396     nz = adiag[i] - adiag[i+1] - 1;
1397     ierr = PetscMemcpy(ls,t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr);
1398     for(m=0;m<nz;m++){
1399       Kernel_v_gets_v_minus_A_times_w(bs,ls,v,t+bs*vi[m]);
1400       v += bs2;
1401     }
1402     Kernel_w_gets_A_times_v(bs,ls,v,t+i*bs); /* *inv(diagonal[i]) */
1403     ierr = PetscMemcpy(x + bs*c[i],t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr);
1404   }
1405   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
1406   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
1407   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
1408   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1409   ierr = PetscLogFlops(2.0*(a->bs2)*(a->nz) - A->rmap->bs*A->cmap->n);CHKERRQ(ierr);
1410   PetscFunctionReturn(0);
1411 }
1412 
1413 #undef __FUNCT__
1414 #define __FUNCT__ "BlockAbs_privat"
1415 PetscErrorCode BlockAbs_private(PetscInt nbs,PetscInt bs2,PetscScalar *blockarray,PetscReal *absarray)
1416 {
1417   PetscErrorCode     ierr;
1418   PetscInt           i,j;
1419   PetscFunctionBegin;
1420   ierr = PetscMemzero(absarray,(nbs+1)*sizeof(PetscReal));CHKERRQ(ierr);
1421   for (i=0; i<nbs; i++){
1422     for (j=0; j<bs2; j++){
1423       if (absarray[i] < PetscAbsScalar(blockarray[i*nbs+j])) absarray[i] = PetscAbsScalar(blockarray[i*nbs+j]);
1424     }
1425   }
1426   PetscFunctionReturn(0);
1427 }
1428 
1429 #undef __FUNCT__
1430 #define __FUNCT__ "MatILUDTFactor_SeqBAIJ"
1431 PetscErrorCode MatILUDTFactor_SeqBAIJ(Mat A,IS isrow,IS iscol,const MatFactorInfo *info,Mat *fact)
1432 {
1433   Mat                B = *fact;
1434   Mat_SeqBAIJ        *a=(Mat_SeqBAIJ*)A->data,*b;
1435   IS                 isicol;
1436   PetscErrorCode     ierr;
1437   const PetscInt     *r,*ic;
1438   PetscInt           i,mbs=a->mbs,bs=A->rmap->bs,bs2=a->bs2,*ai=a->i,*aj=a->j,*ajtmp,*adiag;
1439   PetscInt           *bi,*bj,*bdiag;
1440 
1441   PetscInt           row,nzi,nzi_bl,nzi_bu,*im,dtcount,nzi_al,nzi_au;
1442   PetscInt           nlnk,*lnk;
1443   PetscBT            lnkbt;
1444   PetscTruth         row_identity,icol_identity,both_identity;
1445   MatScalar          *aatmp,*pv,*batmp,*ba,*rtmp,*pc,*multiplier,*vtmp;
1446   const PetscInt     *ics;
1447   PetscInt           j,nz,*pj,*bjtmp,k,ncut,*jtmp;
1448 
1449   PetscReal          dt=info->dt; /* shift=info->shiftinblocks; */
1450   PetscInt           nnz_max;
1451   PetscTruth         missing;
1452   PetscReal          *vtmp_abs;
1453   MatScalar          *v_work;
1454   PetscInt           *v_pivots;
1455 
1456   PetscFunctionBegin;
1457   /* ------- symbolic factorization, can be reused ---------*/
1458   ierr = MatMissingDiagonal(A,&missing,&i);CHKERRQ(ierr);
1459   if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i);
1460   adiag=a->diag;
1461 
1462   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
1463 
1464   /* bdiag is location of diagonal in factor */
1465   ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
1466 
1467   /* allocate row pointers bi */
1468   ierr = PetscMalloc((2*mbs+2)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
1469 
1470   /* allocate bj and ba; max num of nonzero entries is (ai[n]+2*n*dtcount+2) */
1471   dtcount = (PetscInt)info->dtcount;
1472   if (dtcount > mbs-1) dtcount = mbs-1;
1473   nnz_max  = ai[mbs]+2*mbs*dtcount +2;
1474   /* printf("MatILUDTFactor_SeqBAIJ, bs %d, ai[mbs] %d, nnz_max  %d, dtcount %d\n",bs,ai[mbs],nnz_max,dtcount); */
1475   ierr = PetscMalloc(nnz_max*sizeof(PetscInt),&bj);CHKERRQ(ierr);
1476   nnz_max = nnz_max*bs2;
1477   ierr = PetscMalloc(nnz_max*sizeof(MatScalar),&ba);CHKERRQ(ierr);
1478 
1479   /* put together the new matrix */
1480   ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
1481   ierr = PetscLogObjectParent(B,isicol);CHKERRQ(ierr);
1482   b    = (Mat_SeqBAIJ*)(B)->data;
1483   b->free_a       = PETSC_TRUE;
1484   b->free_ij      = PETSC_TRUE;
1485   b->singlemalloc = PETSC_FALSE;
1486   b->a          = ba;
1487   b->j          = bj;
1488   b->i          = bi;
1489   b->diag       = bdiag;
1490   b->ilen       = 0;
1491   b->imax       = 0;
1492   b->row        = isrow;
1493   b->col        = iscol;
1494   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1495   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
1496   b->icol       = isicol;
1497   ierr = PetscMalloc((bs*(mbs+1))*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1498 
1499   ierr = PetscLogObjectMemory(B,nnz_max*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
1500   b->maxnz = nnz_max/bs2;
1501 
1502   (B)->factor                = MAT_FACTOR_ILUDT;
1503   (B)->info.factor_mallocs   = 0;
1504   (B)->info.fill_ratio_given = ((PetscReal)nnz_max)/((PetscReal)(ai[mbs]*bs2));
1505   CHKMEMQ;
1506   /* ------- end of symbolic factorization ---------*/
1507   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
1508   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
1509   ics  = ic;
1510 
1511   /* linked list for storing column indices of the active row */
1512   nlnk = mbs + 1;
1513   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
1514 
1515   /* im: used by PetscLLAddSortedLU(); jtmp: working array for column indices of active row */
1516   ierr = PetscMalloc2(mbs,PetscInt,&im,mbs,PetscInt,&jtmp);CHKERRQ(ierr);
1517   /* rtmp, vtmp: working arrays for sparse and contiguous row entries of active row */
1518   ierr = PetscMalloc2(mbs*bs2,MatScalar,&rtmp,mbs*bs2,MatScalar,&vtmp);CHKERRQ(ierr);
1519   ierr = PetscMalloc((mbs+1)*sizeof(PetscReal),&vtmp_abs);CHKERRQ(ierr);
1520   ierr = PetscMalloc3(bs,MatScalar,&v_work,bs2,MatScalar,&multiplier,bs,PetscInt,&v_pivots);CHKERRQ(ierr);
1521 
1522   bi[0]    = 0;
1523   bdiag[0] = (nnz_max/bs2)-1; /* location of diagonal in factor B */
1524   bi[2*mbs+1] = bdiag[0]+1; /* endof bj and ba array */
1525   for (i=0; i<mbs; i++) {
1526     /* copy initial fill into linked list */
1527     nzi = 0; /* nonzeros for active row i */
1528     nzi = ai[r[i]+1] - ai[r[i]];
1529     if (!nzi) SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i);
1530     nzi_al = adiag[r[i]] - ai[r[i]];
1531     nzi_au = ai[r[i]+1] - adiag[r[i]] -1;
1532     /* printf("row %d, nzi_al/au %d %d\n",i,nzi_al,nzi_au); */
1533 
1534     /* load in initial unfactored row */
1535     ajtmp = aj + ai[r[i]];
1536     ierr = PetscLLAddPerm(nzi,ajtmp,ic,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
1537     ierr = PetscMemzero(rtmp,mbs*bs2*sizeof(PetscScalar));CHKERRQ(ierr);
1538     aatmp = a->a + bs2*ai[r[i]];
1539     for (j=0; j<nzi; j++) {
1540       ierr = PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],aatmp+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1541     }
1542 
1543     /* add pivot rows into linked list */
1544     row = lnk[mbs];
1545     while (row < i) {
1546       nzi_bl = bi[row+1] - bi[row] + 1;
1547       bjtmp = bj + bdiag[row+1]+1; /* points to 1st column next to the diagonal in U */
1548       ierr  = PetscLLAddSortedLU(bjtmp,row,nlnk,lnk,lnkbt,i,nzi_bl,im);CHKERRQ(ierr);
1549       nzi  += nlnk;
1550       row   = lnk[row];
1551     }
1552 
1553     /* copy data from lnk into jtmp, then initialize lnk */
1554     ierr = PetscLLClean(mbs,mbs,nzi,lnk,jtmp,lnkbt);CHKERRQ(ierr);
1555 
1556     /* numerical factorization */
1557     bjtmp = jtmp;
1558     row   = *bjtmp++; /* 1st pivot row */
1559 
1560     while  (row < i) {
1561       pc = rtmp + bs2*row;
1562       pv = ba + bs2*bdiag[row]; /* inv(diag) of the pivot row */
1563       Kernel_A_gets_A_times_B(bs,pc,pv,multiplier); /* pc= multiplier = pc*inv(diag[row]) */
1564       ierr = BlockAbs_private(1,bs2,pc,vtmp_abs);CHKERRQ(ierr);
1565       if (vtmp_abs[0] > dt){ /* apply tolerance dropping rule */
1566         pj         = bj + bdiag[row+1] + 1; /* point to 1st entry of U(row,:) */
1567         pv         = ba + bs2*(bdiag[row+1] + 1);
1568         nz         = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:), excluding diagonal */
1569         for (j=0; j<nz; j++){
1570           Kernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j);
1571         }
1572         /* ierr = PetscLogFlops(bslog*(nz+1.0)-bs);CHKERRQ(ierr); */
1573       }
1574       row = *bjtmp++;
1575     }
1576 
1577     /* copy sparse rtmp into contiguous vtmp; separate L and U part */
1578     nzi_bl = 0; j = 0;
1579     while (jtmp[j] < i){ /* L-part. Note: jtmp is sorted */
1580       ierr = PetscMemcpy(vtmp+bs2*j,rtmp+bs2*jtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
1581       nzi_bl++; j++;
1582     }
1583     nzi_bu = nzi - nzi_bl -1;
1584     /* printf("nzi %d, nzi_bl %d, nzi_bu %d\n",nzi,nzi_bl,nzi_bu); */
1585 
1586     while (j < nzi){ /* U-part */
1587       ierr = PetscMemcpy(vtmp+bs2*j,rtmp+bs2*jtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
1588       /*
1589       printf(" col %d: ",jtmp[j]);
1590       for (j1=0; j1<bs2; j1++) printf(" %g",*(vtmp+bs2*j+j1));
1591       printf(" \n");
1592       */
1593       j++;
1594     }
1595 
1596     ierr = BlockAbs_private(nzi,bs2,vtmp,vtmp_abs);CHKERRQ(ierr);
1597     /*
1598     printf(" row %d, nzi %d, vtmp_abs\n",i,nzi);
1599     for (j1=0; j1<nzi; j1++) printf(" (%d %g),",jtmp[j1],vtmp_abs[j1]);
1600     printf(" \n");
1601     */
1602     bjtmp = bj + bi[i];
1603     batmp = ba + bs2*bi[i];
1604     /* apply level dropping rule to L part */
1605     ncut = nzi_al + dtcount;
1606     if (ncut < nzi_bl){
1607       ierr = PetscSortSplitReal(ncut,nzi_bl,vtmp_abs,jtmp);CHKERRQ(ierr);
1608       ierr = PetscSortIntWithScalarArray(ncut,jtmp,vtmp);CHKERRQ(ierr);
1609     } else {
1610       ncut = nzi_bl;
1611     }
1612     for (j=0; j<ncut; j++){
1613       bjtmp[j] = jtmp[j];
1614       ierr = PetscMemcpy(batmp+bs2*j,rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
1615       /*
1616       printf(" col %d: ",bjtmp[j]);
1617       for (j1=0; j1<bs2; j1++) printf(" %g,",*(batmp+bs2*j+j1));
1618       printf("\n");
1619       */
1620     }
1621     bi[i+1] = bi[i] + ncut;
1622     nzi = ncut + 1;
1623 
1624     /* apply level dropping rule to U part */
1625     ncut = nzi_au + dtcount;
1626     if (ncut < nzi_bu){
1627       ierr = PetscSortSplitReal(ncut,nzi_bu,vtmp_abs+nzi_bl+1,jtmp+nzi_bl+1);CHKERRQ(ierr);
1628       ierr = PetscSortIntWithScalarArray(ncut,jtmp+nzi_bl+1,vtmp+nzi_bl+1);CHKERRQ(ierr);
1629     } else {
1630       ncut = nzi_bu;
1631     }
1632     nzi += ncut;
1633 
1634     /* mark bdiagonal */
1635     bdiag[i+1]    = bdiag[i] - (ncut + 1);
1636     bi[2*mbs - i] = bi[2*mbs - i +1] - (ncut + 1);
1637 
1638     bjtmp = bj + bdiag[i];
1639     batmp = ba + bs2*bdiag[i];
1640     ierr = PetscMemcpy(batmp,rtmp+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1641     *bjtmp = i;
1642     /*
1643     printf(" diag %d: ",*bjtmp);
1644     for (j=0; j<bs2; j++){
1645       printf(" %g,",batmp[j]);
1646     }
1647     printf("\n");
1648     */
1649     bjtmp = bj + bdiag[i+1]+1;
1650     batmp = ba + (bdiag[i+1]+1)*bs2;
1651 
1652     for (k=0; k<ncut; k++){
1653       bjtmp[k] = jtmp[nzi_bl+1+k];
1654       ierr = PetscMemcpy(batmp+bs2*k,rtmp+bs2*bjtmp[k],bs2*sizeof(MatScalar));CHKERRQ(ierr);
1655       /*
1656       printf(" col %d:",bjtmp[k]);
1657       for (j1=0; j1<bs2; j1++) printf(" %g,",*(batmp+bs2*k+j1));
1658       printf("\n");
1659       */
1660     }
1661 
1662     im[i] = nzi; /* used by PetscLLAddSortedLU() */
1663 
1664     /* invert diagonal block for simplier triangular solves - add shift??? */
1665     batmp = ba + bs2*bdiag[i];
1666     ierr = Kernel_A_gets_inverse_A(bs,batmp,v_pivots,v_work);CHKERRQ(ierr);
1667   } /* for (i=0; i<mbs; i++) */
1668   ierr = PetscFree3(v_work,multiplier,v_pivots);CHKERRQ(ierr);
1669 
1670   /* printf("end of L %d, beginning of U %d\n",bi[mbs],bdiag[mbs]); */
1671   if (bi[mbs] >= bdiag[mbs]) SETERRQ2(PETSC_ERR_ARG_SIZ,"end of L array %d cannot >= the beginning of U array %d",bi[mbs],bdiag[mbs]);
1672 
1673   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
1674   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
1675 
1676   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1677 
1678   ierr = PetscFree2(im,jtmp);CHKERRQ(ierr);
1679   ierr = PetscFree2(rtmp,vtmp);CHKERRQ(ierr);
1680 
1681   ierr = PetscLogFlops(bs2*B->cmap->n);CHKERRQ(ierr);
1682   b->maxnz = b->nz = bi[mbs] + bdiag[0] - bdiag[mbs];
1683 
1684   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
1685   ierr = ISIdentity(isicol,&icol_identity);CHKERRQ(ierr);
1686   both_identity = (PetscTruth) (row_identity && icol_identity);
1687   if (row_identity && icol_identity) {
1688     B->ops->solve = MatSolve_SeqBAIJ_N_NaturalOrdering_newdatastruct;
1689   } else {
1690     B->ops->solve = MatSolve_SeqBAIJ_N_newdatastruct;
1691   }
1692 
1693   B->ops->solveadd          = 0;
1694   B->ops->solvetranspose    = 0;
1695   B->ops->solvetransposeadd = 0;
1696   B->ops->matsolve          = 0;
1697   B->assembled              = PETSC_TRUE;
1698   B->preallocated           = PETSC_TRUE;
1699   PetscFunctionReturn(0);
1700 }
1701