xref: /petsc/src/mat/impls/baij/seq/baijfact.c (revision bbd65245a52e23ed6d5d64ba4c4d98fbdb9919b2)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2be1d678aSKris Buschelman 
35d517e7eSBarry Smith /*
4ec3a8b7bSBarry Smith     Factorization code for BAIJ format.
55d517e7eSBarry Smith */
67c4f633dSBarry Smith #include "../src/mat/impls/baij/seq/baij.h"
7c60f0209SBarry Smith #include "../src/mat/blockinvert.h"
8ec3a8b7bSBarry Smith 
9b588c5a2SHong Zhang #undef __FUNCT__
104dd39f65SShri Abhyankar #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2"
114dd39f65SShri Abhyankar PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2(Mat B,Mat A,const MatFactorInfo *info)
12b588c5a2SHong Zhang {
13b588c5a2SHong Zhang   Mat             C=B;
14b588c5a2SHong Zhang   Mat_SeqBAIJ     *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ *)C->data;
15b588c5a2SHong Zhang   IS              isrow = b->row,isicol = b->icol;
16b588c5a2SHong Zhang   PetscErrorCode  ierr;
17b588c5a2SHong Zhang   const PetscInt  *r,*ic,*ics;
18*bbd65245SShri Abhyankar   PetscInt        i,j,k,nz,nzL,row,*pj;
19*bbd65245SShri Abhyankar   const PetscInt  n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,bs2=a->bs2;
20*bbd65245SShri Abhyankar   const PetscInt  *ajtmp,*bjtmp,*bdiag=b->diag;
21*bbd65245SShri Abhyankar   MatScalar       *rtmp,*pc,*mwork,*pv;
22*bbd65245SShri Abhyankar   MatScalar       *aa=a->a,*v;
23*bbd65245SShri Abhyankar   PetscInt        flg;
24b588c5a2SHong Zhang   PetscReal       shift = info->shiftinblocks;
25b588c5a2SHong Zhang 
26b588c5a2SHong Zhang   PetscFunctionBegin;
27b588c5a2SHong Zhang   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
28b588c5a2SHong Zhang   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
29b588c5a2SHong Zhang 
304c000e2eSHong Zhang   /* generate work space needed by the factorization */
31fca92195SBarry Smith   ierr = PetscMalloc2(bs2*n,MatScalar,&rtmp,bs2,MatScalar,&mwork);CHKERRQ(ierr);
324c000e2eSHong Zhang   ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr);
33b588c5a2SHong Zhang   ics  = ic;
34b588c5a2SHong Zhang 
35b588c5a2SHong Zhang   for (i=0; i<n; i++){
36b588c5a2SHong Zhang     /* zero rtmp */
37b588c5a2SHong Zhang     /* L part */
38b588c5a2SHong Zhang     nz    = bi[i+1] - bi[i];
39b588c5a2SHong Zhang     bjtmp = bj + bi[i];
40b588c5a2SHong Zhang     for  (j=0; j<nz; j++){
41b588c5a2SHong Zhang       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
42b588c5a2SHong Zhang     }
43b588c5a2SHong Zhang 
44b588c5a2SHong Zhang     /* U part */
450c4413a7SShri Abhyankar     nz = bdiag[i] - bdiag[i+1];
460c4413a7SShri Abhyankar     bjtmp = bj + bdiag[i+1]+1;
470c4413a7SShri Abhyankar     for  (j=0; j<nz; j++){
480c4413a7SShri Abhyankar       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
490c4413a7SShri Abhyankar     }
500c4413a7SShri Abhyankar 
510c4413a7SShri Abhyankar     /* load in initial (unfactored row) */
520c4413a7SShri Abhyankar     nz    = ai[r[i]+1] - ai[r[i]];
530c4413a7SShri Abhyankar     ajtmp = aj + ai[r[i]];
540c4413a7SShri Abhyankar     v     = aa + bs2*ai[r[i]];
550c4413a7SShri Abhyankar     for (j=0; j<nz; j++) {
560c4413a7SShri Abhyankar       ierr = PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr);
570c4413a7SShri Abhyankar     }
580c4413a7SShri Abhyankar 
590c4413a7SShri Abhyankar     /* elimination */
600c4413a7SShri Abhyankar     bjtmp = bj + bi[i];
610c4413a7SShri Abhyankar     nzL   = bi[i+1] - bi[i];
620c4413a7SShri Abhyankar     for(k=0;k < nzL;k++) {
630c4413a7SShri Abhyankar       row = bjtmp[k];
640c4413a7SShri Abhyankar       pc = rtmp + bs2*row;
650c4413a7SShri Abhyankar       for (flg=0,j=0; j<bs2; j++) { if (pc[j]!=0.0) { flg = 1; break; }}
660c4413a7SShri Abhyankar       if (flg) {
670c4413a7SShri Abhyankar         pv = b->a + bs2*bdiag[row];
680c4413a7SShri Abhyankar         /* Kernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
690c4413a7SShri Abhyankar         ierr = Kernel_A_gets_A_times_B_2(pc,pv,mwork);CHKERRQ(ierr);
700c4413a7SShri Abhyankar 
710c4413a7SShri Abhyankar         pj = b->j + bdiag[row+1]+1; /* begining of U(row,:) */
720c4413a7SShri Abhyankar         pv = b->a + bs2*(bdiag[row+1]+1);
730c4413a7SShri Abhyankar         nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries inU(row,:), excluding diag */
740c4413a7SShri Abhyankar         for (j=0; j<nz; j++) {
750c4413a7SShri Abhyankar           /* Kernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
760c4413a7SShri Abhyankar           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
770c4413a7SShri Abhyankar           v    = rtmp + 4*pj[j];
780c4413a7SShri Abhyankar           ierr = Kernel_A_gets_A_minus_B_times_C_2(v,pc,pv);CHKERRQ(ierr);
790c4413a7SShri Abhyankar           pv  += 4;
800c4413a7SShri Abhyankar         }
810c4413a7SShri Abhyankar         ierr = PetscLogFlops(16*nz+12);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
820c4413a7SShri Abhyankar       }
830c4413a7SShri Abhyankar     }
840c4413a7SShri Abhyankar 
850c4413a7SShri Abhyankar     /* finished row so stick it into b->a */
860c4413a7SShri Abhyankar     /* L part */
870c4413a7SShri Abhyankar     pv   = b->a + bs2*bi[i] ;
880c4413a7SShri Abhyankar     pj   = b->j + bi[i] ;
890c4413a7SShri Abhyankar     nz   = bi[i+1] - bi[i];
900c4413a7SShri Abhyankar     for (j=0; j<nz; j++) {
910c4413a7SShri Abhyankar       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
920c4413a7SShri Abhyankar     }
930c4413a7SShri Abhyankar 
940c4413a7SShri Abhyankar     /* Mark diagonal and invert diagonal for simplier triangular solves */
950c4413a7SShri Abhyankar     pv   = b->a + bs2*bdiag[i];
960c4413a7SShri Abhyankar     pj   = b->j + bdiag[i];
970c4413a7SShri Abhyankar     ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr);
980c4413a7SShri Abhyankar     /* ierr = Kernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); */
990c4413a7SShri Abhyankar     ierr = Kernel_A_gets_inverse_A_2(pv,shift);CHKERRQ(ierr);
1000c4413a7SShri Abhyankar 
1010c4413a7SShri Abhyankar     /* U part */
1020c4413a7SShri Abhyankar     pv = b->a + bs2*(bdiag[i+1]+1);
1030c4413a7SShri Abhyankar     pj = b->j + bdiag[i+1]+1;
1040c4413a7SShri Abhyankar     nz = bdiag[i] - bdiag[i+1] - 1;
1050c4413a7SShri Abhyankar     for (j=0; j<nz; j++){
1060c4413a7SShri Abhyankar       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
1070c4413a7SShri Abhyankar     }
1080c4413a7SShri Abhyankar   }
1090c4413a7SShri Abhyankar 
110fca92195SBarry Smith   ierr = PetscFree2(rtmp,mwork);CHKERRQ(ierr);
1110c4413a7SShri Abhyankar   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
1120c4413a7SShri Abhyankar   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
1134dd39f65SShri Abhyankar   C->ops->solve          = MatSolve_SeqBAIJ_2;
1144dd39f65SShri Abhyankar   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2;
1150c4413a7SShri Abhyankar 
1160c4413a7SShri Abhyankar   C->assembled = PETSC_TRUE;
117766f9fbaSBarry Smith   ierr = PetscLogFlops(1.333333333333*2*2*2*n);CHKERRQ(ierr); /* from inverting diagonal blocks */
1180c4413a7SShri Abhyankar   PetscFunctionReturn(0);
1190c4413a7SShri Abhyankar }
1200c4413a7SShri Abhyankar 
1214c000e2eSHong Zhang #undef __FUNCT__
1224dd39f65SShri Abhyankar #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering"
1234dd39f65SShri Abhyankar PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info)
1244c000e2eSHong Zhang {
1254c000e2eSHong Zhang   Mat             C=B;
1264c000e2eSHong Zhang   Mat_SeqBAIJ     *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ *)C->data;
1274c000e2eSHong Zhang   PetscErrorCode  ierr;
128*bbd65245SShri Abhyankar   PetscInt        i,j,k,nz,nzL,row,*pj;
129*bbd65245SShri Abhyankar   const PetscInt  n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,bs2=a->bs2;
130*bbd65245SShri Abhyankar   const PetscInt  *ajtmp,*bjtmp,*bdiag=b->diag;
131*bbd65245SShri Abhyankar   MatScalar       *rtmp,*pc,*mwork,*pv;
132*bbd65245SShri Abhyankar   MatScalar       *aa=a->a,*v;
133*bbd65245SShri Abhyankar   PetscInt       flg;
1344c000e2eSHong Zhang   PetscReal      shift = info->shiftinblocks;
1354c000e2eSHong Zhang 
1364c000e2eSHong Zhang   PetscFunctionBegin;
1374c000e2eSHong Zhang   /* generate work space needed by the factorization */
138fca92195SBarry Smith   ierr = PetscMalloc2(bs2*n,MatScalar,&rtmp,bs2,MatScalar,&mwork);CHKERRQ(ierr);
1394c000e2eSHong Zhang   ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr);
1404c000e2eSHong Zhang 
1414c000e2eSHong Zhang   for (i=0; i<n; i++){
1424c000e2eSHong Zhang     /* zero rtmp */
1434c000e2eSHong Zhang     /* L part */
1444c000e2eSHong Zhang     nz    = bi[i+1] - bi[i];
1454c000e2eSHong Zhang     bjtmp = bj + bi[i];
1464c000e2eSHong Zhang     for  (j=0; j<nz; j++){
1474c000e2eSHong Zhang       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
1484c000e2eSHong Zhang     }
1494c000e2eSHong Zhang 
1504c000e2eSHong Zhang     /* U part */
151b2b2dd24SShri Abhyankar     nz = bdiag[i] - bdiag[i+1];
152b2b2dd24SShri Abhyankar     bjtmp = bj + bdiag[i+1]+1;
153b2b2dd24SShri Abhyankar     for  (j=0; j<nz; j++){
154b2b2dd24SShri Abhyankar       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
155b2b2dd24SShri Abhyankar     }
156b2b2dd24SShri Abhyankar 
157b2b2dd24SShri Abhyankar     /* load in initial (unfactored row) */
158b2b2dd24SShri Abhyankar     nz    = ai[i+1] - ai[i];
159b2b2dd24SShri Abhyankar     ajtmp = aj + ai[i];
160b2b2dd24SShri Abhyankar     v     = aa + bs2*ai[i];
161b2b2dd24SShri Abhyankar     for (j=0; j<nz; j++) {
162b2b2dd24SShri Abhyankar       ierr = PetscMemcpy(rtmp+bs2*ajtmp[j],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr);
163b2b2dd24SShri Abhyankar     }
164b2b2dd24SShri Abhyankar 
165b2b2dd24SShri Abhyankar     /* elimination */
166b2b2dd24SShri Abhyankar     bjtmp = bj + bi[i];
167b2b2dd24SShri Abhyankar     nzL   = bi[i+1] - bi[i];
168b2b2dd24SShri Abhyankar     for(k=0;k < nzL;k++) {
169b2b2dd24SShri Abhyankar       row = bjtmp[k];
170b2b2dd24SShri Abhyankar       pc = rtmp + bs2*row;
171b2b2dd24SShri Abhyankar       for (flg=0,j=0; j<bs2; j++) { if (pc[j]!=0.0) { flg = 1; break; }}
172b2b2dd24SShri Abhyankar       if (flg) {
173b2b2dd24SShri Abhyankar         pv = b->a + bs2*bdiag[row];
174b2b2dd24SShri Abhyankar         /* Kernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
175b2b2dd24SShri Abhyankar         ierr = Kernel_A_gets_A_times_B_2(pc,pv,mwork);CHKERRQ(ierr);
176b2b2dd24SShri Abhyankar 
177b2b2dd24SShri Abhyankar         pj = b->j + bdiag[row+1]+1; /* beginning of U(row,:) */
178b2b2dd24SShri Abhyankar 	pv = b->a + bs2*(bdiag[row+1]+1);
179b2b2dd24SShri Abhyankar 	nz = bdiag[row]-bdiag[row+1] - 1; /* num of entries in U(row,:) excluding diag */
180b2b2dd24SShri Abhyankar         for (j=0; j<nz; j++) {
181b2b2dd24SShri Abhyankar           /* Kernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
182b2b2dd24SShri Abhyankar           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
183b2b2dd24SShri Abhyankar           v    = rtmp + 4*pj[j];
184b2b2dd24SShri Abhyankar           ierr = Kernel_A_gets_A_minus_B_times_C_2(v,pc,pv);CHKERRQ(ierr);
185b2b2dd24SShri Abhyankar           pv  += 4;
186b2b2dd24SShri Abhyankar         }
187b2b2dd24SShri Abhyankar         ierr = PetscLogFlops(16*nz+12);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
188b2b2dd24SShri Abhyankar       }
189b2b2dd24SShri Abhyankar     }
190b2b2dd24SShri Abhyankar 
191b2b2dd24SShri Abhyankar     /* finished row so stick it into b->a */
192b2b2dd24SShri Abhyankar     /* L part */
193b2b2dd24SShri Abhyankar     pv   = b->a + bs2*bi[i] ;
194b2b2dd24SShri Abhyankar     pj   = b->j + bi[i] ;
195b2b2dd24SShri Abhyankar     nz   = bi[i+1] - bi[i];
196b2b2dd24SShri Abhyankar     for (j=0; j<nz; j++) {
197b2b2dd24SShri Abhyankar       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
198b2b2dd24SShri Abhyankar     }
199b2b2dd24SShri Abhyankar 
200b2b2dd24SShri Abhyankar     /* Mark diagonal and invert diagonal for simplier triangular solves */
201b2b2dd24SShri Abhyankar     pv   = b->a + bs2*bdiag[i];
202b2b2dd24SShri Abhyankar     pj   = b->j + bdiag[i];
203b2b2dd24SShri Abhyankar     ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr);
204b2b2dd24SShri Abhyankar     /* ierr = Kernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); */
205b2b2dd24SShri Abhyankar     ierr = Kernel_A_gets_inverse_A_2(pv,shift);CHKERRQ(ierr);
206b2b2dd24SShri Abhyankar 
207b2b2dd24SShri Abhyankar     /* U part */
208b2b2dd24SShri Abhyankar     /*
209b2b2dd24SShri Abhyankar     pv = b->a + bs2*bi[2*n-i];
210b2b2dd24SShri Abhyankar     pj = b->j + bi[2*n-i];
211b2b2dd24SShri Abhyankar     nz = bi[2*n-i+1] - bi[2*n-i] - 1;
212b2b2dd24SShri Abhyankar     */
213b2b2dd24SShri Abhyankar     pv = b->a + bs2*(bdiag[i+1]+1);
214b2b2dd24SShri Abhyankar     pj = b->j + bdiag[i+1]+1;
215b2b2dd24SShri Abhyankar     nz = bdiag[i] - bdiag[i+1] - 1;
216b2b2dd24SShri Abhyankar     for (j=0; j<nz; j++){
217b2b2dd24SShri Abhyankar       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
218b2b2dd24SShri Abhyankar     }
219b2b2dd24SShri Abhyankar   }
220fca92195SBarry Smith   ierr = PetscFree2(rtmp,mwork);CHKERRQ(ierr);
221ae3d28f0SHong Zhang 
2224dd39f65SShri Abhyankar   C->ops->solve          = MatSolve_SeqBAIJ_2_NaturalOrdering;
2234dd39f65SShri Abhyankar   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_NaturalOrdering;
224b2b2dd24SShri Abhyankar   C->assembled = PETSC_TRUE;
225766f9fbaSBarry Smith   ierr = PetscLogFlops(1.333333333333*2*2*2*n);CHKERRQ(ierr); /* from inverting diagonal blocks */
226b2b2dd24SShri Abhyankar   PetscFunctionReturn(0);
227b2b2dd24SShri Abhyankar }
228b2b2dd24SShri Abhyankar 
229b2b2dd24SShri Abhyankar #undef __FUNCT__
23006e38f1dSHong Zhang #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2_inplace"
23106e38f1dSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_inplace(Mat B,Mat A,const MatFactorInfo *info)
2324eeb42bcSBarry Smith {
233719d5645SBarry Smith   Mat            C = B;
2344eeb42bcSBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data;
2357cdf2dbbSSatish Balay   IS             isrow = b->row,isicol = b->icol;
2366849ba73SBarry Smith   PetscErrorCode ierr;
2375d0c19d7SBarry Smith   const PetscInt *r,*ic;
2385d0c19d7SBarry Smith   PetscInt       i,j,n = a->mbs,*bi = b->i,*bj = b->j;
239690b6cddSBarry Smith   PetscInt       *ajtmpold,*ajtmp,nz,row;
240690b6cddSBarry Smith   PetscInt       *diag_offset=b->diag,idx,*ai=a->i,*aj=a->j,*pj;
241329f5518SBarry Smith   MatScalar      *pv,*v,*rtmp,m1,m2,m3,m4,*pc,*w,*x,x1,x2,x3,x4;
2422515f8d2SSatish Balay   MatScalar      p1,p2,p3,p4;
2433f1db9ecSBarry Smith   MatScalar      *ba = b->a,*aa = a->a;
24462bba022SBarry Smith   PetscReal      shift = info->shiftinblocks;
2454eeb42bcSBarry Smith 
2463a40ed3dSBarry Smith   PetscFunctionBegin;
2474eeb42bcSBarry Smith   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
2484eeb42bcSBarry Smith   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
249b0a32e0cSBarry Smith   ierr  = PetscMalloc(4*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
2504eeb42bcSBarry Smith 
2514eeb42bcSBarry Smith   for (i=0; i<n; i++) {
2524078e994SBarry Smith     nz    = bi[i+1] - bi[i];
2534078e994SBarry Smith     ajtmp = bj + bi[i];
2544eeb42bcSBarry Smith     for  (j=0; j<nz; j++) {
2554eeb42bcSBarry Smith       x = rtmp+4*ajtmp[j]; x[0] = x[1] = x[2] = x[3] = 0.0;
2564eeb42bcSBarry Smith     }
2574eeb42bcSBarry Smith     /* load in initial (unfactored row) */
2584eeb42bcSBarry Smith     idx      = r[i];
2594078e994SBarry Smith     nz       = ai[idx+1] - ai[idx];
2604078e994SBarry Smith     ajtmpold = aj + ai[idx];
2614078e994SBarry Smith     v        = aa + 4*ai[idx];
2624eeb42bcSBarry Smith     for (j=0; j<nz; j++) {
2634eeb42bcSBarry Smith       x    = rtmp+4*ic[ajtmpold[j]];
2644eeb42bcSBarry Smith       x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3];
2654eeb42bcSBarry Smith       v    += 4;
2664eeb42bcSBarry Smith     }
2674eeb42bcSBarry Smith     row = *ajtmp++;
2684eeb42bcSBarry Smith     while (row < i) {
2694eeb42bcSBarry Smith       pc = rtmp + 4*row;
2704eeb42bcSBarry Smith       p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3];
27188685aaeSLois Curfman McInnes       if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0) {
2724078e994SBarry Smith         pv = ba + 4*diag_offset[row];
2734078e994SBarry Smith         pj = bj + diag_offset[row] + 1;
2744eeb42bcSBarry Smith         x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
2754eeb42bcSBarry Smith         pc[0] = m1 = p1*x1 + p3*x2;
2764eeb42bcSBarry Smith         pc[1] = m2 = p2*x1 + p4*x2;
2774eeb42bcSBarry Smith         pc[2] = m3 = p1*x3 + p3*x4;
2784eeb42bcSBarry Smith         pc[3] = m4 = p2*x3 + p4*x4;
2794078e994SBarry Smith         nz = bi[row+1] - diag_offset[row] - 1;
2804eeb42bcSBarry Smith         pv += 4;
2814eeb42bcSBarry Smith         for (j=0; j<nz; j++) {
2824eeb42bcSBarry Smith           x1   = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
2834eeb42bcSBarry Smith           x    = rtmp + 4*pj[j];
2844eeb42bcSBarry Smith           x[0] -= m1*x1 + m3*x2;
2854eeb42bcSBarry Smith           x[1] -= m2*x1 + m4*x2;
2864eeb42bcSBarry Smith           x[2] -= m1*x3 + m3*x4;
2874eeb42bcSBarry Smith           x[3] -= m2*x3 + m4*x4;
2884eeb42bcSBarry Smith           pv   += 4;
2894eeb42bcSBarry Smith         }
290dc0b31edSSatish Balay         ierr = PetscLogFlops(16.0*nz+12.0);CHKERRQ(ierr);
2914eeb42bcSBarry Smith       }
2924eeb42bcSBarry Smith       row = *ajtmp++;
2934eeb42bcSBarry Smith     }
2944eeb42bcSBarry Smith     /* finished row so stick it into b->a */
2954078e994SBarry Smith     pv = ba + 4*bi[i];
2964078e994SBarry Smith     pj = bj + bi[i];
2974078e994SBarry Smith     nz = bi[i+1] - bi[i];
2984eeb42bcSBarry Smith     for (j=0; j<nz; j++) {
2994eeb42bcSBarry Smith       x     = rtmp+4*pj[j];
3004eeb42bcSBarry Smith       pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3];
3014eeb42bcSBarry Smith       pv   += 4;
3024eeb42bcSBarry Smith     }
3034eeb42bcSBarry Smith     /* invert diagonal block */
3044078e994SBarry Smith     w = ba + 4*diag_offset[i];
30562bba022SBarry Smith     ierr = Kernel_A_gets_inverse_A_2(w,shift);CHKERRQ(ierr);
3064eeb42bcSBarry Smith   }
3074eeb42bcSBarry Smith 
308606d414cSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
3094eeb42bcSBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
3104eeb42bcSBarry Smith   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
31106e38f1dSHong Zhang   C->ops->solve          = MatSolve_SeqBAIJ_2_inplace;
31206e38f1dSHong Zhang   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_inplace;
3134eeb42bcSBarry Smith   C->assembled = PETSC_TRUE;
314766f9fbaSBarry Smith   ierr = PetscLogFlops(1.333333333333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
3153a40ed3dSBarry Smith   PetscFunctionReturn(0);
3164eeb42bcSBarry Smith }
3174cd38560SBarry Smith /*
3184cd38560SBarry Smith       Version for when blocks are 2 by 2 Using natural ordering
3194cd38560SBarry Smith */
3204a2ae208SSatish Balay #undef __FUNCT__
32106e38f1dSHong Zhang #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_inplace"
32206e38f1dSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_inplace(Mat C,Mat A,const MatFactorInfo *info)
3234cd38560SBarry Smith {
3244cd38560SBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data;
325dfbe8321SBarry Smith   PetscErrorCode ierr;
326690b6cddSBarry Smith   PetscInt       i,j,n = a->mbs,*bi = b->i,*bj = b->j;
327690b6cddSBarry Smith   PetscInt       *ajtmpold,*ajtmp,nz,row;
328690b6cddSBarry Smith   PetscInt       *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj;
329329f5518SBarry Smith   MatScalar      *pv,*v,*rtmp,*pc,*w,*x;
3302515f8d2SSatish Balay   MatScalar      p1,p2,p3,p4,m1,m2,m3,m4,x1,x2,x3,x4;
3314cd38560SBarry Smith   MatScalar      *ba = b->a,*aa = a->a;
33262bba022SBarry Smith   PetscReal      shift = info->shiftinblocks;
3334cd38560SBarry Smith 
3344cd38560SBarry Smith   PetscFunctionBegin;
335b0a32e0cSBarry Smith   ierr = PetscMalloc(4*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
3364cd38560SBarry Smith   for (i=0; i<n; i++) {
3374cd38560SBarry Smith     nz    = bi[i+1] - bi[i];
3384cd38560SBarry Smith     ajtmp = bj + bi[i];
3394cd38560SBarry Smith     for  (j=0; j<nz; j++) {
3404cd38560SBarry Smith       x = rtmp+4*ajtmp[j];
3414cd38560SBarry Smith       x[0]  = x[1]  = x[2]  = x[3]  = 0.0;
3424cd38560SBarry Smith     }
3434cd38560SBarry Smith     /* load in initial (unfactored row) */
3444cd38560SBarry Smith     nz       = ai[i+1] - ai[i];
3454cd38560SBarry Smith     ajtmpold = aj + ai[i];
3464cd38560SBarry Smith     v        = aa + 4*ai[i];
3474cd38560SBarry Smith     for (j=0; j<nz; j++) {
3484cd38560SBarry Smith       x    = rtmp+4*ajtmpold[j];
3494cd38560SBarry Smith       x[0]  = v[0];  x[1]  = v[1];  x[2]  = v[2];  x[3]  = v[3];
3504cd38560SBarry Smith       v    += 4;
3514cd38560SBarry Smith     }
3524cd38560SBarry Smith     row = *ajtmp++;
3534cd38560SBarry Smith     while (row < i) {
3544cd38560SBarry Smith       pc  = rtmp + 4*row;
3554cd38560SBarry Smith       p1  = pc[0];  p2  = pc[1];  p3  = pc[2];  p4  = pc[3];
3564cd38560SBarry Smith       if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0) {
3574cd38560SBarry Smith         pv = ba + 4*diag_offset[row];
3584cd38560SBarry Smith         pj = bj + diag_offset[row] + 1;
3594cd38560SBarry Smith         x1  = pv[0];  x2  = pv[1];  x3  = pv[2];  x4  = pv[3];
3604cd38560SBarry Smith         pc[0] = m1 = p1*x1 + p3*x2;
3614cd38560SBarry Smith         pc[1] = m2 = p2*x1 + p4*x2;
3624cd38560SBarry Smith         pc[2] = m3 = p1*x3 + p3*x4;
3634cd38560SBarry Smith         pc[3] = m4 = p2*x3 + p4*x4;
3644cd38560SBarry Smith         nz = bi[row+1] - diag_offset[row] - 1;
3654cd38560SBarry Smith         pv += 4;
3664cd38560SBarry Smith         for (j=0; j<nz; j++) {
3674cd38560SBarry Smith           x1   = pv[0];  x2  = pv[1];   x3 = pv[2];  x4  = pv[3];
3684cd38560SBarry Smith           x    = rtmp + 4*pj[j];
3694cd38560SBarry Smith           x[0] -= m1*x1 + m3*x2;
3704cd38560SBarry Smith           x[1] -= m2*x1 + m4*x2;
3714cd38560SBarry Smith           x[2] -= m1*x3 + m3*x4;
3724cd38560SBarry Smith           x[3] -= m2*x3 + m4*x4;
3734cd38560SBarry Smith           pv   += 4;
3744cd38560SBarry Smith         }
375dc0b31edSSatish Balay         ierr = PetscLogFlops(16.0*nz+12.0);CHKERRQ(ierr);
3764cd38560SBarry Smith       }
3774cd38560SBarry Smith       row = *ajtmp++;
3784cd38560SBarry Smith     }
3794cd38560SBarry Smith     /* finished row so stick it into b->a */
3804cd38560SBarry Smith     pv = ba + 4*bi[i];
3814cd38560SBarry Smith     pj = bj + bi[i];
3824cd38560SBarry Smith     nz = bi[i+1] - bi[i];
3834cd38560SBarry Smith     for (j=0; j<nz; j++) {
3844cd38560SBarry Smith       x      = rtmp+4*pj[j];
3854cd38560SBarry Smith       pv[0]  = x[0];  pv[1]  = x[1];  pv[2]  = x[2];  pv[3]  = x[3];
3862efa7f71SHong Zhang       /*
3872efa7f71SHong Zhang       printf(" col %d:",pj[j]);
3882efa7f71SHong Zhang       PetscInt j1;
3892efa7f71SHong Zhang       for (j1=0; j1<4; j1++) printf(" %g,",*(pv+j1));
3902efa7f71SHong Zhang       printf("\n");
3912efa7f71SHong Zhang       */
3924cd38560SBarry Smith       pv   += 4;
3934cd38560SBarry Smith     }
3944cd38560SBarry Smith     /* invert diagonal block */
3954cd38560SBarry Smith     w = ba + 4*diag_offset[i];
3962efa7f71SHong Zhang     /*
3972efa7f71SHong Zhang     printf(" \n%d -th: diag: ",i);
3982efa7f71SHong Zhang     for (j=0; j<4; j++){
3992efa7f71SHong Zhang       printf(" %g,",w[j]);
4002efa7f71SHong Zhang     }
4012efa7f71SHong Zhang     printf("\n----------------------------\n");
4022efa7f71SHong Zhang     */
40362bba022SBarry Smith     ierr = Kernel_A_gets_inverse_A_2(w,shift);CHKERRQ(ierr);
4044cd38560SBarry Smith   }
4054cd38560SBarry Smith 
406606d414cSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
40706e38f1dSHong Zhang   C->ops->solve          = MatSolve_SeqBAIJ_2_NaturalOrdering_inplace;
40806e38f1dSHong Zhang   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_NaturalOrdering_inplace;
4094cd38560SBarry Smith   C->assembled = PETSC_TRUE;
410766f9fbaSBarry Smith   ierr = PetscLogFlops(1.333333333333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
4114cd38560SBarry Smith   PetscFunctionReturn(0);
4124cd38560SBarry Smith }
4137fc0212eSBarry Smith 
4147fc0212eSBarry Smith /* ----------------------------------------------------------- */
4157fc0212eSBarry Smith /*
4167fc0212eSBarry Smith      Version for when blocks are 1 by 1.
4177fc0212eSBarry Smith */
4184a2ae208SSatish Balay #undef __FUNCT__
41906e38f1dSHong Zhang #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_1_inplace"
42006e38f1dSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqBAIJ_1_inplace(Mat C,Mat A,const MatFactorInfo *info)
4217fc0212eSBarry Smith {
4227fc0212eSBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data;
4237cdf2dbbSSatish Balay   IS             isrow = b->row,isicol = b->icol;
4246849ba73SBarry Smith   PetscErrorCode ierr;
4255d0c19d7SBarry Smith   const PetscInt *r,*ic;
4265d0c19d7SBarry Smith   PetscInt       i,j,n = a->mbs,*bi = b->i,*bj = b->j;
427690b6cddSBarry Smith   PetscInt       *ajtmpold,*ajtmp,nz,row,*ai = a->i,*aj = a->j;
428690b6cddSBarry Smith   PetscInt       *diag_offset = b->diag,diag,*pj;
429329f5518SBarry Smith   MatScalar      *pv,*v,*rtmp,multiplier,*pc;
4303f1db9ecSBarry Smith   MatScalar      *ba = b->a,*aa = a->a;
431f58c8c74SBarry Smith   PetscTruth     row_identity, col_identity;
4327fc0212eSBarry Smith 
4333a40ed3dSBarry Smith   PetscFunctionBegin;
4347fc0212eSBarry Smith   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
4357fc0212eSBarry Smith   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
436b0a32e0cSBarry Smith   ierr  = PetscMalloc((n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
4377fc0212eSBarry Smith 
4387fc0212eSBarry Smith   for (i=0; i<n; i++) {
4394078e994SBarry Smith     nz    = bi[i+1] - bi[i];
4404078e994SBarry Smith     ajtmp = bj + bi[i];
4417fc0212eSBarry Smith     for  (j=0; j<nz; j++) rtmp[ajtmp[j]] = 0.0;
4427fc0212eSBarry Smith 
4437fc0212eSBarry Smith     /* load in initial (unfactored row) */
4444078e994SBarry Smith     nz       = ai[r[i]+1] - ai[r[i]];
4454078e994SBarry Smith     ajtmpold = aj + ai[r[i]];
4464078e994SBarry Smith     v        = aa + ai[r[i]];
4477fc0212eSBarry Smith     for (j=0; j<nz; j++) rtmp[ic[ajtmpold[j]]] =  v[j];
4487fc0212eSBarry Smith 
4497fc0212eSBarry Smith     row = *ajtmp++;
4507fc0212eSBarry Smith     while (row < i) {
4517fc0212eSBarry Smith       pc = rtmp + row;
4527fc0212eSBarry Smith       if (*pc != 0.0) {
4534078e994SBarry Smith         pv         = ba + diag_offset[row];
4544078e994SBarry Smith         pj         = bj + diag_offset[row] + 1;
4557fc0212eSBarry Smith         multiplier = *pc * *pv++;
4567fc0212eSBarry Smith         *pc        = multiplier;
4574078e994SBarry Smith         nz         = bi[row+1] - diag_offset[row] - 1;
4587fc0212eSBarry Smith         for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j];
459dc0b31edSSatish Balay         ierr = PetscLogFlops(1.0+2.0*nz);CHKERRQ(ierr);
4607fc0212eSBarry Smith       }
4617fc0212eSBarry Smith       row = *ajtmp++;
4627fc0212eSBarry Smith     }
4637fc0212eSBarry Smith     /* finished row so stick it into b->a */
4644078e994SBarry Smith     pv = ba + bi[i];
4654078e994SBarry Smith     pj = bj + bi[i];
4664078e994SBarry Smith     nz = bi[i+1] - bi[i];
4677fc0212eSBarry Smith     for (j=0; j<nz; j++) {pv[j] = rtmp[pj[j]];}
4684078e994SBarry Smith     diag = diag_offset[i] - bi[i];
4697fc0212eSBarry Smith     /* check pivot entry for current row */
470a73cf429SBarry Smith     if (pv[diag] == 0.0) {
4713b4a8b6dSBarry Smith       SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot: row in original ordering %D in permuted ordering %D",r[i],i);
4727fc0212eSBarry Smith     }
4737fc0212eSBarry Smith     pv[diag] = 1.0/pv[diag];
4747fc0212eSBarry Smith   }
4757fc0212eSBarry Smith 
476606d414cSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
4777fc0212eSBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
4787fc0212eSBarry Smith   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
479f58c8c74SBarry Smith   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
480f58c8c74SBarry Smith   ierr = ISIdentity(isicol,&col_identity);CHKERRQ(ierr);
481f58c8c74SBarry Smith   if (row_identity && col_identity) {
48206e38f1dSHong Zhang     C->ops->solve          = MatSolve_SeqBAIJ_1_NaturalOrdering_inplace;
48306e38f1dSHong Zhang     C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1_NaturalOrdering_inplace;
484f58c8c74SBarry Smith   } else {
48506e38f1dSHong Zhang     C->ops->solve          = MatSolve_SeqBAIJ_1_inplace;
48606e38f1dSHong Zhang     C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1_inplace;
487f58c8c74SBarry Smith   }
4887fc0212eSBarry Smith   C->assembled = PETSC_TRUE;
489d0f46423SBarry Smith   ierr = PetscLogFlops(C->cmap->n);CHKERRQ(ierr);
4903a40ed3dSBarry Smith   PetscFunctionReturn(0);
4917fc0212eSBarry Smith }
4927fc0212eSBarry Smith 
493e631078cSBarry Smith EXTERN_C_BEGIN
494b24902e0SBarry Smith #undef __FUNCT__
495b24902e0SBarry Smith #define __FUNCT__ "MatGetFactor_seqbaij_petsc"
496b24902e0SBarry Smith PetscErrorCode MatGetFactor_seqbaij_petsc(Mat A,MatFactorType ftype,Mat *B)
497b24902e0SBarry Smith {
498d0f46423SBarry Smith   PetscInt           n = A->rmap->n;
499b24902e0SBarry Smith   PetscErrorCode     ierr;
500b24902e0SBarry Smith 
501b24902e0SBarry Smith   PetscFunctionBegin;
502b24902e0SBarry Smith   ierr = MatCreate(((PetscObject)A)->comm,B);CHKERRQ(ierr);
503b24902e0SBarry Smith   ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr);
504c8342467SHong Zhang   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU || ftype == MAT_FACTOR_ILUDT) {
505b24902e0SBarry Smith     ierr = MatSetType(*B,MATSEQBAIJ);CHKERRQ(ierr);
5064dd39f65SShri Abhyankar     (*B)->ops->lufactorsymbolic  = MatLUFactorSymbolic_SeqBAIJ;
5074dd39f65SShri Abhyankar     (*B)->ops->ilufactorsymbolic = MatILUFactorSymbolic_SeqBAIJ;
508b24902e0SBarry Smith   } else if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) {
509b24902e0SBarry Smith     ierr = MatSetType(*B,MATSEQSBAIJ);CHKERRQ(ierr);
5105c9eb25fSBarry Smith     ierr = MatSeqSBAIJSetPreallocation(*B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
5115c9eb25fSBarry Smith     (*B)->ops->iccfactorsymbolic      = MatICCFactorSymbolic_SeqBAIJ;
5125c9eb25fSBarry Smith     (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqBAIJ;
513b24902e0SBarry Smith   } else SETERRQ(PETSC_ERR_SUP,"Factor type not supported");
514719d5645SBarry Smith   (*B)->factor = ftype;
515b24902e0SBarry Smith   PetscFunctionReturn(0);
516b24902e0SBarry Smith }
517e631078cSBarry Smith EXTERN_C_END
518273d9f13SBarry Smith 
519db4efbfdSBarry Smith EXTERN_C_BEGIN
520db4efbfdSBarry Smith #undef __FUNCT__
521db4efbfdSBarry Smith #define __FUNCT__ "MatGetFactorAvailable_seqaij_petsc"
522db4efbfdSBarry Smith PetscErrorCode MatGetFactorAvailable_seqbaij_petsc(Mat A,MatFactorType ftype,PetscTruth *flg)
523db4efbfdSBarry Smith {
524db4efbfdSBarry Smith   PetscFunctionBegin;
525db4efbfdSBarry Smith   *flg = PETSC_TRUE;
526db4efbfdSBarry Smith   PetscFunctionReturn(0);
527db4efbfdSBarry Smith }
528db4efbfdSBarry Smith EXTERN_C_END
529db4efbfdSBarry Smith 
5305d517e7eSBarry Smith /* ----------------------------------------------------------- */
5314a2ae208SSatish Balay #undef __FUNCT__
5324a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqBAIJ"
5330481f469SBarry Smith PetscErrorCode MatLUFactor_SeqBAIJ(Mat A,IS row,IS col,const MatFactorInfo *info)
5345d517e7eSBarry Smith {
535dfbe8321SBarry Smith   PetscErrorCode ierr;
5365d517e7eSBarry Smith   Mat            C;
5375d517e7eSBarry Smith 
5383a40ed3dSBarry Smith   PetscFunctionBegin;
53943244d56SBarry Smith   ierr = MatGetFactor(A,MAT_SOLVER_PETSC,MAT_FACTOR_LU,&C);CHKERRQ(ierr);
540719d5645SBarry Smith   ierr = MatLUFactorSymbolic(C,A,row,col,info);CHKERRQ(ierr);
541719d5645SBarry Smith   ierr = MatLUFactorNumeric(C,A,info);CHKERRQ(ierr);
542db4efbfdSBarry Smith   A->ops->solve            = C->ops->solve;
543db4efbfdSBarry Smith   A->ops->solvetranspose   = C->ops->solvetranspose;
544273d9f13SBarry Smith   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
54552e6d16bSBarry Smith   ierr = PetscLogObjectParent(A,((Mat_SeqBAIJ*)(A->data))->icol);CHKERRQ(ierr);
5463a40ed3dSBarry Smith   PetscFunctionReturn(0);
5475d517e7eSBarry Smith }
5484cd38560SBarry Smith 
5497c4f633dSBarry Smith #include "../src/mat/impls/sbaij/seq/sbaij.h"
550c05c3958SHong Zhang #undef __FUNCT__
551c05c3958SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N"
5520481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N(Mat C,Mat A,const MatFactorInfo *info)
553c05c3958SHong Zhang {
554f3086b4bSHong Zhang   PetscErrorCode ierr;
55578910aadSHong Zhang   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data;
55678910aadSHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
55778910aadSHong Zhang   IS             ip=b->row;
5585d0c19d7SBarry Smith   const PetscInt *rip;
5595d0c19d7SBarry Smith   PetscInt       i,j,mbs=a->mbs,bs=A->rmap->bs,*bi=b->i,*bj=b->j,*bcol;
56078910aadSHong Zhang   PetscInt       *ai=a->i,*aj=a->j;
56178910aadSHong Zhang   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
56278910aadSHong Zhang   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
5633cea4cbeSHong Zhang   PetscReal      zeropivot,rs,shiftnz;
564fbf22428SSatish Balay   PetscReal      shiftpd;
5653cea4cbeSHong Zhang   ChShift_Ctx    sctx;
5663cea4cbeSHong Zhang   PetscInt       newshift;
56778910aadSHong Zhang 
568c05c3958SHong Zhang   PetscFunctionBegin;
5696ad2eaddSHong Zhang   if (bs > 1) {
5706ad2eaddSHong Zhang     if (!a->sbaijMat){
571ceb03754SKris Buschelman       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
5726ad2eaddSHong Zhang     }
573719d5645SBarry Smith     ierr = (a->sbaijMat)->ops->choleskyfactornumeric(C,a->sbaijMat,info);CHKERRQ(ierr);
5746ad2eaddSHong Zhang     ierr = MatDestroy(a->sbaijMat);CHKERRQ(ierr);
5756ad2eaddSHong Zhang     a->sbaijMat = PETSC_NULL;
5766ad2eaddSHong Zhang     PetscFunctionReturn(0);
5776ad2eaddSHong Zhang   }
57878910aadSHong Zhang 
57978910aadSHong Zhang   /* initialization */
5803cea4cbeSHong Zhang   shiftnz   = info->shiftnz;
5813cea4cbeSHong Zhang   shiftpd   = info->shiftpd;
5823cea4cbeSHong Zhang   zeropivot = info->zeropivot;
5833cea4cbeSHong Zhang 
5846ad2eaddSHong Zhang   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
585fca92195SBarry Smith   ierr = PetscMalloc3(mbs,MatScalar,&rtmp,mbs,PetscInt,&il,mbs,PetscInt,&jl);CHKERRQ(ierr);
58678910aadSHong Zhang 
58775567043SBarry Smith   sctx.shift_amount = 0.;
5883cea4cbeSHong Zhang   sctx.nshift       = 0;
58978910aadSHong Zhang   do {
5903cea4cbeSHong Zhang     sctx.chshift = PETSC_FALSE;
59178910aadSHong Zhang     for (i=0; i<mbs; i++) {
592e60cf9a0SBarry Smith       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
59378910aadSHong Zhang     }
59478910aadSHong Zhang 
59578910aadSHong Zhang     for (k = 0; k<mbs; k++){
59678910aadSHong Zhang       bval = ba + bi[k];
59778910aadSHong Zhang       /* initialize k-th row by the perm[k]-th row of A */
59878910aadSHong Zhang       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
59978910aadSHong Zhang       for (j = jmin; j < jmax; j++){
60078910aadSHong Zhang         col = rip[aj[j]];
60178910aadSHong Zhang         if (col >= k){ /* only take upper triangular entry */
60278910aadSHong Zhang           rtmp[col] = aa[j];
60378910aadSHong Zhang           *bval++  = 0.0; /* for in-place factorization */
60478910aadSHong Zhang         }
60578910aadSHong Zhang       }
6063cea4cbeSHong Zhang 
6073cea4cbeSHong Zhang       /* shift the diagonal of the matrix */
6083cea4cbeSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
60978910aadSHong Zhang 
61078910aadSHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
61178910aadSHong Zhang       dk = rtmp[k];
61278910aadSHong Zhang       i = jl[k]; /* first row to be added to k_th row  */
61378910aadSHong Zhang 
61478910aadSHong Zhang       while (i < k){
61578910aadSHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
61678910aadSHong Zhang 
61778910aadSHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
61878910aadSHong Zhang         ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
61978910aadSHong Zhang         uikdi = - ba[ili]*ba[bi[i]];  /* diagonal(k) */
62078910aadSHong Zhang         dk += uikdi*ba[ili];
62178910aadSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
62278910aadSHong Zhang 
62378910aadSHong Zhang         /* add multiple of row i to k-th row */
62478910aadSHong Zhang         jmin = ili + 1; jmax = bi[i+1];
62578910aadSHong Zhang         if (jmin < jmax){
62678910aadSHong Zhang           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
62778910aadSHong Zhang           /* update il and jl for row i */
62878910aadSHong Zhang           il[i] = jmin;
62978910aadSHong Zhang           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
63078910aadSHong Zhang         }
63178910aadSHong Zhang         i = nexti;
63278910aadSHong Zhang       }
63378910aadSHong Zhang 
6343cea4cbeSHong Zhang       /* shift the diagonals when zero pivot is detected */
6353cea4cbeSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
6363cea4cbeSHong Zhang       rs   = 0.0;
63778910aadSHong Zhang       jmin = bi[k]+1;
63878910aadSHong Zhang       nz   = bi[k+1] - jmin;
63978910aadSHong Zhang       if (nz){
64078910aadSHong Zhang         bcol = bj + jmin;
64178910aadSHong Zhang         while (nz--){
64289f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
64389f845c8SHong Zhang           bcol++;
64478910aadSHong Zhang         }
64578910aadSHong Zhang       }
6463cea4cbeSHong Zhang 
6473cea4cbeSHong Zhang       sctx.rs = rs;
6483cea4cbeSHong Zhang       sctx.pv = dk;
64945938b79SHong Zhang       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
65045938b79SHong Zhang       if (newshift == 1) break;
65178910aadSHong Zhang 
65278910aadSHong Zhang       /* copy data into U(k,:) */
65378910aadSHong Zhang       ba[bi[k]] = 1.0/dk; /* U(k,k) */
65478910aadSHong Zhang       jmin = bi[k]+1; jmax = bi[k+1];
65578910aadSHong Zhang       if (jmin < jmax) {
65678910aadSHong Zhang         for (j=jmin; j<jmax; j++){
65778910aadSHong Zhang           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
65878910aadSHong Zhang         }
65978910aadSHong Zhang         /* add the k-th row into il and jl */
66078910aadSHong Zhang         il[k] = jmin;
66178910aadSHong Zhang         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
66278910aadSHong Zhang       }
66378910aadSHong Zhang     }
6643cea4cbeSHong Zhang   } while (sctx.chshift);
665fca92195SBarry Smith   ierr = PetscFree3(rtmp,il,jl);CHKERRQ(ierr);
66678910aadSHong Zhang 
66778910aadSHong Zhang   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
66878910aadSHong Zhang   C->assembled    = PETSC_TRUE;
66978910aadSHong Zhang   C->preallocated = PETSC_TRUE;
670d0f46423SBarry Smith   ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr);
6713cea4cbeSHong Zhang   if (sctx.nshift){
672e738e422SBarry Smith     if (shiftpd) {
6731e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
674e738e422SBarry Smith     } else if (shiftnz) {
675e738e422SBarry Smith       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
67678910aadSHong Zhang     }
67778910aadSHong Zhang   }
678c05c3958SHong Zhang   PetscFunctionReturn(0);
679c05c3958SHong Zhang }
6804cd38560SBarry Smith 
681c05c3958SHong Zhang #undef __FUNCT__
682c05c3958SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering"
6830481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info)
684c05c3958SHong Zhang {
68578910aadSHong Zhang   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data;
68678910aadSHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
68778910aadSHong Zhang   PetscErrorCode ierr;
6883cea4cbeSHong Zhang   PetscInt       i,j,am=a->mbs;
68978910aadSHong Zhang   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
6903cea4cbeSHong Zhang   PetscInt       k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz;
69178910aadSHong Zhang   MatScalar      *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval;
6923cea4cbeSHong Zhang   PetscReal      zeropivot,rs,shiftnz;
693fbf22428SSatish Balay   PetscReal      shiftpd;
6943cea4cbeSHong Zhang   ChShift_Ctx    sctx;
6953cea4cbeSHong Zhang   PetscInt       newshift;
69678910aadSHong Zhang 
697c05c3958SHong Zhang   PetscFunctionBegin;
69878910aadSHong Zhang   /* initialization */
6993cea4cbeSHong Zhang   shiftnz   = info->shiftnz;
7003cea4cbeSHong Zhang   shiftpd   = info->shiftpd;
7013cea4cbeSHong Zhang   zeropivot = info->zeropivot;
7023cea4cbeSHong Zhang 
703fca92195SBarry Smith   ierr = PetscMalloc3(am,MatScalar,&rtmp,am,PetscInt,&il,am,PetscInt,&jl);CHKERRQ(ierr);
70478910aadSHong Zhang 
70575567043SBarry Smith   sctx.shift_amount = 0.;
7063cea4cbeSHong Zhang   sctx.nshift       = 0;
70778910aadSHong Zhang   do {
7083cea4cbeSHong Zhang     sctx.chshift = PETSC_FALSE;
70978910aadSHong Zhang     for (i=0; i<am; i++) {
710e60cf9a0SBarry Smith       rtmp[i] = 0.0; jl[i] = am; il[0] = 0;
71178910aadSHong Zhang     }
71278910aadSHong Zhang 
71378910aadSHong Zhang     for (k = 0; k<am; k++){
71478910aadSHong Zhang     /* initialize k-th row with elements nonzero in row perm(k) of A */
71578910aadSHong Zhang       nz   = ai[k+1] - ai[k];
71678910aadSHong Zhang       acol = aj + ai[k];
71778910aadSHong Zhang       aval = aa + ai[k];
71878910aadSHong Zhang       bval = ba + bi[k];
71978910aadSHong Zhang       while (nz -- ){
72078910aadSHong Zhang         if (*acol < k) { /* skip lower triangular entries */
72178910aadSHong Zhang           acol++; aval++;
72278910aadSHong Zhang         } else {
72378910aadSHong Zhang           rtmp[*acol++] = *aval++;
72478910aadSHong Zhang           *bval++       = 0.0; /* for in-place factorization */
72578910aadSHong Zhang         }
72678910aadSHong Zhang       }
7273cea4cbeSHong Zhang 
7283cea4cbeSHong Zhang       /* shift the diagonal of the matrix */
7293cea4cbeSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
73078910aadSHong Zhang 
73178910aadSHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
73278910aadSHong Zhang       dk = rtmp[k];
73378910aadSHong Zhang       i  = jl[k]; /* first row to be added to k_th row  */
73478910aadSHong Zhang 
73578910aadSHong Zhang       while (i < k){
73678910aadSHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
73778910aadSHong Zhang         /* compute multiplier, update D(k) and U(i,k) */
73878910aadSHong Zhang         ili   = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
73978910aadSHong Zhang         uikdi = - ba[ili]*ba[bi[i]];
74078910aadSHong Zhang         dk   += uikdi*ba[ili];
74178910aadSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
74278910aadSHong Zhang 
74378910aadSHong Zhang         /* add multiple of row i to k-th row ... */
74478910aadSHong Zhang         jmin = ili + 1;
74578910aadSHong Zhang         nz   = bi[i+1] - jmin;
74678910aadSHong Zhang         if (nz > 0){
74778910aadSHong Zhang           bcol = bj + jmin;
74878910aadSHong Zhang           bval = ba + jmin;
74978910aadSHong Zhang           while (nz --) rtmp[*bcol++] += uikdi*(*bval++);
75078910aadSHong Zhang           /* update il and jl for i-th row */
75178910aadSHong Zhang           il[i] = jmin;
75278910aadSHong Zhang           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
75378910aadSHong Zhang         }
75478910aadSHong Zhang         i = nexti;
75578910aadSHong Zhang       }
75678910aadSHong Zhang 
7573cea4cbeSHong Zhang       /* shift the diagonals when zero pivot is detected */
7583cea4cbeSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
7593cea4cbeSHong Zhang       rs   = 0.0;
76078910aadSHong Zhang       jmin = bi[k]+1;
76178910aadSHong Zhang       nz   = bi[k+1] - jmin;
76278910aadSHong Zhang       if (nz){
76378910aadSHong Zhang         bcol = bj + jmin;
76478910aadSHong Zhang         while (nz--){
76589f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
76689f845c8SHong Zhang           bcol++;
76778910aadSHong Zhang         }
76878910aadSHong Zhang       }
7693cea4cbeSHong Zhang 
7703cea4cbeSHong Zhang       sctx.rs = rs;
7713cea4cbeSHong Zhang       sctx.pv = dk;
77245938b79SHong Zhang       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
77345938b79SHong Zhang       if (newshift == 1) break;    /* sctx.shift_amount is updated */
77478910aadSHong Zhang 
77578910aadSHong Zhang       /* copy data into U(k,:) */
77678910aadSHong Zhang       ba[bi[k]] = 1.0/dk;
77778910aadSHong Zhang       jmin      = bi[k]+1;
77878910aadSHong Zhang       nz        = bi[k+1] - jmin;
77978910aadSHong Zhang       if (nz){
78078910aadSHong Zhang         bcol = bj + jmin;
78178910aadSHong Zhang         bval = ba + jmin;
78278910aadSHong Zhang         while (nz--){
78378910aadSHong Zhang           *bval++       = rtmp[*bcol];
78478910aadSHong Zhang           rtmp[*bcol++] = 0.0;
78578910aadSHong Zhang         }
78678910aadSHong Zhang         /* add k-th row into il and jl */
78778910aadSHong Zhang         il[k] = jmin;
78878910aadSHong Zhang         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
78978910aadSHong Zhang       }
79078910aadSHong Zhang     }
7913cea4cbeSHong Zhang   } while (sctx.chshift);
792fca92195SBarry Smith   ierr = PetscFree3(rtmp,il,jl);CHKERRQ(ierr);
79378910aadSHong Zhang 
7940a3c351aSHong Zhang   C->ops->solve                 = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
7950a3c351aSHong Zhang   C->ops->solvetranspose        = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
79678910aadSHong Zhang   C->assembled    = PETSC_TRUE;
79778910aadSHong Zhang   C->preallocated = PETSC_TRUE;
798d0f46423SBarry Smith   ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr);
7993cea4cbeSHong Zhang     if (sctx.nshift){
8003cea4cbeSHong Zhang     if (shiftnz) {
8011e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
8023cea4cbeSHong Zhang     } else if (shiftpd) {
8031e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
80478910aadSHong Zhang     }
80578910aadSHong Zhang   }
806c05c3958SHong Zhang   PetscFunctionReturn(0);
807c05c3958SHong Zhang }
808c05c3958SHong Zhang 
809c05c3958SHong Zhang #include "petscbt.h"
8107c4f633dSBarry Smith #include "../src/mat/utils/freespace.h"
811c05c3958SHong Zhang #undef __FUNCT__
812c05c3958SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqBAIJ"
8130481f469SBarry Smith PetscErrorCode MatICCFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
814c05c3958SHong Zhang {
81578910aadSHong Zhang   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data;
81678910aadSHong Zhang   Mat_SeqSBAIJ       *b;
81778910aadSHong Zhang   Mat                B;
81878910aadSHong Zhang   PetscErrorCode     ierr;
81978910aadSHong Zhang   PetscTruth         perm_identity;
8205d0c19d7SBarry Smith   PetscInt           reallocs=0,i,*ai=a->i,*aj=a->j,am=a->mbs,bs=A->rmap->bs,*ui;
8215d0c19d7SBarry Smith   const PetscInt     *rip;
82278910aadSHong Zhang   PetscInt           jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow;
823cfdb8b8aSHong Zhang   PetscInt           nlnk,*lnk,*lnk_lvl=PETSC_NULL,ncols,ncols_upper,*cols,*cols_lvl,*uj,**uj_ptr,**uj_lvl_ptr;
82478910aadSHong Zhang   PetscReal          fill=info->fill,levels=info->levels;
825a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
826a1a86e44SBarry Smith   PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
82778910aadSHong Zhang   PetscBT            lnkbt;
82878910aadSHong Zhang 
829c05c3958SHong Zhang   PetscFunctionBegin;
8306ad2eaddSHong Zhang   if (bs > 1){
8316ad2eaddSHong Zhang     if (!a->sbaijMat){
832ceb03754SKris Buschelman       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
8336ad2eaddSHong Zhang     }
834719d5645SBarry Smith     (fact)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqSBAIJ;  /* undue the change made in MatGetFactor_seqbaij_petsc */
835719d5645SBarry Smith     ierr = MatICCFactorSymbolic(fact,a->sbaijMat,perm,info);CHKERRQ(ierr);
8366ad2eaddSHong Zhang     PetscFunctionReturn(0);
8376ad2eaddSHong Zhang   }
8386ad2eaddSHong Zhang 
83978910aadSHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
84078910aadSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
84178910aadSHong Zhang 
84278910aadSHong Zhang   /* special case that simply copies fill pattern */
84378910aadSHong Zhang   if (!levels && perm_identity) {
84478910aadSHong Zhang     ierr = MatMarkDiagonal_SeqBAIJ(A);CHKERRQ(ierr);
84578910aadSHong Zhang     ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
84678910aadSHong Zhang     for (i=0; i<am; i++) {
84778910aadSHong Zhang       ui[i] = ai[i+1] - a->diag[i]; /* ui: rowlengths - changes when !perm_identity */
84878910aadSHong Zhang     }
849719d5645SBarry Smith     B = fact;
85078910aadSHong Zhang     ierr = MatSeqSBAIJSetPreallocation(B,1,0,ui);CHKERRQ(ierr);
85178910aadSHong Zhang 
852b24902e0SBarry Smith 
85378910aadSHong Zhang     b  = (Mat_SeqSBAIJ*)B->data;
85478910aadSHong Zhang     uj = b->j;
85578910aadSHong Zhang     for (i=0; i<am; i++) {
85678910aadSHong Zhang       aj = a->j + a->diag[i];
85778910aadSHong Zhang       for (j=0; j<ui[i]; j++){
85878910aadSHong Zhang         *uj++ = *aj++;
85978910aadSHong Zhang       }
86078910aadSHong Zhang       b->ilen[i] = ui[i];
86178910aadSHong Zhang     }
86278910aadSHong Zhang     ierr = PetscFree(ui);CHKERRQ(ierr);
863719d5645SBarry Smith     B->factor = MAT_FACTOR_NONE;
86478910aadSHong Zhang     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
86578910aadSHong Zhang     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
866719d5645SBarry Smith     B->factor = MAT_FACTOR_ICC;
86778910aadSHong Zhang 
86878910aadSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
86978910aadSHong Zhang     PetscFunctionReturn(0);
87078910aadSHong Zhang   }
87178910aadSHong Zhang 
87278910aadSHong Zhang   /* initialization */
87378910aadSHong Zhang   ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
874e60cf9a0SBarry Smith   ui[0] = 0;
87578910aadSHong Zhang   ierr  = PetscMalloc((2*am+1)*sizeof(PetscInt),&cols_lvl);CHKERRQ(ierr);
87678910aadSHong Zhang 
87778910aadSHong Zhang   /* jl: linked list for storing indices of the pivot rows
87878910aadSHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
879fca92195SBarry Smith   ierr = PetscMalloc4(am,PetscInt*,&uj_ptr,am,PetscInt*,&uj_lvl_ptr,am,PetscInt,&il,am,PetscInt,&jl);CHKERRQ(ierr);
88078910aadSHong Zhang   for (i=0; i<am; i++){
881e60cf9a0SBarry Smith     jl[i] = am; il[i] = 0;
88278910aadSHong Zhang   }
88378910aadSHong Zhang 
88478910aadSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
88578910aadSHong Zhang   nlnk = am + 1;
88678910aadSHong Zhang   ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
88778910aadSHong Zhang 
88878910aadSHong Zhang   /* initial FreeSpace size is fill*(ai[am]+1) */
889a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
89078910aadSHong Zhang   current_space = free_space;
891a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr);
89278910aadSHong Zhang   current_space_lvl = free_space_lvl;
89378910aadSHong Zhang 
89478910aadSHong Zhang   for (k=0; k<am; k++){  /* for each active row k */
89578910aadSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
89678910aadSHong Zhang     nzk   = 0;
89778910aadSHong Zhang     ncols = ai[rip[k]+1] - ai[rip[k]];
89878910aadSHong Zhang     ncols_upper = 0;
89978910aadSHong Zhang     cols        = cols_lvl + am;
90078910aadSHong Zhang     for (j=0; j<ncols; j++){
90178910aadSHong Zhang       i = rip[*(aj + ai[rip[k]] + j)];
90278910aadSHong Zhang       if (i >= k){ /* only take upper triangular entry */
90378910aadSHong Zhang         cols[ncols_upper] = i;
90478910aadSHong Zhang         cols_lvl[ncols_upper] = -1;  /* initialize level for nonzero entries */
90578910aadSHong Zhang         ncols_upper++;
90678910aadSHong Zhang       }
90778910aadSHong Zhang     }
90878910aadSHong Zhang     ierr = PetscIncompleteLLAdd(ncols_upper,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
90978910aadSHong Zhang     nzk += nlnk;
91078910aadSHong Zhang 
91178910aadSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
91278910aadSHong Zhang     prow = jl[k]; /* 1st pivot row */
91378910aadSHong Zhang 
91478910aadSHong Zhang     while (prow < k){
91578910aadSHong Zhang       nextprow = jl[prow];
91678910aadSHong Zhang 
91778910aadSHong Zhang       /* merge prow into k-th row */
91878910aadSHong Zhang       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
91978910aadSHong Zhang       jmax = ui[prow+1];
92078910aadSHong Zhang       ncols = jmax-jmin;
92178910aadSHong Zhang       i     = jmin - ui[prow];
92278910aadSHong Zhang       cols = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
92378910aadSHong Zhang       for (j=0; j<ncols; j++) cols_lvl[j] = *(uj_lvl_ptr[prow] + i + j);
9245a8e39fbSHong Zhang       ierr = PetscIncompleteLLAddSorted(ncols,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
92578910aadSHong Zhang       nzk += nlnk;
92678910aadSHong Zhang 
92778910aadSHong Zhang       /* update il and jl for prow */
92878910aadSHong Zhang       if (jmin < jmax){
92978910aadSHong Zhang         il[prow] = jmin;
93078910aadSHong Zhang         j = *cols; jl[prow] = jl[j]; jl[j] = prow;
93178910aadSHong Zhang       }
93278910aadSHong Zhang       prow = nextprow;
93378910aadSHong Zhang     }
93478910aadSHong Zhang 
93578910aadSHong Zhang     /* if free space is not available, make more free space */
93678910aadSHong Zhang     if (current_space->local_remaining<nzk) {
93778910aadSHong Zhang       i = am - k + 1; /* num of unfactored rows */
93878910aadSHong Zhang       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
939a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
940a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(i,&current_space_lvl);CHKERRQ(ierr);
94178910aadSHong Zhang       reallocs++;
94278910aadSHong Zhang     }
94378910aadSHong Zhang 
94478910aadSHong Zhang     /* copy data into free_space and free_space_lvl, then initialize lnk */
94578910aadSHong Zhang     ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
94678910aadSHong Zhang 
94778910aadSHong Zhang     /* add the k-th row into il and jl */
94878910aadSHong Zhang     if (nzk-1 > 0){
94978910aadSHong Zhang       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
95078910aadSHong Zhang       jl[k] = jl[i]; jl[i] = k;
95178910aadSHong Zhang       il[k] = ui[k] + 1;
95278910aadSHong Zhang     }
95378910aadSHong Zhang     uj_ptr[k]     = current_space->array;
95478910aadSHong Zhang     uj_lvl_ptr[k] = current_space_lvl->array;
95578910aadSHong Zhang 
95678910aadSHong Zhang     current_space->array           += nzk;
95778910aadSHong Zhang     current_space->local_used      += nzk;
95878910aadSHong Zhang     current_space->local_remaining -= nzk;
95978910aadSHong Zhang 
96078910aadSHong Zhang     current_space_lvl->array           += nzk;
96178910aadSHong Zhang     current_space_lvl->local_used      += nzk;
96278910aadSHong Zhang     current_space_lvl->local_remaining -= nzk;
96378910aadSHong Zhang 
96478910aadSHong Zhang     ui[k+1] = ui[k] + nzk;
96578910aadSHong Zhang   }
96678910aadSHong Zhang 
9676cf91177SBarry Smith #if defined(PETSC_USE_INFO)
96878910aadSHong Zhang   if (ai[am] != 0) {
96978910aadSHong Zhang     PetscReal af = ((PetscReal)(2*ui[am]-am))/((PetscReal)ai[am]);
970ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
971ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
972ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
97378910aadSHong Zhang   } else {
974ae15b995SBarry Smith     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
97578910aadSHong Zhang   }
97663ba0a88SBarry Smith #endif
97778910aadSHong Zhang 
97878910aadSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
979fca92195SBarry Smith   ierr = PetscFree4(uj_ptr,uj_lvl_ptr,il,jl);CHKERRQ(ierr);
98078910aadSHong Zhang   ierr = PetscFree(cols_lvl);CHKERRQ(ierr);
98178910aadSHong Zhang 
98278910aadSHong Zhang   /* destroy list of free space and other temporary array(s) */
98378910aadSHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
984a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
98578910aadSHong Zhang   ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
986a1a86e44SBarry Smith   ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
98778910aadSHong Zhang 
98878910aadSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
989719d5645SBarry Smith   B = fact;
990ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation(B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
99178910aadSHong Zhang 
99278910aadSHong Zhang   b = (Mat_SeqSBAIJ*)B->data;
99378910aadSHong Zhang   b->singlemalloc = PETSC_FALSE;
994e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
995e6b907acSBarry Smith   b->free_ij       = PETSC_TRUE;
99678910aadSHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
99778910aadSHong Zhang   b->j    = uj;
99878910aadSHong Zhang   b->i    = ui;
99978910aadSHong Zhang   b->diag = 0;
100078910aadSHong Zhang   b->ilen = 0;
100178910aadSHong Zhang   b->imax = 0;
100278910aadSHong Zhang   b->row  = perm;
100378910aadSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
100478910aadSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
100578910aadSHong Zhang   b->icol = perm;
100678910aadSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
100778910aadSHong Zhang   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
100852e6d16bSBarry Smith   ierr    = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
100978910aadSHong Zhang   b->maxnz = b->nz = ui[am];
101078910aadSHong Zhang 
101178910aadSHong Zhang   B->info.factor_mallocs    = reallocs;
101278910aadSHong Zhang   B->info.fill_ratio_given  = fill;
101375567043SBarry Smith   if (ai[am] != 0.) {
101478910aadSHong Zhang     B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
101578910aadSHong Zhang   } else {
101678910aadSHong Zhang     B->info.fill_ratio_needed = 0.0;
101778910aadSHong Zhang   }
101878910aadSHong Zhang   if (perm_identity){
10190a3c351aSHong Zhang     B->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
10200a3c351aSHong Zhang     B->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
102178910aadSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
102278910aadSHong Zhang   } else {
1023719d5645SBarry Smith     (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
102478910aadSHong Zhang   }
1025c05c3958SHong Zhang   PetscFunctionReturn(0);
1026c05c3958SHong Zhang }
1027c05c3958SHong Zhang 
1028c05c3958SHong Zhang #undef __FUNCT__
1029c05c3958SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqBAIJ"
10300481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
1031c05c3958SHong Zhang {
103278910aadSHong Zhang   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data;
103378910aadSHong Zhang   Mat_SeqSBAIJ       *b;
103478910aadSHong Zhang   Mat                B;
103578910aadSHong Zhang   PetscErrorCode     ierr;
103678910aadSHong Zhang   PetscTruth         perm_identity;
103778910aadSHong Zhang   PetscReal          fill = info->fill;
10385d0c19d7SBarry Smith   const PetscInt     *rip;
10395d0c19d7SBarry Smith   PetscInt           i,mbs=a->mbs,bs=A->rmap->bs,*ai=a->i,*aj=a->j,reallocs=0,prow;
104078910aadSHong Zhang   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
104178910aadSHong Zhang   PetscInt           nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr;
1042a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
104378910aadSHong Zhang   PetscBT            lnkbt;
104478910aadSHong Zhang 
1045c05c3958SHong Zhang   PetscFunctionBegin;
10466ad2eaddSHong Zhang   if (bs > 1) { /* convert to seqsbaij */
10476ad2eaddSHong Zhang     if (!a->sbaijMat){
1048ceb03754SKris Buschelman       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
10496ad2eaddSHong Zhang     }
1050719d5645SBarry Smith     (fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqSBAIJ; /* undue the change made in MatGetFactor_seqbaij_petsc */
1051719d5645SBarry Smith     ierr = MatCholeskyFactorSymbolic(fact,a->sbaijMat,perm,info);CHKERRQ(ierr);
10526ad2eaddSHong Zhang     PetscFunctionReturn(0);
10536ad2eaddSHong Zhang   }
10546ad2eaddSHong Zhang 
105578910aadSHong Zhang   /* check whether perm is the identity mapping */
105678910aadSHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
1057c84f5b01SHong Zhang   if (!perm_identity) SETERRQ(PETSC_ERR_SUP,"Matrix reordering is not supported");
105878910aadSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
105978910aadSHong Zhang 
106078910aadSHong Zhang   /* initialization */
106178910aadSHong Zhang   ierr  = PetscMalloc((mbs+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
1062e60cf9a0SBarry Smith   ui[0] = 0;
106378910aadSHong Zhang 
106478910aadSHong Zhang   /* jl: linked list for storing indices of the pivot rows
106578910aadSHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
1066fca92195SBarry Smith   ierr = PetscMalloc4(mbs,PetscInt*,&ui_ptr,mbs,PetscInt,&il,mbs,PetscInt,&jl,mbs,PetscInt,&cols);CHKERRQ(ierr);
106778910aadSHong Zhang   for (i=0; i<mbs; i++){
1068e60cf9a0SBarry Smith     jl[i] = mbs; il[i] = 0;
106978910aadSHong Zhang   }
107078910aadSHong Zhang 
107178910aadSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
107278910aadSHong Zhang   nlnk = mbs + 1;
107378910aadSHong Zhang   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
107478910aadSHong Zhang 
107578910aadSHong Zhang   /* initial FreeSpace size is fill*(ai[mbs]+1) */
1076a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[mbs]+1)),&free_space);CHKERRQ(ierr);
107778910aadSHong Zhang   current_space = free_space;
107878910aadSHong Zhang 
107978910aadSHong Zhang   for (k=0; k<mbs; k++){  /* for each active row k */
108078910aadSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
108178910aadSHong Zhang     nzk   = 0;
108278910aadSHong Zhang     ncols = ai[rip[k]+1] - ai[rip[k]];
1083e60cf9a0SBarry Smith     ncols_upper = 0;
108478910aadSHong Zhang     for (j=0; j<ncols; j++){
108578910aadSHong Zhang       i = rip[*(aj + ai[rip[k]] + j)];
108678910aadSHong Zhang       if (i >= k){ /* only take upper triangular entry */
108778910aadSHong Zhang         cols[ncols_upper] = i;
108878910aadSHong Zhang         ncols_upper++;
108978910aadSHong Zhang       }
109078910aadSHong Zhang     }
109178910aadSHong Zhang     ierr = PetscLLAdd(ncols_upper,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
109278910aadSHong Zhang     nzk += nlnk;
109378910aadSHong Zhang 
109478910aadSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
109578910aadSHong Zhang     prow = jl[k]; /* 1st pivot row */
109678910aadSHong Zhang 
109778910aadSHong Zhang     while (prow < k){
109878910aadSHong Zhang       nextprow = jl[prow];
109978910aadSHong Zhang       /* merge prow into k-th row */
110078910aadSHong Zhang       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
110178910aadSHong Zhang       jmax = ui[prow+1];
110278910aadSHong Zhang       ncols = jmax-jmin;
110378910aadSHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
11045a8e39fbSHong Zhang       ierr = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
110578910aadSHong Zhang       nzk += nlnk;
110678910aadSHong Zhang 
110778910aadSHong Zhang       /* update il and jl for prow */
110878910aadSHong Zhang       if (jmin < jmax){
110978910aadSHong Zhang         il[prow] = jmin;
111078910aadSHong Zhang         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
111178910aadSHong Zhang       }
111278910aadSHong Zhang       prow = nextprow;
111378910aadSHong Zhang     }
111478910aadSHong Zhang 
111578910aadSHong Zhang     /* if free space is not available, make more free space */
111678910aadSHong Zhang     if (current_space->local_remaining<nzk) {
111778910aadSHong Zhang       i = mbs - k + 1; /* num of unfactored rows */
111878910aadSHong Zhang       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
1119a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
112078910aadSHong Zhang       reallocs++;
112178910aadSHong Zhang     }
112278910aadSHong Zhang 
112378910aadSHong Zhang     /* copy data into free space, then initialize lnk */
112478910aadSHong Zhang     ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
112578910aadSHong Zhang 
112678910aadSHong Zhang     /* add the k-th row into il and jl */
112778910aadSHong Zhang     if (nzk-1 > 0){
112878910aadSHong Zhang       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
112978910aadSHong Zhang       jl[k] = jl[i]; jl[i] = k;
113078910aadSHong Zhang       il[k] = ui[k] + 1;
113178910aadSHong Zhang     }
113278910aadSHong Zhang     ui_ptr[k] = current_space->array;
113378910aadSHong Zhang     current_space->array           += nzk;
113478910aadSHong Zhang     current_space->local_used      += nzk;
113578910aadSHong Zhang     current_space->local_remaining -= nzk;
113678910aadSHong Zhang 
113778910aadSHong Zhang     ui[k+1] = ui[k] + nzk;
113878910aadSHong Zhang   }
113978910aadSHong Zhang 
11406cf91177SBarry Smith #if defined(PETSC_USE_INFO)
114175567043SBarry Smith   if (ai[mbs] != 0.) {
114278910aadSHong Zhang     PetscReal af = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
1143ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
1144ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
1145ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
114678910aadSHong Zhang   } else {
1147ae15b995SBarry Smith     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
114878910aadSHong Zhang   }
114963ba0a88SBarry Smith #endif
115078910aadSHong Zhang 
115178910aadSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
1152fca92195SBarry Smith   ierr = PetscFree4(ui_ptr,il,jl,cols);CHKERRQ(ierr);
115378910aadSHong Zhang 
115478910aadSHong Zhang   /* destroy list of free space and other temporary array(s) */
115578910aadSHong Zhang   ierr = PetscMalloc((ui[mbs]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
1156a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
115778910aadSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
115878910aadSHong Zhang 
115978910aadSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
1160719d5645SBarry Smith   B    = fact;
1161ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
116278910aadSHong Zhang 
116378910aadSHong Zhang   b = (Mat_SeqSBAIJ*)B->data;
116478910aadSHong Zhang   b->singlemalloc = PETSC_FALSE;
1165e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
1166e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
116778910aadSHong Zhang   ierr = PetscMalloc((ui[mbs]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
116878910aadSHong Zhang   b->j    = uj;
116978910aadSHong Zhang   b->i    = ui;
117078910aadSHong Zhang   b->diag = 0;
117178910aadSHong Zhang   b->ilen = 0;
117278910aadSHong Zhang   b->imax = 0;
117378910aadSHong Zhang   b->row  = perm;
117478910aadSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
117578910aadSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
117678910aadSHong Zhang   b->icol = perm;
117778910aadSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
117878910aadSHong Zhang   ierr    = PetscMalloc((mbs+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
117952e6d16bSBarry Smith   ierr    = PetscLogObjectMemory(B,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
118078910aadSHong Zhang   b->maxnz = b->nz = ui[mbs];
118178910aadSHong Zhang 
118278910aadSHong Zhang   B->info.factor_mallocs    = reallocs;
118378910aadSHong Zhang   B->info.fill_ratio_given  = fill;
118475567043SBarry Smith   if (ai[mbs] != 0.) {
118578910aadSHong Zhang     B->info.fill_ratio_needed = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
118678910aadSHong Zhang   } else {
118778910aadSHong Zhang     B->info.fill_ratio_needed = 0.0;
118878910aadSHong Zhang   }
118978910aadSHong Zhang   if (perm_identity){
11906ad2eaddSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
119178910aadSHong Zhang   } else {
11926ad2eaddSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
119378910aadSHong Zhang   }
1194c05c3958SHong Zhang   PetscFunctionReturn(0);
1195c05c3958SHong Zhang }
1196c8342467SHong Zhang 
11972efa7f71SHong Zhang #undef __FUNCT__
11984dd39f65SShri Abhyankar #define __FUNCT__ "MatSolve_SeqBAIJ_N_NaturalOrdering"
11994dd39f65SShri Abhyankar PetscErrorCode MatSolve_SeqBAIJ_N_NaturalOrdering(Mat A,Vec bb,Vec xx)
12001a83e813SShri Abhyankar {
12011a83e813SShri Abhyankar   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ *)A->data;
12021a83e813SShri Abhyankar   PetscErrorCode ierr;
12031a83e813SShri Abhyankar   const PetscInt *ai=a->i,*aj=a->j,*adiag=a->diag,*vi;
12041a83e813SShri Abhyankar   PetscInt       i,k,n=a->mbs;
12051a83e813SShri Abhyankar   PetscInt       nz,bs=A->rmap->bs,bs2=a->bs2;
12061a83e813SShri Abhyankar   MatScalar      *aa=a->a,*v;
12071a83e813SShri Abhyankar   PetscScalar    *x,*b,*s,*t,*ls;
12081a83e813SShri Abhyankar 
12091a83e813SShri Abhyankar   PetscFunctionBegin;
12101a83e813SShri Abhyankar   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
12111a83e813SShri Abhyankar   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
12121a83e813SShri Abhyankar   t  = a->solve_work;
12131a83e813SShri Abhyankar 
12141a83e813SShri Abhyankar   /* forward solve the lower triangular */
12151a83e813SShri Abhyankar   ierr = PetscMemcpy(t,b,bs*sizeof(PetscScalar));CHKERRQ(ierr); /* copy 1st block of b to t */
12161a83e813SShri Abhyankar 
12171a83e813SShri Abhyankar   for (i=1; i<n; i++) {
12181a83e813SShri Abhyankar     v   = aa + bs2*ai[i];
12191a83e813SShri Abhyankar     vi  = aj + ai[i];
12201a83e813SShri Abhyankar     nz = ai[i+1] - ai[i];
12211a83e813SShri Abhyankar     s = t + bs*i;
12221a83e813SShri Abhyankar     ierr = PetscMemcpy(s,b+bs*i,bs*sizeof(PetscScalar));CHKERRQ(ierr); /* copy i_th block of b to t */
12231a83e813SShri Abhyankar     for(k=0;k<nz;k++){
12241a83e813SShri Abhyankar       Kernel_v_gets_v_minus_A_times_w(bs,s,v,t+bs*vi[k]);
12251a83e813SShri Abhyankar       v += bs2;
12261a83e813SShri Abhyankar     }
12271a83e813SShri Abhyankar   }
12281a83e813SShri Abhyankar 
12291a83e813SShri Abhyankar   /* backward solve the upper triangular */
12301a83e813SShri Abhyankar   ls = a->solve_work + A->cmap->n;
12311a83e813SShri Abhyankar   for (i=n-1; i>=0; i--){
12321a83e813SShri Abhyankar     v  = aa + bs2*(adiag[i+1]+1);
12331a83e813SShri Abhyankar     vi = aj + adiag[i+1]+1;
12341a83e813SShri Abhyankar     nz = adiag[i] - adiag[i+1]-1;
12351a83e813SShri Abhyankar     ierr = PetscMemcpy(ls,t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr);
12361a83e813SShri Abhyankar     for(k=0;k<nz;k++){
12371a83e813SShri Abhyankar       Kernel_v_gets_v_minus_A_times_w(bs,ls,v,t+bs*vi[k]);
12381a83e813SShri Abhyankar       v += bs2;
12391a83e813SShri Abhyankar     }
12401a83e813SShri Abhyankar     Kernel_w_gets_A_times_v(bs,ls,aa+bs2*adiag[i],t+i*bs); /* *inv(diagonal[i]) */
12411a83e813SShri Abhyankar     ierr = PetscMemcpy(x+i*bs,t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr);
12421a83e813SShri Abhyankar   }
12431a83e813SShri Abhyankar 
12441a83e813SShri Abhyankar   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
12451a83e813SShri Abhyankar   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
12461a83e813SShri Abhyankar   ierr = PetscLogFlops(2.0*(a->bs2)*(a->nz) - A->rmap->bs*A->cmap->n);CHKERRQ(ierr);
12471a83e813SShri Abhyankar   PetscFunctionReturn(0);
12481a83e813SShri Abhyankar }
12491a83e813SShri Abhyankar 
12502efa7f71SHong Zhang #undef __FUNCT__
12514dd39f65SShri Abhyankar #define __FUNCT__ "MatSolve_SeqBAIJ_N"
12524dd39f65SShri Abhyankar PetscErrorCode MatSolve_SeqBAIJ_N(Mat A,Vec bb,Vec xx)
125335aa4fcfSShri Abhyankar {
125435aa4fcfSShri Abhyankar   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ *)A->data;
125535aa4fcfSShri Abhyankar   IS             iscol=a->col,isrow=a->row;
125635aa4fcfSShri Abhyankar   PetscErrorCode ierr;
125735aa4fcfSShri Abhyankar   const PetscInt *r,*c,*rout,*cout,*ai=a->i,*aj=a->j,*adiag=a->diag,*vi;
125835aa4fcfSShri Abhyankar   PetscInt       i,m,n=a->mbs;
125935aa4fcfSShri Abhyankar   PetscInt       nz,bs=A->rmap->bs,bs2=a->bs2;
126035aa4fcfSShri Abhyankar   MatScalar      *aa=a->a,*v;
126135aa4fcfSShri Abhyankar   PetscScalar    *x,*b,*s,*t,*ls;
126235aa4fcfSShri Abhyankar 
126335aa4fcfSShri Abhyankar   PetscFunctionBegin;
126435aa4fcfSShri Abhyankar   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
126535aa4fcfSShri Abhyankar   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
126635aa4fcfSShri Abhyankar   t  = a->solve_work;
126735aa4fcfSShri Abhyankar 
126835aa4fcfSShri Abhyankar   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
126935aa4fcfSShri Abhyankar   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
127035aa4fcfSShri Abhyankar 
127135aa4fcfSShri Abhyankar   /* forward solve the lower triangular */
127235aa4fcfSShri Abhyankar   ierr = PetscMemcpy(t,b+bs*r[0],bs*sizeof(PetscScalar));CHKERRQ(ierr);
127335aa4fcfSShri Abhyankar   for (i=1; i<n; i++) {
127435aa4fcfSShri Abhyankar     v   = aa + bs2*ai[i];
127535aa4fcfSShri Abhyankar     vi  = aj + ai[i];
127635aa4fcfSShri Abhyankar     nz = ai[i+1] - ai[i];
127735aa4fcfSShri Abhyankar     s = t + bs*i;
127835aa4fcfSShri Abhyankar     ierr = PetscMemcpy(s,b+bs*r[i],bs*sizeof(PetscScalar));CHKERRQ(ierr);
127935aa4fcfSShri Abhyankar     for(m=0;m<nz;m++){
128035aa4fcfSShri Abhyankar       Kernel_v_gets_v_minus_A_times_w(bs,s,v,t+bs*vi[m]);
128135aa4fcfSShri Abhyankar       v += bs2;
128235aa4fcfSShri Abhyankar     }
128335aa4fcfSShri Abhyankar   }
128435aa4fcfSShri Abhyankar 
128535aa4fcfSShri Abhyankar   /* backward solve the upper triangular */
128635aa4fcfSShri Abhyankar   ls = a->solve_work + A->cmap->n;
128735aa4fcfSShri Abhyankar   for (i=n-1; i>=0; i--){
128835aa4fcfSShri Abhyankar     v  = aa + bs2*(adiag[i+1]+1);
128935aa4fcfSShri Abhyankar     vi = aj + adiag[i+1]+1;
129035aa4fcfSShri Abhyankar     nz = adiag[i] - adiag[i+1] - 1;
129135aa4fcfSShri Abhyankar     ierr = PetscMemcpy(ls,t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr);
129235aa4fcfSShri Abhyankar     for(m=0;m<nz;m++){
129335aa4fcfSShri Abhyankar       Kernel_v_gets_v_minus_A_times_w(bs,ls,v,t+bs*vi[m]);
129435aa4fcfSShri Abhyankar       v += bs2;
129535aa4fcfSShri Abhyankar     }
129635aa4fcfSShri Abhyankar     Kernel_w_gets_A_times_v(bs,ls,v,t+i*bs); /* *inv(diagonal[i]) */
129735aa4fcfSShri Abhyankar     ierr = PetscMemcpy(x + bs*c[i],t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr);
129835aa4fcfSShri Abhyankar   }
129935aa4fcfSShri Abhyankar   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
130035aa4fcfSShri Abhyankar   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
130135aa4fcfSShri Abhyankar   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
130235aa4fcfSShri Abhyankar   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
130335aa4fcfSShri Abhyankar   ierr = PetscLogFlops(2.0*(a->bs2)*(a->nz) - A->rmap->bs*A->cmap->n);CHKERRQ(ierr);
130435aa4fcfSShri Abhyankar   PetscFunctionReturn(0);
130535aa4fcfSShri Abhyankar }
130635aa4fcfSShri Abhyankar 
130735aa4fcfSShri Abhyankar #undef __FUNCT__
13082efa7f71SHong Zhang #define __FUNCT__ "BlockAbs_privat"
130991d91c02SMatthew Knepley PetscErrorCode BlockAbs_private(PetscInt nbs,PetscInt bs2,PetscScalar *blockarray,PetscReal *absarray)
13102efa7f71SHong Zhang {
13112efa7f71SHong Zhang   PetscErrorCode     ierr;
13122efa7f71SHong Zhang   PetscInt           i,j;
13132efa7f71SHong Zhang   PetscFunctionBegin;
13142efa7f71SHong Zhang   ierr = PetscMemzero(absarray,(nbs+1)*sizeof(PetscReal));CHKERRQ(ierr);
13152efa7f71SHong Zhang   for (i=0; i<nbs; i++){
13162efa7f71SHong Zhang     for (j=0; j<bs2; j++){
13172efa7f71SHong Zhang       if (absarray[i] < PetscAbsScalar(blockarray[i*nbs+j])) absarray[i] = PetscAbsScalar(blockarray[i*nbs+j]);
13182efa7f71SHong Zhang     }
13192efa7f71SHong Zhang   }
13202efa7f71SHong Zhang   PetscFunctionReturn(0);
13212efa7f71SHong Zhang }
13222efa7f71SHong Zhang 
1323c8342467SHong Zhang #undef __FUNCT__
1324c8342467SHong Zhang #define __FUNCT__ "MatILUDTFactor_SeqBAIJ"
1325fe97e370SBarry Smith /*
1326fe97e370SBarry Smith      This needs to be renamed and called by the regular MatILUFactor_SeqBAIJ when drop tolerance is used
1327fe97e370SBarry Smith */
1328c8342467SHong Zhang PetscErrorCode MatILUDTFactor_SeqBAIJ(Mat A,IS isrow,IS iscol,const MatFactorInfo *info,Mat *fact)
1329c8342467SHong Zhang {
13302efa7f71SHong Zhang   Mat                B = *fact;
13312efa7f71SHong Zhang   Mat_SeqBAIJ        *a=(Mat_SeqBAIJ*)A->data,*b;
13322efa7f71SHong Zhang   IS                 isicol;
13332efa7f71SHong Zhang   PetscErrorCode     ierr;
13342efa7f71SHong Zhang   const PetscInt     *r,*ic;
13352efa7f71SHong Zhang   PetscInt           i,mbs=a->mbs,bs=A->rmap->bs,bs2=a->bs2,*ai=a->i,*aj=a->j,*ajtmp,*adiag;
13362efa7f71SHong Zhang   PetscInt           *bi,*bj,*bdiag;
13372efa7f71SHong Zhang 
13382efa7f71SHong Zhang   PetscInt           row,nzi,nzi_bl,nzi_bu,*im,dtcount,nzi_al,nzi_au;
13392efa7f71SHong Zhang   PetscInt           nlnk,*lnk;
13402efa7f71SHong Zhang   PetscBT            lnkbt;
13412efa7f71SHong Zhang   PetscTruth         row_identity,icol_identity,both_identity;
13422efa7f71SHong Zhang   MatScalar          *aatmp,*pv,*batmp,*ba,*rtmp,*pc,*multiplier,*vtmp;
13432efa7f71SHong Zhang   const PetscInt     *ics;
13442efa7f71SHong Zhang   PetscInt           j,nz,*pj,*bjtmp,k,ncut,*jtmp;
13452efa7f71SHong Zhang 
13466da515a1SHong Zhang   PetscReal          dt=info->dt; /* shift=info->shiftinblocks; */
13472efa7f71SHong Zhang   PetscInt           nnz_max;
13482efa7f71SHong Zhang   PetscTruth         missing;
1349021822bcSHong Zhang   PetscReal          *vtmp_abs;
135097ef94ebSSatish Balay   MatScalar          *v_work;
135197ef94ebSSatish Balay   PetscInt           *v_pivots;
13522efa7f71SHong Zhang 
1353c8342467SHong Zhang   PetscFunctionBegin;
13542efa7f71SHong Zhang   /* ------- symbolic factorization, can be reused ---------*/
13552efa7f71SHong Zhang   ierr = MatMissingDiagonal(A,&missing,&i);CHKERRQ(ierr);
13562efa7f71SHong Zhang   if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i);
13572efa7f71SHong Zhang   adiag=a->diag;
13582efa7f71SHong Zhang 
13592efa7f71SHong Zhang   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
13602efa7f71SHong Zhang 
13612efa7f71SHong Zhang   /* bdiag is location of diagonal in factor */
13622efa7f71SHong Zhang   ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
13632efa7f71SHong Zhang 
13642efa7f71SHong Zhang   /* allocate row pointers bi */
13656bce7ff8SHong Zhang   ierr = PetscMalloc((2*mbs+2)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
13662efa7f71SHong Zhang 
13672efa7f71SHong Zhang   /* allocate bj and ba; max num of nonzero entries is (ai[n]+2*n*dtcount+2) */
13682efa7f71SHong Zhang   dtcount = (PetscInt)info->dtcount;
13696bce7ff8SHong Zhang   if (dtcount > mbs-1) dtcount = mbs-1;
13706bce7ff8SHong Zhang   nnz_max  = ai[mbs]+2*mbs*dtcount +2;
13716da515a1SHong Zhang   /* printf("MatILUDTFactor_SeqBAIJ, bs %d, ai[mbs] %d, nnz_max  %d, dtcount %d\n",bs,ai[mbs],nnz_max,dtcount); */
13722efa7f71SHong Zhang   ierr = PetscMalloc(nnz_max*sizeof(PetscInt),&bj);CHKERRQ(ierr);
13736bce7ff8SHong Zhang   nnz_max = nnz_max*bs2;
13742efa7f71SHong Zhang   ierr = PetscMalloc(nnz_max*sizeof(MatScalar),&ba);CHKERRQ(ierr);
13752efa7f71SHong Zhang 
13762efa7f71SHong Zhang   /* put together the new matrix */
13772efa7f71SHong Zhang   ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
13782efa7f71SHong Zhang   ierr = PetscLogObjectParent(B,isicol);CHKERRQ(ierr);
13792efa7f71SHong Zhang   b    = (Mat_SeqBAIJ*)(B)->data;
13802efa7f71SHong Zhang   b->free_a       = PETSC_TRUE;
13812efa7f71SHong Zhang   b->free_ij      = PETSC_TRUE;
13822efa7f71SHong Zhang   b->singlemalloc = PETSC_FALSE;
13832efa7f71SHong Zhang   b->a          = ba;
13842efa7f71SHong Zhang   b->j          = bj;
13852efa7f71SHong Zhang   b->i          = bi;
13862efa7f71SHong Zhang   b->diag       = bdiag;
13872efa7f71SHong Zhang   b->ilen       = 0;
13882efa7f71SHong Zhang   b->imax       = 0;
13892efa7f71SHong Zhang   b->row        = isrow;
13902efa7f71SHong Zhang   b->col        = iscol;
13912efa7f71SHong Zhang   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
13922efa7f71SHong Zhang   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
13932efa7f71SHong Zhang   b->icol       = isicol;
13942efa7f71SHong Zhang   ierr = PetscMalloc((bs*(mbs+1))*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
13952efa7f71SHong Zhang 
13962efa7f71SHong Zhang   ierr = PetscLogObjectMemory(B,nnz_max*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
13972efa7f71SHong Zhang   b->maxnz = nnz_max/bs2;
13982efa7f71SHong Zhang 
13992efa7f71SHong Zhang   (B)->factor                = MAT_FACTOR_ILUDT;
14002efa7f71SHong Zhang   (B)->info.factor_mallocs   = 0;
14012efa7f71SHong Zhang   (B)->info.fill_ratio_given = ((PetscReal)nnz_max)/((PetscReal)(ai[mbs]*bs2));
14022efa7f71SHong Zhang   CHKMEMQ;
14032efa7f71SHong Zhang   /* ------- end of symbolic factorization ---------*/
14042efa7f71SHong Zhang   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
14052efa7f71SHong Zhang   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
14062efa7f71SHong Zhang   ics  = ic;
14072efa7f71SHong Zhang 
14082efa7f71SHong Zhang   /* linked list for storing column indices of the active row */
14092efa7f71SHong Zhang   nlnk = mbs + 1;
14102efa7f71SHong Zhang   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
14112efa7f71SHong Zhang 
14122efa7f71SHong Zhang   /* im: used by PetscLLAddSortedLU(); jtmp: working array for column indices of active row */
1413fca92195SBarry Smith   ierr = PetscMalloc2(mbs,PetscInt,&im,mbs,PetscInt,&jtmp);CHKERRQ(ierr);
14142efa7f71SHong Zhang   /* rtmp, vtmp: working arrays for sparse and contiguous row entries of active row */
1415fca92195SBarry Smith   ierr = PetscMalloc2(mbs*bs2,MatScalar,&rtmp,mbs*bs2,MatScalar,&vtmp);CHKERRQ(ierr);
1416021822bcSHong Zhang   ierr = PetscMalloc((mbs+1)*sizeof(PetscReal),&vtmp_abs);CHKERRQ(ierr);
1417fca92195SBarry Smith   ierr = PetscMalloc3(bs,MatScalar,&v_work,bs2,MatScalar,&multiplier,bs,PetscInt,&v_pivots);CHKERRQ(ierr);
14182efa7f71SHong Zhang 
1419e60cf9a0SBarry Smith   bi[0]    = 0;
14202efa7f71SHong Zhang   bdiag[0] = (nnz_max/bs2)-1; /* location of diagonal in factor B */
14216bce7ff8SHong Zhang   bi[2*mbs+1] = bdiag[0]+1; /* endof bj and ba array */
14222efa7f71SHong Zhang   for (i=0; i<mbs; i++) {
14232efa7f71SHong Zhang     /* copy initial fill into linked list */
14242efa7f71SHong Zhang     nzi = 0; /* nonzeros for active row i */
14252efa7f71SHong Zhang     nzi = ai[r[i]+1] - ai[r[i]];
14262efa7f71SHong Zhang     if (!nzi) SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i);
14272efa7f71SHong Zhang     nzi_al = adiag[r[i]] - ai[r[i]];
14282efa7f71SHong Zhang     nzi_au = ai[r[i]+1] - adiag[r[i]] -1;
14296da515a1SHong Zhang     /* printf("row %d, nzi_al/au %d %d\n",i,nzi_al,nzi_au); */
14302efa7f71SHong Zhang 
14312efa7f71SHong Zhang     /* load in initial unfactored row */
14322efa7f71SHong Zhang     ajtmp = aj + ai[r[i]];
14332efa7f71SHong Zhang     ierr = PetscLLAddPerm(nzi,ajtmp,ic,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
1434fca92195SBarry Smith     ierr = PetscMemzero(rtmp,mbs*bs2*sizeof(PetscScalar));CHKERRQ(ierr);
14352efa7f71SHong Zhang     aatmp = a->a + bs2*ai[r[i]];
14362efa7f71SHong Zhang     for (j=0; j<nzi; j++) {
14372efa7f71SHong Zhang       ierr = PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],aatmp+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr);
14382efa7f71SHong Zhang     }
14392efa7f71SHong Zhang 
14402efa7f71SHong Zhang     /* add pivot rows into linked list */
14412efa7f71SHong Zhang     row = lnk[mbs];
14422efa7f71SHong Zhang     while (row < i) {
14432efa7f71SHong Zhang       nzi_bl = bi[row+1] - bi[row] + 1;
14442efa7f71SHong Zhang       bjtmp = bj + bdiag[row+1]+1; /* points to 1st column next to the diagonal in U */
14452efa7f71SHong Zhang       ierr  = PetscLLAddSortedLU(bjtmp,row,nlnk,lnk,lnkbt,i,nzi_bl,im);CHKERRQ(ierr);
14462efa7f71SHong Zhang       nzi  += nlnk;
14472efa7f71SHong Zhang       row   = lnk[row];
14482efa7f71SHong Zhang     }
14492efa7f71SHong Zhang 
14502efa7f71SHong Zhang     /* copy data from lnk into jtmp, then initialize lnk */
14512efa7f71SHong Zhang     ierr = PetscLLClean(mbs,mbs,nzi,lnk,jtmp,lnkbt);CHKERRQ(ierr);
14522efa7f71SHong Zhang 
14532efa7f71SHong Zhang     /* numerical factorization */
14542efa7f71SHong Zhang     bjtmp = jtmp;
14552efa7f71SHong Zhang     row   = *bjtmp++; /* 1st pivot row */
14562efa7f71SHong Zhang 
14572efa7f71SHong Zhang     while  (row < i) {
14582efa7f71SHong Zhang       pc = rtmp + bs2*row;
14592efa7f71SHong Zhang       pv = ba + bs2*bdiag[row]; /* inv(diag) of the pivot row */
14602efa7f71SHong Zhang       Kernel_A_gets_A_times_B(bs,pc,pv,multiplier); /* pc= multiplier = pc*inv(diag[row]) */
14612efa7f71SHong Zhang       ierr = BlockAbs_private(1,bs2,pc,vtmp_abs);CHKERRQ(ierr);
14622efa7f71SHong Zhang       if (vtmp_abs[0] > dt){ /* apply tolerance dropping rule */
14632efa7f71SHong Zhang         pj         = bj + bdiag[row+1] + 1; /* point to 1st entry of U(row,:) */
14642efa7f71SHong Zhang         pv         = ba + bs2*(bdiag[row+1] + 1);
14652efa7f71SHong Zhang         nz         = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:), excluding diagonal */
14662efa7f71SHong Zhang         for (j=0; j<nz; j++){
14672efa7f71SHong Zhang           Kernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j);
14682efa7f71SHong Zhang         }
14696da515a1SHong Zhang         /* ierr = PetscLogFlops(bslog*(nz+1.0)-bs);CHKERRQ(ierr); */
14702efa7f71SHong Zhang       }
14712efa7f71SHong Zhang       row = *bjtmp++;
14722efa7f71SHong Zhang     }
14732efa7f71SHong Zhang 
14742efa7f71SHong Zhang     /* copy sparse rtmp into contiguous vtmp; separate L and U part */
14752efa7f71SHong Zhang     nzi_bl = 0; j = 0;
14762efa7f71SHong Zhang     while (jtmp[j] < i){ /* L-part. Note: jtmp is sorted */
14772efa7f71SHong Zhang       ierr = PetscMemcpy(vtmp+bs2*j,rtmp+bs2*jtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
14782efa7f71SHong Zhang       nzi_bl++; j++;
14792efa7f71SHong Zhang     }
14802efa7f71SHong Zhang     nzi_bu = nzi - nzi_bl -1;
14816da515a1SHong Zhang     /* printf("nzi %d, nzi_bl %d, nzi_bu %d\n",nzi,nzi_bl,nzi_bu); */
14822efa7f71SHong Zhang 
14832efa7f71SHong Zhang     while (j < nzi){ /* U-part */
14842efa7f71SHong Zhang       ierr = PetscMemcpy(vtmp+bs2*j,rtmp+bs2*jtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
14852efa7f71SHong Zhang       /*
14862efa7f71SHong Zhang       printf(" col %d: ",jtmp[j]);
14872efa7f71SHong Zhang       for (j1=0; j1<bs2; j1++) printf(" %g",*(vtmp+bs2*j+j1));
14882efa7f71SHong Zhang       printf(" \n");
14892efa7f71SHong Zhang       */
14902efa7f71SHong Zhang       j++;
14912efa7f71SHong Zhang     }
14922efa7f71SHong Zhang 
14932efa7f71SHong Zhang     ierr = BlockAbs_private(nzi,bs2,vtmp,vtmp_abs);CHKERRQ(ierr);
14942efa7f71SHong Zhang     /*
14952efa7f71SHong Zhang     printf(" row %d, nzi %d, vtmp_abs\n",i,nzi);
14962efa7f71SHong Zhang     for (j1=0; j1<nzi; j1++) printf(" (%d %g),",jtmp[j1],vtmp_abs[j1]);
14972efa7f71SHong Zhang     printf(" \n");
14982efa7f71SHong Zhang     */
14992efa7f71SHong Zhang     bjtmp = bj + bi[i];
15002efa7f71SHong Zhang     batmp = ba + bs2*bi[i];
15012efa7f71SHong Zhang     /* apply level dropping rule to L part */
15022efa7f71SHong Zhang     ncut = nzi_al + dtcount;
15032efa7f71SHong Zhang     if (ncut < nzi_bl){
1504021822bcSHong Zhang       ierr = PetscSortSplitReal(ncut,nzi_bl,vtmp_abs,jtmp);CHKERRQ(ierr);
15052efa7f71SHong Zhang       ierr = PetscSortIntWithScalarArray(ncut,jtmp,vtmp);CHKERRQ(ierr);
15062efa7f71SHong Zhang     } else {
15072efa7f71SHong Zhang       ncut = nzi_bl;
15082efa7f71SHong Zhang     }
15092efa7f71SHong Zhang     for (j=0; j<ncut; j++){
15102efa7f71SHong Zhang       bjtmp[j] = jtmp[j];
15112efa7f71SHong Zhang       ierr = PetscMemcpy(batmp+bs2*j,rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
15122efa7f71SHong Zhang       /*
15132efa7f71SHong Zhang       printf(" col %d: ",bjtmp[j]);
15142efa7f71SHong Zhang       for (j1=0; j1<bs2; j1++) printf(" %g,",*(batmp+bs2*j+j1));
15152efa7f71SHong Zhang       printf("\n");
15162efa7f71SHong Zhang       */
15172efa7f71SHong Zhang     }
15182efa7f71SHong Zhang     bi[i+1] = bi[i] + ncut;
15192efa7f71SHong Zhang     nzi = ncut + 1;
15202efa7f71SHong Zhang 
15212efa7f71SHong Zhang     /* apply level dropping rule to U part */
15222efa7f71SHong Zhang     ncut = nzi_au + dtcount;
15232efa7f71SHong Zhang     if (ncut < nzi_bu){
1524021822bcSHong Zhang       ierr = PetscSortSplitReal(ncut,nzi_bu,vtmp_abs+nzi_bl+1,jtmp+nzi_bl+1);CHKERRQ(ierr);
15252efa7f71SHong Zhang       ierr = PetscSortIntWithScalarArray(ncut,jtmp+nzi_bl+1,vtmp+nzi_bl+1);CHKERRQ(ierr);
15262efa7f71SHong Zhang     } else {
15272efa7f71SHong Zhang       ncut = nzi_bu;
15282efa7f71SHong Zhang     }
15292efa7f71SHong Zhang     nzi += ncut;
15302efa7f71SHong Zhang 
15312efa7f71SHong Zhang     /* mark bdiagonal */
15322efa7f71SHong Zhang     bdiag[i+1]    = bdiag[i] - (ncut + 1);
15336bce7ff8SHong Zhang     bi[2*mbs - i] = bi[2*mbs - i +1] - (ncut + 1);
15346bce7ff8SHong Zhang 
15352efa7f71SHong Zhang     bjtmp = bj + bdiag[i];
15362efa7f71SHong Zhang     batmp = ba + bs2*bdiag[i];
15372efa7f71SHong Zhang     ierr = PetscMemcpy(batmp,rtmp+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
15382efa7f71SHong Zhang     *bjtmp = i;
15392efa7f71SHong Zhang     /*
15402efa7f71SHong Zhang     printf(" diag %d: ",*bjtmp);
15412efa7f71SHong Zhang     for (j=0; j<bs2; j++){
15422efa7f71SHong Zhang       printf(" %g,",batmp[j]);
15432efa7f71SHong Zhang     }
15442efa7f71SHong Zhang     printf("\n");
15452efa7f71SHong Zhang     */
15462efa7f71SHong Zhang     bjtmp = bj + bdiag[i+1]+1;
15472efa7f71SHong Zhang     batmp = ba + (bdiag[i+1]+1)*bs2;
15482efa7f71SHong Zhang 
15492efa7f71SHong Zhang     for (k=0; k<ncut; k++){
15502efa7f71SHong Zhang       bjtmp[k] = jtmp[nzi_bl+1+k];
15512efa7f71SHong Zhang       ierr = PetscMemcpy(batmp+bs2*k,rtmp+bs2*bjtmp[k],bs2*sizeof(MatScalar));CHKERRQ(ierr);
15522efa7f71SHong Zhang       /*
15532efa7f71SHong Zhang       printf(" col %d:",bjtmp[k]);
15542efa7f71SHong Zhang       for (j1=0; j1<bs2; j1++) printf(" %g,",*(batmp+bs2*k+j1));
15552efa7f71SHong Zhang       printf("\n");
15562efa7f71SHong Zhang       */
15572efa7f71SHong Zhang     }
15582efa7f71SHong Zhang 
15592efa7f71SHong Zhang     im[i] = nzi; /* used by PetscLLAddSortedLU() */
15602efa7f71SHong Zhang 
15612efa7f71SHong Zhang     /* invert diagonal block for simplier triangular solves - add shift??? */
15622efa7f71SHong Zhang     batmp = ba + bs2*bdiag[i];
15632efa7f71SHong Zhang     ierr = Kernel_A_gets_inverse_A(bs,batmp,v_pivots,v_work);CHKERRQ(ierr);
15642efa7f71SHong Zhang   } /* for (i=0; i<mbs; i++) */
1565fca92195SBarry Smith   ierr = PetscFree3(v_work,multiplier,v_pivots);CHKERRQ(ierr);
15662efa7f71SHong Zhang 
15676da515a1SHong Zhang   /* printf("end of L %d, beginning of U %d\n",bi[mbs],bdiag[mbs]); */
15682efa7f71SHong Zhang   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]);
15692efa7f71SHong Zhang 
15702efa7f71SHong Zhang   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
15712efa7f71SHong Zhang   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
15722efa7f71SHong Zhang 
15732efa7f71SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
15742efa7f71SHong Zhang 
1575fca92195SBarry Smith   ierr = PetscFree2(im,jtmp);CHKERRQ(ierr);
1576fca92195SBarry Smith   ierr = PetscFree2(rtmp,vtmp);CHKERRQ(ierr);
15772efa7f71SHong Zhang 
15782efa7f71SHong Zhang   ierr = PetscLogFlops(bs2*B->cmap->n);CHKERRQ(ierr);
15792efa7f71SHong Zhang   b->maxnz = b->nz = bi[mbs] + bdiag[0] - bdiag[mbs];
15802efa7f71SHong Zhang 
15812efa7f71SHong Zhang   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
15822efa7f71SHong Zhang   ierr = ISIdentity(isicol,&icol_identity);CHKERRQ(ierr);
15832efa7f71SHong Zhang   both_identity = (PetscTruth) (row_identity && icol_identity);
15842efa7f71SHong Zhang   if (row_identity && icol_identity) {
15854dd39f65SShri Abhyankar     B->ops->solve = MatSolve_SeqBAIJ_N_NaturalOrdering;
15862efa7f71SHong Zhang   } else {
15874dd39f65SShri Abhyankar     B->ops->solve = MatSolve_SeqBAIJ_N;
15882efa7f71SHong Zhang   }
15892efa7f71SHong Zhang 
15902efa7f71SHong Zhang   B->ops->solveadd          = 0;
15912efa7f71SHong Zhang   B->ops->solvetranspose    = 0;
15922efa7f71SHong Zhang   B->ops->solvetransposeadd = 0;
15932efa7f71SHong Zhang   B->ops->matsolve          = 0;
15942efa7f71SHong Zhang   B->assembled              = PETSC_TRUE;
15952efa7f71SHong Zhang   B->preallocated           = PETSC_TRUE;
1596c8342467SHong Zhang   PetscFunctionReturn(0);
1597c8342467SHong Zhang }
1598