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