xref: /petsc/src/mat/impls/baij/seq/baijfact.c (revision 0a3c351a51b8afeba0e28ebe856db6cb3362db0e)
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__
10b588c5a2SHong Zhang #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2_newdatastruct"
11b588c5a2SHong Zhang PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_newdatastruct(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;
18b588c5a2SHong Zhang   PetscInt       i,j,k,n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
19b588c5a2SHong Zhang   PetscInt       *ajtmp,*bjtmp,nz,nzL,row,*bdiag=b->diag,*pj;
20b588c5a2SHong Zhang   MatScalar      *rtmp,*pc,*mwork,*v,*pv,*aa=a->a;
214c000e2eSHong Zhang   PetscInt       bs2 = a->bs2,flg;
22b588c5a2SHong Zhang   PetscReal      shift = info->shiftinblocks;
23b588c5a2SHong Zhang 
24b588c5a2SHong Zhang   PetscFunctionBegin;
25b588c5a2SHong Zhang   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
26b588c5a2SHong Zhang   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
27b588c5a2SHong Zhang 
284c000e2eSHong Zhang   /* generate work space needed by the factorization */
29fca92195SBarry Smith   ierr = PetscMalloc2(bs2*n,MatScalar,&rtmp,bs2,MatScalar,&mwork);CHKERRQ(ierr);
304c000e2eSHong Zhang   ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr);
31b588c5a2SHong Zhang   ics  = ic;
32b588c5a2SHong Zhang 
33b588c5a2SHong Zhang   for (i=0; i<n; i++){
34b588c5a2SHong Zhang     /* zero rtmp */
35b588c5a2SHong Zhang     /* L part */
36b588c5a2SHong Zhang     nz    = bi[i+1] - bi[i];
37b588c5a2SHong Zhang     bjtmp = bj + bi[i];
38b588c5a2SHong Zhang     for  (j=0; j<nz; j++){
39b588c5a2SHong Zhang       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
40b588c5a2SHong Zhang     }
41b588c5a2SHong Zhang 
42b588c5a2SHong Zhang     /* U part */
430c4413a7SShri Abhyankar     nz = bdiag[i] - bdiag[i+1];
440c4413a7SShri Abhyankar     bjtmp = bj + bdiag[i+1]+1;
450c4413a7SShri Abhyankar     for  (j=0; j<nz; j++){
460c4413a7SShri Abhyankar       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
470c4413a7SShri Abhyankar     }
480c4413a7SShri Abhyankar 
490c4413a7SShri Abhyankar     /* load in initial (unfactored row) */
500c4413a7SShri Abhyankar     nz    = ai[r[i]+1] - ai[r[i]];
510c4413a7SShri Abhyankar     ajtmp = aj + ai[r[i]];
520c4413a7SShri Abhyankar     v     = aa + bs2*ai[r[i]];
530c4413a7SShri Abhyankar     for (j=0; j<nz; j++) {
540c4413a7SShri Abhyankar       ierr = PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr);
550c4413a7SShri Abhyankar     }
560c4413a7SShri Abhyankar 
570c4413a7SShri Abhyankar     /* elimination */
580c4413a7SShri Abhyankar     bjtmp = bj + bi[i];
590c4413a7SShri Abhyankar     nzL   = bi[i+1] - bi[i];
600c4413a7SShri Abhyankar     for(k=0;k < nzL;k++) {
610c4413a7SShri Abhyankar       row = bjtmp[k];
620c4413a7SShri Abhyankar       pc = rtmp + bs2*row;
630c4413a7SShri Abhyankar       for (flg=0,j=0; j<bs2; j++) { if (pc[j]!=0.0) { flg = 1; break; }}
640c4413a7SShri Abhyankar       if (flg) {
650c4413a7SShri Abhyankar         pv = b->a + bs2*bdiag[row];
660c4413a7SShri Abhyankar         /* Kernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
670c4413a7SShri Abhyankar         ierr = Kernel_A_gets_A_times_B_2(pc,pv,mwork);CHKERRQ(ierr);
680c4413a7SShri Abhyankar 
690c4413a7SShri Abhyankar         pj = b->j + bdiag[row+1]+1; /* begining of U(row,:) */
700c4413a7SShri Abhyankar         pv = b->a + bs2*(bdiag[row+1]+1);
710c4413a7SShri Abhyankar         nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries inU(row,:), excluding diag */
720c4413a7SShri Abhyankar         for (j=0; j<nz; j++) {
730c4413a7SShri Abhyankar           /* Kernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
740c4413a7SShri Abhyankar           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
750c4413a7SShri Abhyankar           v    = rtmp + 4*pj[j];
760c4413a7SShri Abhyankar           ierr = Kernel_A_gets_A_minus_B_times_C_2(v,pc,pv);CHKERRQ(ierr);
770c4413a7SShri Abhyankar           pv  += 4;
780c4413a7SShri Abhyankar         }
790c4413a7SShri Abhyankar         ierr = PetscLogFlops(16*nz+12);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
800c4413a7SShri Abhyankar       }
810c4413a7SShri Abhyankar     }
820c4413a7SShri Abhyankar 
830c4413a7SShri Abhyankar     /* finished row so stick it into b->a */
840c4413a7SShri Abhyankar     /* L part */
850c4413a7SShri Abhyankar     pv   = b->a + bs2*bi[i] ;
860c4413a7SShri Abhyankar     pj   = b->j + bi[i] ;
870c4413a7SShri Abhyankar     nz   = bi[i+1] - bi[i];
880c4413a7SShri Abhyankar     for (j=0; j<nz; j++) {
890c4413a7SShri Abhyankar       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
900c4413a7SShri Abhyankar     }
910c4413a7SShri Abhyankar 
920c4413a7SShri Abhyankar     /* Mark diagonal and invert diagonal for simplier triangular solves */
930c4413a7SShri Abhyankar     pv   = b->a + bs2*bdiag[i];
940c4413a7SShri Abhyankar     pj   = b->j + bdiag[i];
950c4413a7SShri Abhyankar     ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr);
960c4413a7SShri Abhyankar     /* ierr = Kernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); */
970c4413a7SShri Abhyankar     ierr = Kernel_A_gets_inverse_A_2(pv,shift);CHKERRQ(ierr);
980c4413a7SShri Abhyankar 
990c4413a7SShri Abhyankar     /* U part */
1000c4413a7SShri Abhyankar     pv = b->a + bs2*(bdiag[i+1]+1);
1010c4413a7SShri Abhyankar     pj = b->j + bdiag[i+1]+1;
1020c4413a7SShri Abhyankar     nz = bdiag[i] - bdiag[i+1] - 1;
1030c4413a7SShri Abhyankar     for (j=0; j<nz; j++){
1040c4413a7SShri Abhyankar       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
1050c4413a7SShri Abhyankar     }
1060c4413a7SShri Abhyankar   }
1070c4413a7SShri Abhyankar 
108fca92195SBarry Smith   ierr = PetscFree2(rtmp,mwork);CHKERRQ(ierr);
1090c4413a7SShri Abhyankar   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
1100c4413a7SShri Abhyankar   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
111a2d6a19aSShri Abhyankar   C->ops->solve          = MatSolve_SeqBAIJ_2_newdatastruct;
11232121132SShri Abhyankar   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_newdatastruct;
1130c4413a7SShri Abhyankar 
1140c4413a7SShri Abhyankar   C->assembled = PETSC_TRUE;
1150c4413a7SShri Abhyankar   ierr = PetscLogFlops(1.3333*bs2*n);CHKERRQ(ierr); /* from inverting diagonal blocks */
1160c4413a7SShri Abhyankar   PetscFunctionReturn(0);
1170c4413a7SShri Abhyankar }
1180c4413a7SShri Abhyankar 
1194c000e2eSHong Zhang #undef __FUNCT__
1204c000e2eSHong Zhang #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_newdatastruct"
1214c000e2eSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_newdatastruct(Mat B,Mat A,const MatFactorInfo *info)
1224c000e2eSHong Zhang {
1234c000e2eSHong Zhang   Mat            C=B;
1244c000e2eSHong Zhang   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ *)C->data;
1254c000e2eSHong Zhang   PetscErrorCode ierr;
1264c000e2eSHong Zhang   PetscInt       i,j,k,n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
1274c000e2eSHong Zhang   PetscInt       *ajtmp,*bjtmp,nz,nzL,row,*bdiag=b->diag,*pj;
1284c000e2eSHong Zhang   MatScalar      *rtmp,*pc,*mwork,*v,*pv,*aa=a->a;
1294c000e2eSHong Zhang   PetscInt       bs2 = a->bs2,flg;
1304c000e2eSHong Zhang   PetscReal      shift = info->shiftinblocks;
1314c000e2eSHong Zhang 
1324c000e2eSHong Zhang   PetscFunctionBegin;
1334c000e2eSHong Zhang   /* generate work space needed by the factorization */
134fca92195SBarry Smith   ierr = PetscMalloc2(bs2*n,MatScalar,&rtmp,bs2,MatScalar,&mwork);CHKERRQ(ierr);
1354c000e2eSHong Zhang   ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr);
1364c000e2eSHong Zhang 
1374c000e2eSHong Zhang   for (i=0; i<n; i++){
1384c000e2eSHong Zhang     /* zero rtmp */
1394c000e2eSHong Zhang     /* L part */
1404c000e2eSHong Zhang     nz    = bi[i+1] - bi[i];
1414c000e2eSHong Zhang     bjtmp = bj + bi[i];
1424c000e2eSHong Zhang     for  (j=0; j<nz; j++){
1434c000e2eSHong Zhang       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
1444c000e2eSHong Zhang     }
1454c000e2eSHong Zhang 
1464c000e2eSHong Zhang     /* U part */
147b2b2dd24SShri Abhyankar     nz = bdiag[i] - bdiag[i+1];
148b2b2dd24SShri Abhyankar     bjtmp = bj + bdiag[i+1]+1;
149b2b2dd24SShri Abhyankar     for  (j=0; j<nz; j++){
150b2b2dd24SShri Abhyankar       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
151b2b2dd24SShri Abhyankar     }
152b2b2dd24SShri Abhyankar 
153b2b2dd24SShri Abhyankar     /* load in initial (unfactored row) */
154b2b2dd24SShri Abhyankar     nz    = ai[i+1] - ai[i];
155b2b2dd24SShri Abhyankar     ajtmp = aj + ai[i];
156b2b2dd24SShri Abhyankar     v     = aa + bs2*ai[i];
157b2b2dd24SShri Abhyankar     for (j=0; j<nz; j++) {
158b2b2dd24SShri Abhyankar       ierr = PetscMemcpy(rtmp+bs2*ajtmp[j],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr);
159b2b2dd24SShri Abhyankar     }
160b2b2dd24SShri Abhyankar 
161b2b2dd24SShri Abhyankar     /* elimination */
162b2b2dd24SShri Abhyankar     bjtmp = bj + bi[i];
163b2b2dd24SShri Abhyankar     nzL   = bi[i+1] - bi[i];
164b2b2dd24SShri Abhyankar     for(k=0;k < nzL;k++) {
165b2b2dd24SShri Abhyankar       row = bjtmp[k];
166b2b2dd24SShri Abhyankar       pc = rtmp + bs2*row;
167b2b2dd24SShri Abhyankar       for (flg=0,j=0; j<bs2; j++) { if (pc[j]!=0.0) { flg = 1; break; }}
168b2b2dd24SShri Abhyankar       if (flg) {
169b2b2dd24SShri Abhyankar         pv = b->a + bs2*bdiag[row];
170b2b2dd24SShri Abhyankar         /* Kernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
171b2b2dd24SShri Abhyankar         ierr = Kernel_A_gets_A_times_B_2(pc,pv,mwork);CHKERRQ(ierr);
172b2b2dd24SShri Abhyankar 
173b2b2dd24SShri Abhyankar         pj = b->j + bdiag[row+1]+1; /* beginning of U(row,:) */
174b2b2dd24SShri Abhyankar 	pv = b->a + bs2*(bdiag[row+1]+1);
175b2b2dd24SShri Abhyankar 	nz = bdiag[row]-bdiag[row+1] - 1; /* num of entries in U(row,:) excluding diag */
176b2b2dd24SShri Abhyankar         for (j=0; j<nz; j++) {
177b2b2dd24SShri Abhyankar           /* Kernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
178b2b2dd24SShri Abhyankar           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
179b2b2dd24SShri Abhyankar           v    = rtmp + 4*pj[j];
180b2b2dd24SShri Abhyankar           ierr = Kernel_A_gets_A_minus_B_times_C_2(v,pc,pv);CHKERRQ(ierr);
181b2b2dd24SShri Abhyankar           pv  += 4;
182b2b2dd24SShri Abhyankar         }
183b2b2dd24SShri Abhyankar         ierr = PetscLogFlops(16*nz+12);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
184b2b2dd24SShri Abhyankar       }
185b2b2dd24SShri Abhyankar     }
186b2b2dd24SShri Abhyankar 
187b2b2dd24SShri Abhyankar     /* finished row so stick it into b->a */
188b2b2dd24SShri Abhyankar     /* L part */
189b2b2dd24SShri Abhyankar     pv   = b->a + bs2*bi[i] ;
190b2b2dd24SShri Abhyankar     pj   = b->j + bi[i] ;
191b2b2dd24SShri Abhyankar     nz   = bi[i+1] - bi[i];
192b2b2dd24SShri Abhyankar     for (j=0; j<nz; j++) {
193b2b2dd24SShri Abhyankar       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
194b2b2dd24SShri Abhyankar     }
195b2b2dd24SShri Abhyankar 
196b2b2dd24SShri Abhyankar     /* Mark diagonal and invert diagonal for simplier triangular solves */
197b2b2dd24SShri Abhyankar     pv   = b->a + bs2*bdiag[i];
198b2b2dd24SShri Abhyankar     pj   = b->j + bdiag[i];
199b2b2dd24SShri Abhyankar     ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr);
200b2b2dd24SShri Abhyankar     /* ierr = Kernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); */
201b2b2dd24SShri Abhyankar     ierr = Kernel_A_gets_inverse_A_2(pv,shift);CHKERRQ(ierr);
202b2b2dd24SShri Abhyankar 
203b2b2dd24SShri Abhyankar     /* U part */
204b2b2dd24SShri Abhyankar     /*
205b2b2dd24SShri Abhyankar     pv = b->a + bs2*bi[2*n-i];
206b2b2dd24SShri Abhyankar     pj = b->j + bi[2*n-i];
207b2b2dd24SShri Abhyankar     nz = bi[2*n-i+1] - bi[2*n-i] - 1;
208b2b2dd24SShri Abhyankar     */
209b2b2dd24SShri Abhyankar     pv = b->a + bs2*(bdiag[i+1]+1);
210b2b2dd24SShri Abhyankar     pj = b->j + bdiag[i+1]+1;
211b2b2dd24SShri Abhyankar     nz = bdiag[i] - bdiag[i+1] - 1;
212b2b2dd24SShri Abhyankar     for (j=0; j<nz; j++){
213b2b2dd24SShri Abhyankar       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
214b2b2dd24SShri Abhyankar     }
215b2b2dd24SShri Abhyankar   }
216fca92195SBarry Smith   ierr = PetscFree2(rtmp,mwork);CHKERRQ(ierr);
217ae3d28f0SHong Zhang 
218a2d6a19aSShri Abhyankar   C->ops->solve          = MatSolve_SeqBAIJ_2_NaturalOrdering_newdatastruct;
2196929473cSShri Abhyankar   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_NaturalOrdering_newdatastruct;
220b2b2dd24SShri Abhyankar   C->assembled = PETSC_TRUE;
221b2b2dd24SShri Abhyankar   ierr = PetscLogFlops(1.3333*bs2*n);CHKERRQ(ierr); /* from inverting diagonal blocks */
222b2b2dd24SShri Abhyankar   PetscFunctionReturn(0);
223b2b2dd24SShri Abhyankar }
224b2b2dd24SShri Abhyankar 
225b2b2dd24SShri Abhyankar #undef __FUNCT__
22606e38f1dSHong Zhang #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2_inplace"
22706e38f1dSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_inplace(Mat B,Mat A,const MatFactorInfo *info)
2284eeb42bcSBarry Smith {
229719d5645SBarry Smith   Mat            C = B;
2304eeb42bcSBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data;
2317cdf2dbbSSatish Balay   IS             isrow = b->row,isicol = b->icol;
2326849ba73SBarry Smith   PetscErrorCode ierr;
2335d0c19d7SBarry Smith   const PetscInt *r,*ic;
2345d0c19d7SBarry Smith   PetscInt       i,j,n = a->mbs,*bi = b->i,*bj = b->j;
235690b6cddSBarry Smith   PetscInt       *ajtmpold,*ajtmp,nz,row;
236690b6cddSBarry Smith   PetscInt       *diag_offset=b->diag,idx,*ai=a->i,*aj=a->j,*pj;
237329f5518SBarry Smith   MatScalar      *pv,*v,*rtmp,m1,m2,m3,m4,*pc,*w,*x,x1,x2,x3,x4;
2382515f8d2SSatish Balay   MatScalar      p1,p2,p3,p4;
2393f1db9ecSBarry Smith   MatScalar      *ba = b->a,*aa = a->a;
24062bba022SBarry Smith   PetscReal      shift = info->shiftinblocks;
2414eeb42bcSBarry Smith 
2423a40ed3dSBarry Smith   PetscFunctionBegin;
2434eeb42bcSBarry Smith   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
2444eeb42bcSBarry Smith   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
245b0a32e0cSBarry Smith   ierr  = PetscMalloc(4*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
2464eeb42bcSBarry Smith 
2474eeb42bcSBarry Smith   for (i=0; i<n; i++) {
2484078e994SBarry Smith     nz    = bi[i+1] - bi[i];
2494078e994SBarry Smith     ajtmp = bj + bi[i];
2504eeb42bcSBarry Smith     for  (j=0; j<nz; j++) {
2514eeb42bcSBarry Smith       x = rtmp+4*ajtmp[j]; x[0] = x[1] = x[2] = x[3] = 0.0;
2524eeb42bcSBarry Smith     }
2534eeb42bcSBarry Smith     /* load in initial (unfactored row) */
2544eeb42bcSBarry Smith     idx      = r[i];
2554078e994SBarry Smith     nz       = ai[idx+1] - ai[idx];
2564078e994SBarry Smith     ajtmpold = aj + ai[idx];
2574078e994SBarry Smith     v        = aa + 4*ai[idx];
2584eeb42bcSBarry Smith     for (j=0; j<nz; j++) {
2594eeb42bcSBarry Smith       x    = rtmp+4*ic[ajtmpold[j]];
2604eeb42bcSBarry Smith       x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3];
2614eeb42bcSBarry Smith       v    += 4;
2624eeb42bcSBarry Smith     }
2634eeb42bcSBarry Smith     row = *ajtmp++;
2644eeb42bcSBarry Smith     while (row < i) {
2654eeb42bcSBarry Smith       pc = rtmp + 4*row;
2664eeb42bcSBarry Smith       p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3];
26788685aaeSLois Curfman McInnes       if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0) {
2684078e994SBarry Smith         pv = ba + 4*diag_offset[row];
2694078e994SBarry Smith         pj = bj + diag_offset[row] + 1;
2704eeb42bcSBarry Smith         x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
2714eeb42bcSBarry Smith         pc[0] = m1 = p1*x1 + p3*x2;
2724eeb42bcSBarry Smith         pc[1] = m2 = p2*x1 + p4*x2;
2734eeb42bcSBarry Smith         pc[2] = m3 = p1*x3 + p3*x4;
2744eeb42bcSBarry Smith         pc[3] = m4 = p2*x3 + p4*x4;
2754078e994SBarry Smith         nz = bi[row+1] - diag_offset[row] - 1;
2764eeb42bcSBarry Smith         pv += 4;
2774eeb42bcSBarry Smith         for (j=0; j<nz; j++) {
2784eeb42bcSBarry Smith           x1   = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
2794eeb42bcSBarry Smith           x    = rtmp + 4*pj[j];
2804eeb42bcSBarry Smith           x[0] -= m1*x1 + m3*x2;
2814eeb42bcSBarry Smith           x[1] -= m2*x1 + m4*x2;
2824eeb42bcSBarry Smith           x[2] -= m1*x3 + m3*x4;
2834eeb42bcSBarry Smith           x[3] -= m2*x3 + m4*x4;
2844eeb42bcSBarry Smith           pv   += 4;
2854eeb42bcSBarry Smith         }
286dc0b31edSSatish Balay         ierr = PetscLogFlops(16.0*nz+12.0);CHKERRQ(ierr);
2874eeb42bcSBarry Smith       }
2884eeb42bcSBarry Smith       row = *ajtmp++;
2894eeb42bcSBarry Smith     }
2904eeb42bcSBarry Smith     /* finished row so stick it into b->a */
2914078e994SBarry Smith     pv = ba + 4*bi[i];
2924078e994SBarry Smith     pj = bj + bi[i];
2934078e994SBarry Smith     nz = bi[i+1] - bi[i];
2944eeb42bcSBarry Smith     for (j=0; j<nz; j++) {
2954eeb42bcSBarry Smith       x     = rtmp+4*pj[j];
2964eeb42bcSBarry Smith       pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3];
2974eeb42bcSBarry Smith       pv   += 4;
2984eeb42bcSBarry Smith     }
2994eeb42bcSBarry Smith     /* invert diagonal block */
3004078e994SBarry Smith     w = ba + 4*diag_offset[i];
30162bba022SBarry Smith     ierr = Kernel_A_gets_inverse_A_2(w,shift);CHKERRQ(ierr);
3024eeb42bcSBarry Smith   }
3034eeb42bcSBarry Smith 
304606d414cSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
3054eeb42bcSBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
3064eeb42bcSBarry Smith   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
30706e38f1dSHong Zhang   C->ops->solve          = MatSolve_SeqBAIJ_2_inplace;
30806e38f1dSHong Zhang   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_inplace;
3094eeb42bcSBarry Smith   C->assembled = PETSC_TRUE;
310efee365bSSatish Balay   ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
3113a40ed3dSBarry Smith   PetscFunctionReturn(0);
3124eeb42bcSBarry Smith }
3134cd38560SBarry Smith /*
3144cd38560SBarry Smith       Version for when blocks are 2 by 2 Using natural ordering
3154cd38560SBarry Smith */
3164a2ae208SSatish Balay #undef __FUNCT__
31706e38f1dSHong Zhang #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_inplace"
31806e38f1dSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_inplace(Mat C,Mat A,const MatFactorInfo *info)
3194cd38560SBarry Smith {
3204cd38560SBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data;
321dfbe8321SBarry Smith   PetscErrorCode ierr;
322690b6cddSBarry Smith   PetscInt       i,j,n = a->mbs,*bi = b->i,*bj = b->j;
323690b6cddSBarry Smith   PetscInt       *ajtmpold,*ajtmp,nz,row;
324690b6cddSBarry Smith   PetscInt       *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj;
325329f5518SBarry Smith   MatScalar      *pv,*v,*rtmp,*pc,*w,*x;
3262515f8d2SSatish Balay   MatScalar      p1,p2,p3,p4,m1,m2,m3,m4,x1,x2,x3,x4;
3274cd38560SBarry Smith   MatScalar      *ba = b->a,*aa = a->a;
32862bba022SBarry Smith   PetscReal      shift = info->shiftinblocks;
3294cd38560SBarry Smith 
3304cd38560SBarry Smith   PetscFunctionBegin;
331b0a32e0cSBarry Smith   ierr = PetscMalloc(4*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
3324cd38560SBarry Smith   for (i=0; i<n; i++) {
3334cd38560SBarry Smith     nz    = bi[i+1] - bi[i];
3344cd38560SBarry Smith     ajtmp = bj + bi[i];
3354cd38560SBarry Smith     for  (j=0; j<nz; j++) {
3364cd38560SBarry Smith       x = rtmp+4*ajtmp[j];
3374cd38560SBarry Smith       x[0]  = x[1]  = x[2]  = x[3]  = 0.0;
3384cd38560SBarry Smith     }
3394cd38560SBarry Smith     /* load in initial (unfactored row) */
3404cd38560SBarry Smith     nz       = ai[i+1] - ai[i];
3414cd38560SBarry Smith     ajtmpold = aj + ai[i];
3424cd38560SBarry Smith     v        = aa + 4*ai[i];
3434cd38560SBarry Smith     for (j=0; j<nz; j++) {
3444cd38560SBarry Smith       x    = rtmp+4*ajtmpold[j];
3454cd38560SBarry Smith       x[0]  = v[0];  x[1]  = v[1];  x[2]  = v[2];  x[3]  = v[3];
3464cd38560SBarry Smith       v    += 4;
3474cd38560SBarry Smith     }
3484cd38560SBarry Smith     row = *ajtmp++;
3494cd38560SBarry Smith     while (row < i) {
3504cd38560SBarry Smith       pc  = rtmp + 4*row;
3514cd38560SBarry Smith       p1  = pc[0];  p2  = pc[1];  p3  = pc[2];  p4  = pc[3];
3524cd38560SBarry Smith       if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0) {
3534cd38560SBarry Smith         pv = ba + 4*diag_offset[row];
3544cd38560SBarry Smith         pj = bj + diag_offset[row] + 1;
3554cd38560SBarry Smith         x1  = pv[0];  x2  = pv[1];  x3  = pv[2];  x4  = pv[3];
3564cd38560SBarry Smith         pc[0] = m1 = p1*x1 + p3*x2;
3574cd38560SBarry Smith         pc[1] = m2 = p2*x1 + p4*x2;
3584cd38560SBarry Smith         pc[2] = m3 = p1*x3 + p3*x4;
3594cd38560SBarry Smith         pc[3] = m4 = p2*x3 + p4*x4;
3604cd38560SBarry Smith         nz = bi[row+1] - diag_offset[row] - 1;
3614cd38560SBarry Smith         pv += 4;
3624cd38560SBarry Smith         for (j=0; j<nz; j++) {
3634cd38560SBarry Smith           x1   = pv[0];  x2  = pv[1];   x3 = pv[2];  x4  = pv[3];
3644cd38560SBarry Smith           x    = rtmp + 4*pj[j];
3654cd38560SBarry Smith           x[0] -= m1*x1 + m3*x2;
3664cd38560SBarry Smith           x[1] -= m2*x1 + m4*x2;
3674cd38560SBarry Smith           x[2] -= m1*x3 + m3*x4;
3684cd38560SBarry Smith           x[3] -= m2*x3 + m4*x4;
3694cd38560SBarry Smith           pv   += 4;
3704cd38560SBarry Smith         }
371dc0b31edSSatish Balay         ierr = PetscLogFlops(16.0*nz+12.0);CHKERRQ(ierr);
3724cd38560SBarry Smith       }
3734cd38560SBarry Smith       row = *ajtmp++;
3744cd38560SBarry Smith     }
3754cd38560SBarry Smith     /* finished row so stick it into b->a */
3764cd38560SBarry Smith     pv = ba + 4*bi[i];
3774cd38560SBarry Smith     pj = bj + bi[i];
3784cd38560SBarry Smith     nz = bi[i+1] - bi[i];
3794cd38560SBarry Smith     for (j=0; j<nz; j++) {
3804cd38560SBarry Smith       x      = rtmp+4*pj[j];
3814cd38560SBarry Smith       pv[0]  = x[0];  pv[1]  = x[1];  pv[2]  = x[2];  pv[3]  = x[3];
3822efa7f71SHong Zhang       /*
3832efa7f71SHong Zhang       printf(" col %d:",pj[j]);
3842efa7f71SHong Zhang       PetscInt j1;
3852efa7f71SHong Zhang       for (j1=0; j1<4; j1++) printf(" %g,",*(pv+j1));
3862efa7f71SHong Zhang       printf("\n");
3872efa7f71SHong Zhang       */
3884cd38560SBarry Smith       pv   += 4;
3894cd38560SBarry Smith     }
3904cd38560SBarry Smith     /* invert diagonal block */
3914cd38560SBarry Smith     w = ba + 4*diag_offset[i];
3922efa7f71SHong Zhang     /*
3932efa7f71SHong Zhang     printf(" \n%d -th: diag: ",i);
3942efa7f71SHong Zhang     for (j=0; j<4; j++){
3952efa7f71SHong Zhang       printf(" %g,",w[j]);
3962efa7f71SHong Zhang     }
3972efa7f71SHong Zhang     printf("\n----------------------------\n");
3982efa7f71SHong Zhang     */
39962bba022SBarry Smith     ierr = Kernel_A_gets_inverse_A_2(w,shift);CHKERRQ(ierr);
4004cd38560SBarry Smith   }
4014cd38560SBarry Smith 
402606d414cSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
40306e38f1dSHong Zhang   C->ops->solve          = MatSolve_SeqBAIJ_2_NaturalOrdering_inplace;
40406e38f1dSHong Zhang   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_NaturalOrdering_inplace;
4054cd38560SBarry Smith   C->assembled = PETSC_TRUE;
406efee365bSSatish Balay   ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
4074cd38560SBarry Smith   PetscFunctionReturn(0);
4084cd38560SBarry Smith }
4097fc0212eSBarry Smith 
4107fc0212eSBarry Smith /* ----------------------------------------------------------- */
4117fc0212eSBarry Smith /*
4127fc0212eSBarry Smith      Version for when blocks are 1 by 1.
4137fc0212eSBarry Smith */
4144a2ae208SSatish Balay #undef __FUNCT__
41506e38f1dSHong Zhang #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_1_inplace"
41606e38f1dSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqBAIJ_1_inplace(Mat C,Mat A,const MatFactorInfo *info)
4177fc0212eSBarry Smith {
4187fc0212eSBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data;
4197cdf2dbbSSatish Balay   IS             isrow = b->row,isicol = b->icol;
4206849ba73SBarry Smith   PetscErrorCode ierr;
4215d0c19d7SBarry Smith   const PetscInt *r,*ic;
4225d0c19d7SBarry Smith   PetscInt       i,j,n = a->mbs,*bi = b->i,*bj = b->j;
423690b6cddSBarry Smith   PetscInt       *ajtmpold,*ajtmp,nz,row,*ai = a->i,*aj = a->j;
424690b6cddSBarry Smith   PetscInt       *diag_offset = b->diag,diag,*pj;
425329f5518SBarry Smith   MatScalar      *pv,*v,*rtmp,multiplier,*pc;
4263f1db9ecSBarry Smith   MatScalar      *ba = b->a,*aa = a->a;
427f58c8c74SBarry Smith   PetscTruth     row_identity, col_identity;
4287fc0212eSBarry Smith 
4293a40ed3dSBarry Smith   PetscFunctionBegin;
4307fc0212eSBarry Smith   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
4317fc0212eSBarry Smith   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
432b0a32e0cSBarry Smith   ierr  = PetscMalloc((n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
4337fc0212eSBarry Smith 
4347fc0212eSBarry Smith   for (i=0; i<n; i++) {
4354078e994SBarry Smith     nz    = bi[i+1] - bi[i];
4364078e994SBarry Smith     ajtmp = bj + bi[i];
4377fc0212eSBarry Smith     for  (j=0; j<nz; j++) rtmp[ajtmp[j]] = 0.0;
4387fc0212eSBarry Smith 
4397fc0212eSBarry Smith     /* load in initial (unfactored row) */
4404078e994SBarry Smith     nz       = ai[r[i]+1] - ai[r[i]];
4414078e994SBarry Smith     ajtmpold = aj + ai[r[i]];
4424078e994SBarry Smith     v        = aa + ai[r[i]];
4437fc0212eSBarry Smith     for (j=0; j<nz; j++) rtmp[ic[ajtmpold[j]]] =  v[j];
4447fc0212eSBarry Smith 
4457fc0212eSBarry Smith     row = *ajtmp++;
4467fc0212eSBarry Smith     while (row < i) {
4477fc0212eSBarry Smith       pc = rtmp + row;
4487fc0212eSBarry Smith       if (*pc != 0.0) {
4494078e994SBarry Smith         pv         = ba + diag_offset[row];
4504078e994SBarry Smith         pj         = bj + diag_offset[row] + 1;
4517fc0212eSBarry Smith         multiplier = *pc * *pv++;
4527fc0212eSBarry Smith         *pc        = multiplier;
4534078e994SBarry Smith         nz         = bi[row+1] - diag_offset[row] - 1;
4547fc0212eSBarry Smith         for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j];
455dc0b31edSSatish Balay         ierr = PetscLogFlops(1.0+2.0*nz);CHKERRQ(ierr);
4567fc0212eSBarry Smith       }
4577fc0212eSBarry Smith       row = *ajtmp++;
4587fc0212eSBarry Smith     }
4597fc0212eSBarry Smith     /* finished row so stick it into b->a */
4604078e994SBarry Smith     pv = ba + bi[i];
4614078e994SBarry Smith     pj = bj + bi[i];
4624078e994SBarry Smith     nz = bi[i+1] - bi[i];
4637fc0212eSBarry Smith     for (j=0; j<nz; j++) {pv[j] = rtmp[pj[j]];}
4644078e994SBarry Smith     diag = diag_offset[i] - bi[i];
4657fc0212eSBarry Smith     /* check pivot entry for current row */
466a73cf429SBarry Smith     if (pv[diag] == 0.0) {
4673b4a8b6dSBarry Smith       SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot: row in original ordering %D in permuted ordering %D",r[i],i);
4687fc0212eSBarry Smith     }
4697fc0212eSBarry Smith     pv[diag] = 1.0/pv[diag];
4707fc0212eSBarry Smith   }
4717fc0212eSBarry Smith 
472606d414cSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
4737fc0212eSBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
4747fc0212eSBarry Smith   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
475f58c8c74SBarry Smith   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
476f58c8c74SBarry Smith   ierr = ISIdentity(isicol,&col_identity);CHKERRQ(ierr);
477f58c8c74SBarry Smith   if (row_identity && col_identity) {
47806e38f1dSHong Zhang     C->ops->solve          = MatSolve_SeqBAIJ_1_NaturalOrdering_inplace;
47906e38f1dSHong Zhang     C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1_NaturalOrdering_inplace;
480f58c8c74SBarry Smith   } else {
48106e38f1dSHong Zhang     C->ops->solve          = MatSolve_SeqBAIJ_1_inplace;
48206e38f1dSHong Zhang     C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1_inplace;
483f58c8c74SBarry Smith   }
4847fc0212eSBarry Smith   C->assembled = PETSC_TRUE;
485d0f46423SBarry Smith   ierr = PetscLogFlops(C->cmap->n);CHKERRQ(ierr);
4863a40ed3dSBarry Smith   PetscFunctionReturn(0);
4877fc0212eSBarry Smith }
4887fc0212eSBarry Smith 
489e631078cSBarry Smith EXTERN_C_BEGIN
490b24902e0SBarry Smith #undef __FUNCT__
491b24902e0SBarry Smith #define __FUNCT__ "MatGetFactor_seqbaij_petsc"
492b24902e0SBarry Smith PetscErrorCode MatGetFactor_seqbaij_petsc(Mat A,MatFactorType ftype,Mat *B)
493b24902e0SBarry Smith {
494d0f46423SBarry Smith   PetscInt           n = A->rmap->n;
495b24902e0SBarry Smith   PetscErrorCode     ierr;
496b24902e0SBarry Smith 
497b24902e0SBarry Smith   PetscFunctionBegin;
498b24902e0SBarry Smith   ierr = MatCreate(((PetscObject)A)->comm,B);CHKERRQ(ierr);
499b24902e0SBarry Smith   ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr);
500c8342467SHong Zhang   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU || ftype == MAT_FACTOR_ILUDT) {
501b24902e0SBarry Smith     ierr = MatSetType(*B,MATSEQBAIJ);CHKERRQ(ierr);
50206e38f1dSHong Zhang     (*B)->ops->lufactorsymbolic  = MatLUFactorSymbolic_SeqBAIJ_newdatastruct;
50306e38f1dSHong Zhang     (*B)->ops->ilufactorsymbolic = MatILUFactorSymbolic_SeqBAIJ_newdatastruct;
504c8342467SHong Zhang     (*B)->ops->iludtfactor       = MatILUDTFactor_SeqBAIJ;
505b24902e0SBarry Smith   } else if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) {
506b24902e0SBarry Smith     ierr = MatSetType(*B,MATSEQSBAIJ);CHKERRQ(ierr);
5075c9eb25fSBarry Smith     ierr = MatSeqSBAIJSetPreallocation(*B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
5085c9eb25fSBarry Smith     (*B)->ops->iccfactorsymbolic      = MatICCFactorSymbolic_SeqBAIJ;
5095c9eb25fSBarry Smith     (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqBAIJ;
510b24902e0SBarry Smith   } else SETERRQ(PETSC_ERR_SUP,"Factor type not supported");
511719d5645SBarry Smith   (*B)->factor = ftype;
512b24902e0SBarry Smith   PetscFunctionReturn(0);
513b24902e0SBarry Smith }
514e631078cSBarry Smith EXTERN_C_END
515273d9f13SBarry Smith 
516db4efbfdSBarry Smith EXTERN_C_BEGIN
517db4efbfdSBarry Smith #undef __FUNCT__
518db4efbfdSBarry Smith #define __FUNCT__ "MatGetFactorAvailable_seqaij_petsc"
519db4efbfdSBarry Smith PetscErrorCode MatGetFactorAvailable_seqbaij_petsc(Mat A,MatFactorType ftype,PetscTruth *flg)
520db4efbfdSBarry Smith {
521db4efbfdSBarry Smith   PetscFunctionBegin;
522db4efbfdSBarry Smith   *flg = PETSC_TRUE;
523db4efbfdSBarry Smith   PetscFunctionReturn(0);
524db4efbfdSBarry Smith }
525db4efbfdSBarry Smith EXTERN_C_END
526db4efbfdSBarry Smith 
5275d517e7eSBarry Smith /* ----------------------------------------------------------- */
5284a2ae208SSatish Balay #undef __FUNCT__
5294a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqBAIJ"
5300481f469SBarry Smith PetscErrorCode MatLUFactor_SeqBAIJ(Mat A,IS row,IS col,const MatFactorInfo *info)
5315d517e7eSBarry Smith {
532dfbe8321SBarry Smith   PetscErrorCode ierr;
5335d517e7eSBarry Smith   Mat            C;
5345d517e7eSBarry Smith 
5353a40ed3dSBarry Smith   PetscFunctionBegin;
53643244d56SBarry Smith   ierr = MatGetFactor(A,MAT_SOLVER_PETSC,MAT_FACTOR_LU,&C);CHKERRQ(ierr);
537719d5645SBarry Smith   ierr = MatLUFactorSymbolic(C,A,row,col,info);CHKERRQ(ierr);
538719d5645SBarry Smith   ierr = MatLUFactorNumeric(C,A,info);CHKERRQ(ierr);
539db4efbfdSBarry Smith   A->ops->solve            = C->ops->solve;
540db4efbfdSBarry Smith   A->ops->solvetranspose   = C->ops->solvetranspose;
541273d9f13SBarry Smith   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
54252e6d16bSBarry Smith   ierr = PetscLogObjectParent(A,((Mat_SeqBAIJ*)(A->data))->icol);CHKERRQ(ierr);
5433a40ed3dSBarry Smith   PetscFunctionReturn(0);
5445d517e7eSBarry Smith }
5454cd38560SBarry Smith 
5467c4f633dSBarry Smith #include "../src/mat/impls/sbaij/seq/sbaij.h"
547c05c3958SHong Zhang #undef __FUNCT__
548c05c3958SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N"
5490481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N(Mat C,Mat A,const MatFactorInfo *info)
550c05c3958SHong Zhang {
551f3086b4bSHong Zhang   PetscErrorCode ierr;
55278910aadSHong Zhang   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data;
55378910aadSHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
55478910aadSHong Zhang   IS             ip=b->row;
5555d0c19d7SBarry Smith   const PetscInt *rip;
5565d0c19d7SBarry Smith   PetscInt       i,j,mbs=a->mbs,bs=A->rmap->bs,*bi=b->i,*bj=b->j,*bcol;
55778910aadSHong Zhang   PetscInt       *ai=a->i,*aj=a->j;
55878910aadSHong Zhang   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
55978910aadSHong Zhang   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
5603cea4cbeSHong Zhang   PetscReal      zeropivot,rs,shiftnz;
561fbf22428SSatish Balay   PetscReal      shiftpd;
5623cea4cbeSHong Zhang   ChShift_Ctx    sctx;
5633cea4cbeSHong Zhang   PetscInt       newshift;
56478910aadSHong Zhang 
565c05c3958SHong Zhang   PetscFunctionBegin;
5666ad2eaddSHong Zhang   if (bs > 1) {
5676ad2eaddSHong Zhang     if (!a->sbaijMat){
568ceb03754SKris Buschelman       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
5696ad2eaddSHong Zhang     }
570719d5645SBarry Smith     ierr = (a->sbaijMat)->ops->choleskyfactornumeric(C,a->sbaijMat,info);CHKERRQ(ierr);
5716ad2eaddSHong Zhang     ierr = MatDestroy(a->sbaijMat);CHKERRQ(ierr);
5726ad2eaddSHong Zhang     a->sbaijMat = PETSC_NULL;
5736ad2eaddSHong Zhang     PetscFunctionReturn(0);
5746ad2eaddSHong Zhang   }
57578910aadSHong Zhang 
57678910aadSHong Zhang   /* initialization */
5773cea4cbeSHong Zhang   shiftnz   = info->shiftnz;
5783cea4cbeSHong Zhang   shiftpd   = info->shiftpd;
5793cea4cbeSHong Zhang   zeropivot = info->zeropivot;
5803cea4cbeSHong Zhang 
5816ad2eaddSHong Zhang   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
582fca92195SBarry Smith   ierr = PetscMalloc3(mbs,MatScalar,&rtmp,mbs,PetscInt,&il,mbs,PetscInt,&jl);CHKERRQ(ierr);
58378910aadSHong Zhang 
58475567043SBarry Smith   sctx.shift_amount = 0.;
5853cea4cbeSHong Zhang   sctx.nshift       = 0;
58678910aadSHong Zhang   do {
5873cea4cbeSHong Zhang     sctx.chshift = PETSC_FALSE;
58878910aadSHong Zhang     for (i=0; i<mbs; i++) {
589e60cf9a0SBarry Smith       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
59078910aadSHong Zhang     }
59178910aadSHong Zhang 
59278910aadSHong Zhang     for (k = 0; k<mbs; k++){
59378910aadSHong Zhang       bval = ba + bi[k];
59478910aadSHong Zhang       /* initialize k-th row by the perm[k]-th row of A */
59578910aadSHong Zhang       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
59678910aadSHong Zhang       for (j = jmin; j < jmax; j++){
59778910aadSHong Zhang         col = rip[aj[j]];
59878910aadSHong Zhang         if (col >= k){ /* only take upper triangular entry */
59978910aadSHong Zhang           rtmp[col] = aa[j];
60078910aadSHong Zhang           *bval++  = 0.0; /* for in-place factorization */
60178910aadSHong Zhang         }
60278910aadSHong Zhang       }
6033cea4cbeSHong Zhang 
6043cea4cbeSHong Zhang       /* shift the diagonal of the matrix */
6053cea4cbeSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
60678910aadSHong Zhang 
60778910aadSHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
60878910aadSHong Zhang       dk = rtmp[k];
60978910aadSHong Zhang       i = jl[k]; /* first row to be added to k_th row  */
61078910aadSHong Zhang 
61178910aadSHong Zhang       while (i < k){
61278910aadSHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
61378910aadSHong Zhang 
61478910aadSHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
61578910aadSHong Zhang         ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
61678910aadSHong Zhang         uikdi = - ba[ili]*ba[bi[i]];  /* diagonal(k) */
61778910aadSHong Zhang         dk += uikdi*ba[ili];
61878910aadSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
61978910aadSHong Zhang 
62078910aadSHong Zhang         /* add multiple of row i to k-th row */
62178910aadSHong Zhang         jmin = ili + 1; jmax = bi[i+1];
62278910aadSHong Zhang         if (jmin < jmax){
62378910aadSHong Zhang           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
62478910aadSHong Zhang           /* update il and jl for row i */
62578910aadSHong Zhang           il[i] = jmin;
62678910aadSHong Zhang           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
62778910aadSHong Zhang         }
62878910aadSHong Zhang         i = nexti;
62978910aadSHong Zhang       }
63078910aadSHong Zhang 
6313cea4cbeSHong Zhang       /* shift the diagonals when zero pivot is detected */
6323cea4cbeSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
6333cea4cbeSHong Zhang       rs   = 0.0;
63478910aadSHong Zhang       jmin = bi[k]+1;
63578910aadSHong Zhang       nz   = bi[k+1] - jmin;
63678910aadSHong Zhang       if (nz){
63778910aadSHong Zhang         bcol = bj + jmin;
63878910aadSHong Zhang         while (nz--){
63989f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
64089f845c8SHong Zhang           bcol++;
64178910aadSHong Zhang         }
64278910aadSHong Zhang       }
6433cea4cbeSHong Zhang 
6443cea4cbeSHong Zhang       sctx.rs = rs;
6453cea4cbeSHong Zhang       sctx.pv = dk;
64645938b79SHong Zhang       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
64745938b79SHong Zhang       if (newshift == 1) break;
64878910aadSHong Zhang 
64978910aadSHong Zhang       /* copy data into U(k,:) */
65078910aadSHong Zhang       ba[bi[k]] = 1.0/dk; /* U(k,k) */
65178910aadSHong Zhang       jmin = bi[k]+1; jmax = bi[k+1];
65278910aadSHong Zhang       if (jmin < jmax) {
65378910aadSHong Zhang         for (j=jmin; j<jmax; j++){
65478910aadSHong Zhang           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
65578910aadSHong Zhang         }
65678910aadSHong Zhang         /* add the k-th row into il and jl */
65778910aadSHong Zhang         il[k] = jmin;
65878910aadSHong Zhang         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
65978910aadSHong Zhang       }
66078910aadSHong Zhang     }
6613cea4cbeSHong Zhang   } while (sctx.chshift);
662fca92195SBarry Smith   ierr = PetscFree3(rtmp,il,jl);CHKERRQ(ierr);
66378910aadSHong Zhang 
66478910aadSHong Zhang   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
66578910aadSHong Zhang   C->assembled    = PETSC_TRUE;
66678910aadSHong Zhang   C->preallocated = PETSC_TRUE;
667d0f46423SBarry Smith   ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr);
6683cea4cbeSHong Zhang   if (sctx.nshift){
669e738e422SBarry Smith     if (shiftpd) {
6701e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
671e738e422SBarry Smith     } else if (shiftnz) {
672e738e422SBarry Smith       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
67378910aadSHong Zhang     }
67478910aadSHong Zhang   }
675c05c3958SHong Zhang   PetscFunctionReturn(0);
676c05c3958SHong Zhang }
6774cd38560SBarry Smith 
678c05c3958SHong Zhang #undef __FUNCT__
679c05c3958SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering"
6800481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info)
681c05c3958SHong Zhang {
68278910aadSHong Zhang   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data;
68378910aadSHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
68478910aadSHong Zhang   PetscErrorCode ierr;
6853cea4cbeSHong Zhang   PetscInt       i,j,am=a->mbs;
68678910aadSHong Zhang   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
6873cea4cbeSHong Zhang   PetscInt       k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz;
68878910aadSHong Zhang   MatScalar      *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval;
6893cea4cbeSHong Zhang   PetscReal      zeropivot,rs,shiftnz;
690fbf22428SSatish Balay   PetscReal      shiftpd;
6913cea4cbeSHong Zhang   ChShift_Ctx    sctx;
6923cea4cbeSHong Zhang   PetscInt       newshift;
69378910aadSHong Zhang 
694c05c3958SHong Zhang   PetscFunctionBegin;
69578910aadSHong Zhang   /* initialization */
6963cea4cbeSHong Zhang   shiftnz   = info->shiftnz;
6973cea4cbeSHong Zhang   shiftpd   = info->shiftpd;
6983cea4cbeSHong Zhang   zeropivot = info->zeropivot;
6993cea4cbeSHong Zhang 
700fca92195SBarry Smith   ierr = PetscMalloc3(am,MatScalar,&rtmp,am,PetscInt,&il,am,PetscInt,&jl);CHKERRQ(ierr);
70178910aadSHong Zhang 
70275567043SBarry Smith   sctx.shift_amount = 0.;
7033cea4cbeSHong Zhang   sctx.nshift       = 0;
70478910aadSHong Zhang   do {
7053cea4cbeSHong Zhang     sctx.chshift = PETSC_FALSE;
70678910aadSHong Zhang     for (i=0; i<am; i++) {
707e60cf9a0SBarry Smith       rtmp[i] = 0.0; jl[i] = am; il[0] = 0;
70878910aadSHong Zhang     }
70978910aadSHong Zhang 
71078910aadSHong Zhang     for (k = 0; k<am; k++){
71178910aadSHong Zhang     /* initialize k-th row with elements nonzero in row perm(k) of A */
71278910aadSHong Zhang       nz   = ai[k+1] - ai[k];
71378910aadSHong Zhang       acol = aj + ai[k];
71478910aadSHong Zhang       aval = aa + ai[k];
71578910aadSHong Zhang       bval = ba + bi[k];
71678910aadSHong Zhang       while (nz -- ){
71778910aadSHong Zhang         if (*acol < k) { /* skip lower triangular entries */
71878910aadSHong Zhang           acol++; aval++;
71978910aadSHong Zhang         } else {
72078910aadSHong Zhang           rtmp[*acol++] = *aval++;
72178910aadSHong Zhang           *bval++       = 0.0; /* for in-place factorization */
72278910aadSHong Zhang         }
72378910aadSHong Zhang       }
7243cea4cbeSHong Zhang 
7253cea4cbeSHong Zhang       /* shift the diagonal of the matrix */
7263cea4cbeSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
72778910aadSHong Zhang 
72878910aadSHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
72978910aadSHong Zhang       dk = rtmp[k];
73078910aadSHong Zhang       i  = jl[k]; /* first row to be added to k_th row  */
73178910aadSHong Zhang 
73278910aadSHong Zhang       while (i < k){
73378910aadSHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
73478910aadSHong Zhang         /* compute multiplier, update D(k) and U(i,k) */
73578910aadSHong Zhang         ili   = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
73678910aadSHong Zhang         uikdi = - ba[ili]*ba[bi[i]];
73778910aadSHong Zhang         dk   += uikdi*ba[ili];
73878910aadSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
73978910aadSHong Zhang 
74078910aadSHong Zhang         /* add multiple of row i to k-th row ... */
74178910aadSHong Zhang         jmin = ili + 1;
74278910aadSHong Zhang         nz   = bi[i+1] - jmin;
74378910aadSHong Zhang         if (nz > 0){
74478910aadSHong Zhang           bcol = bj + jmin;
74578910aadSHong Zhang           bval = ba + jmin;
74678910aadSHong Zhang           while (nz --) rtmp[*bcol++] += uikdi*(*bval++);
74778910aadSHong Zhang           /* update il and jl for i-th row */
74878910aadSHong Zhang           il[i] = jmin;
74978910aadSHong Zhang           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
75078910aadSHong Zhang         }
75178910aadSHong Zhang         i = nexti;
75278910aadSHong Zhang       }
75378910aadSHong Zhang 
7543cea4cbeSHong Zhang       /* shift the diagonals when zero pivot is detected */
7553cea4cbeSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
7563cea4cbeSHong Zhang       rs   = 0.0;
75778910aadSHong Zhang       jmin = bi[k]+1;
75878910aadSHong Zhang       nz   = bi[k+1] - jmin;
75978910aadSHong Zhang       if (nz){
76078910aadSHong Zhang         bcol = bj + jmin;
76178910aadSHong Zhang         while (nz--){
76289f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
76389f845c8SHong Zhang           bcol++;
76478910aadSHong Zhang         }
76578910aadSHong Zhang       }
7663cea4cbeSHong Zhang 
7673cea4cbeSHong Zhang       sctx.rs = rs;
7683cea4cbeSHong Zhang       sctx.pv = dk;
76945938b79SHong Zhang       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
77045938b79SHong Zhang       if (newshift == 1) break;    /* sctx.shift_amount is updated */
77178910aadSHong Zhang 
77278910aadSHong Zhang       /* copy data into U(k,:) */
77378910aadSHong Zhang       ba[bi[k]] = 1.0/dk;
77478910aadSHong Zhang       jmin      = bi[k]+1;
77578910aadSHong Zhang       nz        = bi[k+1] - jmin;
77678910aadSHong Zhang       if (nz){
77778910aadSHong Zhang         bcol = bj + jmin;
77878910aadSHong Zhang         bval = ba + jmin;
77978910aadSHong Zhang         while (nz--){
78078910aadSHong Zhang           *bval++       = rtmp[*bcol];
78178910aadSHong Zhang           rtmp[*bcol++] = 0.0;
78278910aadSHong Zhang         }
78378910aadSHong Zhang         /* add k-th row into il and jl */
78478910aadSHong Zhang         il[k] = jmin;
78578910aadSHong Zhang         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
78678910aadSHong Zhang       }
78778910aadSHong Zhang     }
7883cea4cbeSHong Zhang   } while (sctx.chshift);
789fca92195SBarry Smith   ierr = PetscFree3(rtmp,il,jl);CHKERRQ(ierr);
79078910aadSHong Zhang 
791*0a3c351aSHong Zhang   C->ops->solve                 = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
792*0a3c351aSHong Zhang   C->ops->solvetranspose        = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
79378910aadSHong Zhang   C->assembled    = PETSC_TRUE;
79478910aadSHong Zhang   C->preallocated = PETSC_TRUE;
795d0f46423SBarry Smith   ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr);
7963cea4cbeSHong Zhang     if (sctx.nshift){
7973cea4cbeSHong Zhang     if (shiftnz) {
7981e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
7993cea4cbeSHong Zhang     } else if (shiftpd) {
8001e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
80178910aadSHong Zhang     }
80278910aadSHong Zhang   }
803c05c3958SHong Zhang   PetscFunctionReturn(0);
804c05c3958SHong Zhang }
805c05c3958SHong Zhang 
806c05c3958SHong Zhang #include "petscbt.h"
8077c4f633dSBarry Smith #include "../src/mat/utils/freespace.h"
808c05c3958SHong Zhang #undef __FUNCT__
809c05c3958SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqBAIJ"
8100481f469SBarry Smith PetscErrorCode MatICCFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
811c05c3958SHong Zhang {
81278910aadSHong Zhang   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data;
81378910aadSHong Zhang   Mat_SeqSBAIJ       *b;
81478910aadSHong Zhang   Mat                B;
81578910aadSHong Zhang   PetscErrorCode     ierr;
81678910aadSHong Zhang   PetscTruth         perm_identity;
8175d0c19d7SBarry Smith   PetscInt           reallocs=0,i,*ai=a->i,*aj=a->j,am=a->mbs,bs=A->rmap->bs,*ui;
8185d0c19d7SBarry Smith   const PetscInt     *rip;
81978910aadSHong Zhang   PetscInt           jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow;
820cfdb8b8aSHong Zhang   PetscInt           nlnk,*lnk,*lnk_lvl=PETSC_NULL,ncols,ncols_upper,*cols,*cols_lvl,*uj,**uj_ptr,**uj_lvl_ptr;
82178910aadSHong Zhang   PetscReal          fill=info->fill,levels=info->levels;
822a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
823a1a86e44SBarry Smith   PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
82478910aadSHong Zhang   PetscBT            lnkbt;
82578910aadSHong Zhang 
826c05c3958SHong Zhang   PetscFunctionBegin;
8276ad2eaddSHong Zhang   if (bs > 1){
8286ad2eaddSHong Zhang     if (!a->sbaijMat){
829ceb03754SKris Buschelman       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
8306ad2eaddSHong Zhang     }
831719d5645SBarry Smith     (fact)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqSBAIJ;  /* undue the change made in MatGetFactor_seqbaij_petsc */
832719d5645SBarry Smith     ierr = MatICCFactorSymbolic(fact,a->sbaijMat,perm,info);CHKERRQ(ierr);
8336ad2eaddSHong Zhang     PetscFunctionReturn(0);
8346ad2eaddSHong Zhang   }
8356ad2eaddSHong Zhang 
83678910aadSHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
83778910aadSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
83878910aadSHong Zhang 
83978910aadSHong Zhang   /* special case that simply copies fill pattern */
84078910aadSHong Zhang   if (!levels && perm_identity) {
84178910aadSHong Zhang     ierr = MatMarkDiagonal_SeqBAIJ(A);CHKERRQ(ierr);
84278910aadSHong Zhang     ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
84378910aadSHong Zhang     for (i=0; i<am; i++) {
84478910aadSHong Zhang       ui[i] = ai[i+1] - a->diag[i]; /* ui: rowlengths - changes when !perm_identity */
84578910aadSHong Zhang     }
846719d5645SBarry Smith     B = fact;
84778910aadSHong Zhang     ierr = MatSeqSBAIJSetPreallocation(B,1,0,ui);CHKERRQ(ierr);
84878910aadSHong Zhang 
849b24902e0SBarry Smith 
85078910aadSHong Zhang     b  = (Mat_SeqSBAIJ*)B->data;
85178910aadSHong Zhang     uj = b->j;
85278910aadSHong Zhang     for (i=0; i<am; i++) {
85378910aadSHong Zhang       aj = a->j + a->diag[i];
85478910aadSHong Zhang       for (j=0; j<ui[i]; j++){
85578910aadSHong Zhang         *uj++ = *aj++;
85678910aadSHong Zhang       }
85778910aadSHong Zhang       b->ilen[i] = ui[i];
85878910aadSHong Zhang     }
85978910aadSHong Zhang     ierr = PetscFree(ui);CHKERRQ(ierr);
860719d5645SBarry Smith     B->factor = MAT_FACTOR_NONE;
86178910aadSHong Zhang     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
86278910aadSHong Zhang     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
863719d5645SBarry Smith     B->factor = MAT_FACTOR_ICC;
86478910aadSHong Zhang 
86578910aadSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
86678910aadSHong Zhang     PetscFunctionReturn(0);
86778910aadSHong Zhang   }
86878910aadSHong Zhang 
86978910aadSHong Zhang   /* initialization */
87078910aadSHong Zhang   ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
871e60cf9a0SBarry Smith   ui[0] = 0;
87278910aadSHong Zhang   ierr  = PetscMalloc((2*am+1)*sizeof(PetscInt),&cols_lvl);CHKERRQ(ierr);
87378910aadSHong Zhang 
87478910aadSHong Zhang   /* jl: linked list for storing indices of the pivot rows
87578910aadSHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
876fca92195SBarry Smith   ierr = PetscMalloc4(am,PetscInt*,&uj_ptr,am,PetscInt*,&uj_lvl_ptr,am,PetscInt,&il,am,PetscInt,&jl);CHKERRQ(ierr);
87778910aadSHong Zhang   for (i=0; i<am; i++){
878e60cf9a0SBarry Smith     jl[i] = am; il[i] = 0;
87978910aadSHong Zhang   }
88078910aadSHong Zhang 
88178910aadSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
88278910aadSHong Zhang   nlnk = am + 1;
88378910aadSHong Zhang   ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
88478910aadSHong Zhang 
88578910aadSHong Zhang   /* initial FreeSpace size is fill*(ai[am]+1) */
886a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
88778910aadSHong Zhang   current_space = free_space;
888a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr);
88978910aadSHong Zhang   current_space_lvl = free_space_lvl;
89078910aadSHong Zhang 
89178910aadSHong Zhang   for (k=0; k<am; k++){  /* for each active row k */
89278910aadSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
89378910aadSHong Zhang     nzk   = 0;
89478910aadSHong Zhang     ncols = ai[rip[k]+1] - ai[rip[k]];
89578910aadSHong Zhang     ncols_upper = 0;
89678910aadSHong Zhang     cols        = cols_lvl + am;
89778910aadSHong Zhang     for (j=0; j<ncols; j++){
89878910aadSHong Zhang       i = rip[*(aj + ai[rip[k]] + j)];
89978910aadSHong Zhang       if (i >= k){ /* only take upper triangular entry */
90078910aadSHong Zhang         cols[ncols_upper] = i;
90178910aadSHong Zhang         cols_lvl[ncols_upper] = -1;  /* initialize level for nonzero entries */
90278910aadSHong Zhang         ncols_upper++;
90378910aadSHong Zhang       }
90478910aadSHong Zhang     }
90578910aadSHong Zhang     ierr = PetscIncompleteLLAdd(ncols_upper,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
90678910aadSHong Zhang     nzk += nlnk;
90778910aadSHong Zhang 
90878910aadSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
90978910aadSHong Zhang     prow = jl[k]; /* 1st pivot row */
91078910aadSHong Zhang 
91178910aadSHong Zhang     while (prow < k){
91278910aadSHong Zhang       nextprow = jl[prow];
91378910aadSHong Zhang 
91478910aadSHong Zhang       /* merge prow into k-th row */
91578910aadSHong Zhang       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
91678910aadSHong Zhang       jmax = ui[prow+1];
91778910aadSHong Zhang       ncols = jmax-jmin;
91878910aadSHong Zhang       i     = jmin - ui[prow];
91978910aadSHong Zhang       cols = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
92078910aadSHong Zhang       for (j=0; j<ncols; j++) cols_lvl[j] = *(uj_lvl_ptr[prow] + i + j);
9215a8e39fbSHong Zhang       ierr = PetscIncompleteLLAddSorted(ncols,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
92278910aadSHong Zhang       nzk += nlnk;
92378910aadSHong Zhang 
92478910aadSHong Zhang       /* update il and jl for prow */
92578910aadSHong Zhang       if (jmin < jmax){
92678910aadSHong Zhang         il[prow] = jmin;
92778910aadSHong Zhang         j = *cols; jl[prow] = jl[j]; jl[j] = prow;
92878910aadSHong Zhang       }
92978910aadSHong Zhang       prow = nextprow;
93078910aadSHong Zhang     }
93178910aadSHong Zhang 
93278910aadSHong Zhang     /* if free space is not available, make more free space */
93378910aadSHong Zhang     if (current_space->local_remaining<nzk) {
93478910aadSHong Zhang       i = am - k + 1; /* num of unfactored rows */
93578910aadSHong Zhang       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
936a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
937a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(i,&current_space_lvl);CHKERRQ(ierr);
93878910aadSHong Zhang       reallocs++;
93978910aadSHong Zhang     }
94078910aadSHong Zhang 
94178910aadSHong Zhang     /* copy data into free_space and free_space_lvl, then initialize lnk */
94278910aadSHong Zhang     ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
94378910aadSHong Zhang 
94478910aadSHong Zhang     /* add the k-th row into il and jl */
94578910aadSHong Zhang     if (nzk-1 > 0){
94678910aadSHong Zhang       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
94778910aadSHong Zhang       jl[k] = jl[i]; jl[i] = k;
94878910aadSHong Zhang       il[k] = ui[k] + 1;
94978910aadSHong Zhang     }
95078910aadSHong Zhang     uj_ptr[k]     = current_space->array;
95178910aadSHong Zhang     uj_lvl_ptr[k] = current_space_lvl->array;
95278910aadSHong Zhang 
95378910aadSHong Zhang     current_space->array           += nzk;
95478910aadSHong Zhang     current_space->local_used      += nzk;
95578910aadSHong Zhang     current_space->local_remaining -= nzk;
95678910aadSHong Zhang 
95778910aadSHong Zhang     current_space_lvl->array           += nzk;
95878910aadSHong Zhang     current_space_lvl->local_used      += nzk;
95978910aadSHong Zhang     current_space_lvl->local_remaining -= nzk;
96078910aadSHong Zhang 
96178910aadSHong Zhang     ui[k+1] = ui[k] + nzk;
96278910aadSHong Zhang   }
96378910aadSHong Zhang 
9646cf91177SBarry Smith #if defined(PETSC_USE_INFO)
96578910aadSHong Zhang   if (ai[am] != 0) {
96678910aadSHong Zhang     PetscReal af = ((PetscReal)(2*ui[am]-am))/((PetscReal)ai[am]);
967ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
968ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
969ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
97078910aadSHong Zhang   } else {
971ae15b995SBarry Smith     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
97278910aadSHong Zhang   }
97363ba0a88SBarry Smith #endif
97478910aadSHong Zhang 
97578910aadSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
976fca92195SBarry Smith   ierr = PetscFree4(uj_ptr,uj_lvl_ptr,il,jl);CHKERRQ(ierr);
97778910aadSHong Zhang   ierr = PetscFree(cols_lvl);CHKERRQ(ierr);
97878910aadSHong Zhang 
97978910aadSHong Zhang   /* destroy list of free space and other temporary array(s) */
98078910aadSHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
981a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
98278910aadSHong Zhang   ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
983a1a86e44SBarry Smith   ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
98478910aadSHong Zhang 
98578910aadSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
986719d5645SBarry Smith   B = fact;
987ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation(B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
98878910aadSHong Zhang 
98978910aadSHong Zhang   b = (Mat_SeqSBAIJ*)B->data;
99078910aadSHong Zhang   b->singlemalloc = PETSC_FALSE;
991e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
992e6b907acSBarry Smith   b->free_ij       = PETSC_TRUE;
99378910aadSHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
99478910aadSHong Zhang   b->j    = uj;
99578910aadSHong Zhang   b->i    = ui;
99678910aadSHong Zhang   b->diag = 0;
99778910aadSHong Zhang   b->ilen = 0;
99878910aadSHong Zhang   b->imax = 0;
99978910aadSHong Zhang   b->row  = perm;
100078910aadSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
100178910aadSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
100278910aadSHong Zhang   b->icol = perm;
100378910aadSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
100478910aadSHong Zhang   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
100552e6d16bSBarry Smith   ierr    = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
100678910aadSHong Zhang   b->maxnz = b->nz = ui[am];
100778910aadSHong Zhang 
100878910aadSHong Zhang   B->info.factor_mallocs    = reallocs;
100978910aadSHong Zhang   B->info.fill_ratio_given  = fill;
101075567043SBarry Smith   if (ai[am] != 0.) {
101178910aadSHong Zhang     B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
101278910aadSHong Zhang   } else {
101378910aadSHong Zhang     B->info.fill_ratio_needed = 0.0;
101478910aadSHong Zhang   }
101578910aadSHong Zhang   if (perm_identity){
1016*0a3c351aSHong Zhang     B->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
1017*0a3c351aSHong Zhang     B->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
101878910aadSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
101978910aadSHong Zhang   } else {
1020719d5645SBarry Smith     (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
102178910aadSHong Zhang   }
1022c05c3958SHong Zhang   PetscFunctionReturn(0);
1023c05c3958SHong Zhang }
1024c05c3958SHong Zhang 
1025c05c3958SHong Zhang #undef __FUNCT__
1026c05c3958SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqBAIJ"
10270481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
1028c05c3958SHong Zhang {
102978910aadSHong Zhang   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data;
103078910aadSHong Zhang   Mat_SeqSBAIJ       *b;
103178910aadSHong Zhang   Mat                B;
103278910aadSHong Zhang   PetscErrorCode     ierr;
103378910aadSHong Zhang   PetscTruth         perm_identity;
103478910aadSHong Zhang   PetscReal          fill = info->fill;
10355d0c19d7SBarry Smith   const PetscInt     *rip;
10365d0c19d7SBarry Smith   PetscInt           i,mbs=a->mbs,bs=A->rmap->bs,*ai=a->i,*aj=a->j,reallocs=0,prow;
103778910aadSHong Zhang   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
103878910aadSHong Zhang   PetscInt           nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr;
1039a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
104078910aadSHong Zhang   PetscBT            lnkbt;
104178910aadSHong Zhang 
1042c05c3958SHong Zhang   PetscFunctionBegin;
10436ad2eaddSHong Zhang   if (bs > 1) { /* convert to seqsbaij */
10446ad2eaddSHong Zhang     if (!a->sbaijMat){
1045ceb03754SKris Buschelman       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
10466ad2eaddSHong Zhang     }
1047719d5645SBarry Smith     (fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqSBAIJ; /* undue the change made in MatGetFactor_seqbaij_petsc */
1048719d5645SBarry Smith     ierr = MatCholeskyFactorSymbolic(fact,a->sbaijMat,perm,info);CHKERRQ(ierr);
10496ad2eaddSHong Zhang     PetscFunctionReturn(0);
10506ad2eaddSHong Zhang   }
10516ad2eaddSHong Zhang 
105278910aadSHong Zhang   /* check whether perm is the identity mapping */
105378910aadSHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
1054c84f5b01SHong Zhang   if (!perm_identity) SETERRQ(PETSC_ERR_SUP,"Matrix reordering is not supported");
105578910aadSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
105678910aadSHong Zhang 
105778910aadSHong Zhang   /* initialization */
105878910aadSHong Zhang   ierr  = PetscMalloc((mbs+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
1059e60cf9a0SBarry Smith   ui[0] = 0;
106078910aadSHong Zhang 
106178910aadSHong Zhang   /* jl: linked list for storing indices of the pivot rows
106278910aadSHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
1063fca92195SBarry Smith   ierr = PetscMalloc4(mbs,PetscInt*,&ui_ptr,mbs,PetscInt,&il,mbs,PetscInt,&jl,mbs,PetscInt,&cols);CHKERRQ(ierr);
106478910aadSHong Zhang   for (i=0; i<mbs; i++){
1065e60cf9a0SBarry Smith     jl[i] = mbs; il[i] = 0;
106678910aadSHong Zhang   }
106778910aadSHong Zhang 
106878910aadSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
106978910aadSHong Zhang   nlnk = mbs + 1;
107078910aadSHong Zhang   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
107178910aadSHong Zhang 
107278910aadSHong Zhang   /* initial FreeSpace size is fill*(ai[mbs]+1) */
1073a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[mbs]+1)),&free_space);CHKERRQ(ierr);
107478910aadSHong Zhang   current_space = free_space;
107578910aadSHong Zhang 
107678910aadSHong Zhang   for (k=0; k<mbs; k++){  /* for each active row k */
107778910aadSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
107878910aadSHong Zhang     nzk   = 0;
107978910aadSHong Zhang     ncols = ai[rip[k]+1] - ai[rip[k]];
1080e60cf9a0SBarry Smith     ncols_upper = 0;
108178910aadSHong Zhang     for (j=0; j<ncols; j++){
108278910aadSHong Zhang       i = rip[*(aj + ai[rip[k]] + j)];
108378910aadSHong Zhang       if (i >= k){ /* only take upper triangular entry */
108478910aadSHong Zhang         cols[ncols_upper] = i;
108578910aadSHong Zhang         ncols_upper++;
108678910aadSHong Zhang       }
108778910aadSHong Zhang     }
108878910aadSHong Zhang     ierr = PetscLLAdd(ncols_upper,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
108978910aadSHong Zhang     nzk += nlnk;
109078910aadSHong Zhang 
109178910aadSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
109278910aadSHong Zhang     prow = jl[k]; /* 1st pivot row */
109378910aadSHong Zhang 
109478910aadSHong Zhang     while (prow < k){
109578910aadSHong Zhang       nextprow = jl[prow];
109678910aadSHong Zhang       /* merge prow into k-th row */
109778910aadSHong Zhang       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
109878910aadSHong Zhang       jmax = ui[prow+1];
109978910aadSHong Zhang       ncols = jmax-jmin;
110078910aadSHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
11015a8e39fbSHong Zhang       ierr = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
110278910aadSHong Zhang       nzk += nlnk;
110378910aadSHong Zhang 
110478910aadSHong Zhang       /* update il and jl for prow */
110578910aadSHong Zhang       if (jmin < jmax){
110678910aadSHong Zhang         il[prow] = jmin;
110778910aadSHong Zhang         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
110878910aadSHong Zhang       }
110978910aadSHong Zhang       prow = nextprow;
111078910aadSHong Zhang     }
111178910aadSHong Zhang 
111278910aadSHong Zhang     /* if free space is not available, make more free space */
111378910aadSHong Zhang     if (current_space->local_remaining<nzk) {
111478910aadSHong Zhang       i = mbs - k + 1; /* num of unfactored rows */
111578910aadSHong Zhang       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
1116a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
111778910aadSHong Zhang       reallocs++;
111878910aadSHong Zhang     }
111978910aadSHong Zhang 
112078910aadSHong Zhang     /* copy data into free space, then initialize lnk */
112178910aadSHong Zhang     ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
112278910aadSHong Zhang 
112378910aadSHong Zhang     /* add the k-th row into il and jl */
112478910aadSHong Zhang     if (nzk-1 > 0){
112578910aadSHong Zhang       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
112678910aadSHong Zhang       jl[k] = jl[i]; jl[i] = k;
112778910aadSHong Zhang       il[k] = ui[k] + 1;
112878910aadSHong Zhang     }
112978910aadSHong Zhang     ui_ptr[k] = current_space->array;
113078910aadSHong Zhang     current_space->array           += nzk;
113178910aadSHong Zhang     current_space->local_used      += nzk;
113278910aadSHong Zhang     current_space->local_remaining -= nzk;
113378910aadSHong Zhang 
113478910aadSHong Zhang     ui[k+1] = ui[k] + nzk;
113578910aadSHong Zhang   }
113678910aadSHong Zhang 
11376cf91177SBarry Smith #if defined(PETSC_USE_INFO)
113875567043SBarry Smith   if (ai[mbs] != 0.) {
113978910aadSHong Zhang     PetscReal af = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
1140ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
1141ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
1142ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
114378910aadSHong Zhang   } else {
1144ae15b995SBarry Smith     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
114578910aadSHong Zhang   }
114663ba0a88SBarry Smith #endif
114778910aadSHong Zhang 
114878910aadSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
1149fca92195SBarry Smith   ierr = PetscFree4(ui_ptr,il,jl,cols);CHKERRQ(ierr);
115078910aadSHong Zhang 
115178910aadSHong Zhang   /* destroy list of free space and other temporary array(s) */
115278910aadSHong Zhang   ierr = PetscMalloc((ui[mbs]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
1153a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
115478910aadSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
115578910aadSHong Zhang 
115678910aadSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
1157719d5645SBarry Smith   B    = fact;
1158ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
115978910aadSHong Zhang 
116078910aadSHong Zhang   b = (Mat_SeqSBAIJ*)B->data;
116178910aadSHong Zhang   b->singlemalloc = PETSC_FALSE;
1162e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
1163e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
116478910aadSHong Zhang   ierr = PetscMalloc((ui[mbs]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
116578910aadSHong Zhang   b->j    = uj;
116678910aadSHong Zhang   b->i    = ui;
116778910aadSHong Zhang   b->diag = 0;
116878910aadSHong Zhang   b->ilen = 0;
116978910aadSHong Zhang   b->imax = 0;
117078910aadSHong Zhang   b->row  = perm;
117178910aadSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
117278910aadSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
117378910aadSHong Zhang   b->icol = perm;
117478910aadSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
117578910aadSHong Zhang   ierr    = PetscMalloc((mbs+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
117652e6d16bSBarry Smith   ierr    = PetscLogObjectMemory(B,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
117778910aadSHong Zhang   b->maxnz = b->nz = ui[mbs];
117878910aadSHong Zhang 
117978910aadSHong Zhang   B->info.factor_mallocs    = reallocs;
118078910aadSHong Zhang   B->info.fill_ratio_given  = fill;
118175567043SBarry Smith   if (ai[mbs] != 0.) {
118278910aadSHong Zhang     B->info.fill_ratio_needed = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
118378910aadSHong Zhang   } else {
118478910aadSHong Zhang     B->info.fill_ratio_needed = 0.0;
118578910aadSHong Zhang   }
118678910aadSHong Zhang   if (perm_identity){
11876ad2eaddSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
118878910aadSHong Zhang   } else {
11896ad2eaddSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
119078910aadSHong Zhang   }
1191c05c3958SHong Zhang   PetscFunctionReturn(0);
1192c05c3958SHong Zhang }
1193c8342467SHong Zhang 
11942efa7f71SHong Zhang #undef __FUNCT__
1195a2d6a19aSShri Abhyankar #define __FUNCT__ "MatSolve_SeqBAIJ_N_NaturalOrdering_newdatastruct"
1196a2d6a19aSShri Abhyankar PetscErrorCode MatSolve_SeqBAIJ_N_NaturalOrdering_newdatastruct(Mat A,Vec bb,Vec xx)
11971a83e813SShri Abhyankar {
11981a83e813SShri Abhyankar   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ *)A->data;
11991a83e813SShri Abhyankar   PetscErrorCode ierr;
12001a83e813SShri Abhyankar   const PetscInt *ai=a->i,*aj=a->j,*adiag=a->diag,*vi;
12011a83e813SShri Abhyankar   PetscInt       i,k,n=a->mbs;
12021a83e813SShri Abhyankar   PetscInt       nz,bs=A->rmap->bs,bs2=a->bs2;
12031a83e813SShri Abhyankar   MatScalar      *aa=a->a,*v;
12041a83e813SShri Abhyankar   PetscScalar    *x,*b,*s,*t,*ls;
12051a83e813SShri Abhyankar 
12061a83e813SShri Abhyankar   PetscFunctionBegin;
12071a83e813SShri Abhyankar   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
12081a83e813SShri Abhyankar   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
12091a83e813SShri Abhyankar   t  = a->solve_work;
12101a83e813SShri Abhyankar 
12111a83e813SShri Abhyankar   /* forward solve the lower triangular */
12121a83e813SShri Abhyankar   ierr = PetscMemcpy(t,b,bs*sizeof(PetscScalar));CHKERRQ(ierr); /* copy 1st block of b to t */
12131a83e813SShri Abhyankar 
12141a83e813SShri Abhyankar   for (i=1; i<n; i++) {
12151a83e813SShri Abhyankar     v   = aa + bs2*ai[i];
12161a83e813SShri Abhyankar     vi  = aj + ai[i];
12171a83e813SShri Abhyankar     nz = ai[i+1] - ai[i];
12181a83e813SShri Abhyankar     s = t + bs*i;
12191a83e813SShri Abhyankar     ierr = PetscMemcpy(s,b+bs*i,bs*sizeof(PetscScalar));CHKERRQ(ierr); /* copy i_th block of b to t */
12201a83e813SShri Abhyankar     for(k=0;k<nz;k++){
12211a83e813SShri Abhyankar       Kernel_v_gets_v_minus_A_times_w(bs,s,v,t+bs*vi[k]);
12221a83e813SShri Abhyankar       v += bs2;
12231a83e813SShri Abhyankar     }
12241a83e813SShri Abhyankar   }
12251a83e813SShri Abhyankar 
12261a83e813SShri Abhyankar   /* backward solve the upper triangular */
12271a83e813SShri Abhyankar   ls = a->solve_work + A->cmap->n;
12281a83e813SShri Abhyankar   for (i=n-1; i>=0; i--){
12291a83e813SShri Abhyankar     v  = aa + bs2*(adiag[i+1]+1);
12301a83e813SShri Abhyankar     vi = aj + adiag[i+1]+1;
12311a83e813SShri Abhyankar     nz = adiag[i] - adiag[i+1]-1;
12321a83e813SShri Abhyankar     ierr = PetscMemcpy(ls,t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr);
12331a83e813SShri Abhyankar     for(k=0;k<nz;k++){
12341a83e813SShri Abhyankar       Kernel_v_gets_v_minus_A_times_w(bs,ls,v,t+bs*vi[k]);
12351a83e813SShri Abhyankar       v += bs2;
12361a83e813SShri Abhyankar     }
12371a83e813SShri Abhyankar     Kernel_w_gets_A_times_v(bs,ls,aa+bs2*adiag[i],t+i*bs); /* *inv(diagonal[i]) */
12381a83e813SShri Abhyankar     ierr = PetscMemcpy(x+i*bs,t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr);
12391a83e813SShri Abhyankar   }
12401a83e813SShri Abhyankar 
12411a83e813SShri Abhyankar   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
12421a83e813SShri Abhyankar   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
12431a83e813SShri Abhyankar   ierr = PetscLogFlops(2.0*(a->bs2)*(a->nz) - A->rmap->bs*A->cmap->n);CHKERRQ(ierr);
12441a83e813SShri Abhyankar   PetscFunctionReturn(0);
12451a83e813SShri Abhyankar }
12461a83e813SShri Abhyankar 
12472efa7f71SHong Zhang #undef __FUNCT__
1248a2d6a19aSShri Abhyankar #define __FUNCT__ "MatSolve_SeqBAIJ_N_newdatastruct"
1249a2d6a19aSShri Abhyankar PetscErrorCode MatSolve_SeqBAIJ_N_newdatastruct(Mat A,Vec bb,Vec xx)
125035aa4fcfSShri Abhyankar {
125135aa4fcfSShri Abhyankar   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ *)A->data;
125235aa4fcfSShri Abhyankar   IS             iscol=a->col,isrow=a->row;
125335aa4fcfSShri Abhyankar   PetscErrorCode ierr;
125435aa4fcfSShri Abhyankar   const PetscInt *r,*c,*rout,*cout,*ai=a->i,*aj=a->j,*adiag=a->diag,*vi;
125535aa4fcfSShri Abhyankar   PetscInt       i,m,n=a->mbs;
125635aa4fcfSShri Abhyankar   PetscInt       nz,bs=A->rmap->bs,bs2=a->bs2;
125735aa4fcfSShri Abhyankar   MatScalar      *aa=a->a,*v;
125835aa4fcfSShri Abhyankar   PetscScalar    *x,*b,*s,*t,*ls;
125935aa4fcfSShri Abhyankar 
126035aa4fcfSShri Abhyankar   PetscFunctionBegin;
126135aa4fcfSShri Abhyankar   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
126235aa4fcfSShri Abhyankar   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
126335aa4fcfSShri Abhyankar   t  = a->solve_work;
126435aa4fcfSShri Abhyankar 
126535aa4fcfSShri Abhyankar   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
126635aa4fcfSShri Abhyankar   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
126735aa4fcfSShri Abhyankar 
126835aa4fcfSShri Abhyankar   /* forward solve the lower triangular */
126935aa4fcfSShri Abhyankar   ierr = PetscMemcpy(t,b+bs*r[0],bs*sizeof(PetscScalar));CHKERRQ(ierr);
127035aa4fcfSShri Abhyankar   for (i=1; i<n; i++) {
127135aa4fcfSShri Abhyankar     v   = aa + bs2*ai[i];
127235aa4fcfSShri Abhyankar     vi  = aj + ai[i];
127335aa4fcfSShri Abhyankar     nz = ai[i+1] - ai[i];
127435aa4fcfSShri Abhyankar     s = t + bs*i;
127535aa4fcfSShri Abhyankar     ierr = PetscMemcpy(s,b+bs*r[i],bs*sizeof(PetscScalar));CHKERRQ(ierr);
127635aa4fcfSShri Abhyankar     for(m=0;m<nz;m++){
127735aa4fcfSShri Abhyankar       Kernel_v_gets_v_minus_A_times_w(bs,s,v,t+bs*vi[m]);
127835aa4fcfSShri Abhyankar       v += bs2;
127935aa4fcfSShri Abhyankar     }
128035aa4fcfSShri Abhyankar   }
128135aa4fcfSShri Abhyankar 
128235aa4fcfSShri Abhyankar   /* backward solve the upper triangular */
128335aa4fcfSShri Abhyankar   ls = a->solve_work + A->cmap->n;
128435aa4fcfSShri Abhyankar   for (i=n-1; i>=0; i--){
128535aa4fcfSShri Abhyankar     v  = aa + bs2*(adiag[i+1]+1);
128635aa4fcfSShri Abhyankar     vi = aj + adiag[i+1]+1;
128735aa4fcfSShri Abhyankar     nz = adiag[i] - adiag[i+1] - 1;
128835aa4fcfSShri Abhyankar     ierr = PetscMemcpy(ls,t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr);
128935aa4fcfSShri Abhyankar     for(m=0;m<nz;m++){
129035aa4fcfSShri Abhyankar       Kernel_v_gets_v_minus_A_times_w(bs,ls,v,t+bs*vi[m]);
129135aa4fcfSShri Abhyankar       v += bs2;
129235aa4fcfSShri Abhyankar     }
129335aa4fcfSShri Abhyankar     Kernel_w_gets_A_times_v(bs,ls,v,t+i*bs); /* *inv(diagonal[i]) */
129435aa4fcfSShri Abhyankar     ierr = PetscMemcpy(x + bs*c[i],t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr);
129535aa4fcfSShri Abhyankar   }
129635aa4fcfSShri Abhyankar   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
129735aa4fcfSShri Abhyankar   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
129835aa4fcfSShri Abhyankar   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
129935aa4fcfSShri Abhyankar   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
130035aa4fcfSShri Abhyankar   ierr = PetscLogFlops(2.0*(a->bs2)*(a->nz) - A->rmap->bs*A->cmap->n);CHKERRQ(ierr);
130135aa4fcfSShri Abhyankar   PetscFunctionReturn(0);
130235aa4fcfSShri Abhyankar }
130335aa4fcfSShri Abhyankar 
130435aa4fcfSShri Abhyankar #undef __FUNCT__
13052efa7f71SHong Zhang #define __FUNCT__ "BlockAbs_privat"
130691d91c02SMatthew Knepley PetscErrorCode BlockAbs_private(PetscInt nbs,PetscInt bs2,PetscScalar *blockarray,PetscReal *absarray)
13072efa7f71SHong Zhang {
13082efa7f71SHong Zhang   PetscErrorCode     ierr;
13092efa7f71SHong Zhang   PetscInt           i,j;
13102efa7f71SHong Zhang   PetscFunctionBegin;
13112efa7f71SHong Zhang   ierr = PetscMemzero(absarray,(nbs+1)*sizeof(PetscReal));CHKERRQ(ierr);
13122efa7f71SHong Zhang   for (i=0; i<nbs; i++){
13132efa7f71SHong Zhang     for (j=0; j<bs2; j++){
13142efa7f71SHong Zhang       if (absarray[i] < PetscAbsScalar(blockarray[i*nbs+j])) absarray[i] = PetscAbsScalar(blockarray[i*nbs+j]);
13152efa7f71SHong Zhang     }
13162efa7f71SHong Zhang   }
13172efa7f71SHong Zhang   PetscFunctionReturn(0);
13182efa7f71SHong Zhang }
13192efa7f71SHong Zhang 
1320c8342467SHong Zhang #undef __FUNCT__
1321c8342467SHong Zhang #define __FUNCT__ "MatILUDTFactor_SeqBAIJ"
1322c8342467SHong Zhang PetscErrorCode MatILUDTFactor_SeqBAIJ(Mat A,IS isrow,IS iscol,const MatFactorInfo *info,Mat *fact)
1323c8342467SHong Zhang {
13242efa7f71SHong Zhang   Mat                B = *fact;
13252efa7f71SHong Zhang   Mat_SeqBAIJ        *a=(Mat_SeqBAIJ*)A->data,*b;
13262efa7f71SHong Zhang   IS                 isicol;
13272efa7f71SHong Zhang   PetscErrorCode     ierr;
13282efa7f71SHong Zhang   const PetscInt     *r,*ic;
13292efa7f71SHong Zhang   PetscInt           i,mbs=a->mbs,bs=A->rmap->bs,bs2=a->bs2,*ai=a->i,*aj=a->j,*ajtmp,*adiag;
13302efa7f71SHong Zhang   PetscInt           *bi,*bj,*bdiag;
13312efa7f71SHong Zhang 
13322efa7f71SHong Zhang   PetscInt           row,nzi,nzi_bl,nzi_bu,*im,dtcount,nzi_al,nzi_au;
13332efa7f71SHong Zhang   PetscInt           nlnk,*lnk;
13342efa7f71SHong Zhang   PetscBT            lnkbt;
13352efa7f71SHong Zhang   PetscTruth         row_identity,icol_identity,both_identity;
13362efa7f71SHong Zhang   MatScalar          *aatmp,*pv,*batmp,*ba,*rtmp,*pc,*multiplier,*vtmp;
13372efa7f71SHong Zhang   const PetscInt     *ics;
13382efa7f71SHong Zhang   PetscInt           j,nz,*pj,*bjtmp,k,ncut,*jtmp;
13392efa7f71SHong Zhang 
13406da515a1SHong Zhang   PetscReal          dt=info->dt; /* shift=info->shiftinblocks; */
13412efa7f71SHong Zhang   PetscInt           nnz_max;
13422efa7f71SHong Zhang   PetscTruth         missing;
1343021822bcSHong Zhang   PetscReal          *vtmp_abs;
134497ef94ebSSatish Balay   MatScalar          *v_work;
134597ef94ebSSatish Balay   PetscInt           *v_pivots;
13462efa7f71SHong Zhang 
1347c8342467SHong Zhang   PetscFunctionBegin;
13482efa7f71SHong Zhang   /* ------- symbolic factorization, can be reused ---------*/
13492efa7f71SHong Zhang   ierr = MatMissingDiagonal(A,&missing,&i);CHKERRQ(ierr);
13502efa7f71SHong Zhang   if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i);
13512efa7f71SHong Zhang   adiag=a->diag;
13522efa7f71SHong Zhang 
13532efa7f71SHong Zhang   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
13542efa7f71SHong Zhang 
13552efa7f71SHong Zhang   /* bdiag is location of diagonal in factor */
13562efa7f71SHong Zhang   ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
13572efa7f71SHong Zhang 
13582efa7f71SHong Zhang   /* allocate row pointers bi */
13596bce7ff8SHong Zhang   ierr = PetscMalloc((2*mbs+2)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
13602efa7f71SHong Zhang 
13612efa7f71SHong Zhang   /* allocate bj and ba; max num of nonzero entries is (ai[n]+2*n*dtcount+2) */
13622efa7f71SHong Zhang   dtcount = (PetscInt)info->dtcount;
13636bce7ff8SHong Zhang   if (dtcount > mbs-1) dtcount = mbs-1;
13646bce7ff8SHong Zhang   nnz_max  = ai[mbs]+2*mbs*dtcount +2;
13656da515a1SHong Zhang   /* printf("MatILUDTFactor_SeqBAIJ, bs %d, ai[mbs] %d, nnz_max  %d, dtcount %d\n",bs,ai[mbs],nnz_max,dtcount); */
13662efa7f71SHong Zhang   ierr = PetscMalloc(nnz_max*sizeof(PetscInt),&bj);CHKERRQ(ierr);
13676bce7ff8SHong Zhang   nnz_max = nnz_max*bs2;
13682efa7f71SHong Zhang   ierr = PetscMalloc(nnz_max*sizeof(MatScalar),&ba);CHKERRQ(ierr);
13692efa7f71SHong Zhang 
13702efa7f71SHong Zhang   /* put together the new matrix */
13712efa7f71SHong Zhang   ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
13722efa7f71SHong Zhang   ierr = PetscLogObjectParent(B,isicol);CHKERRQ(ierr);
13732efa7f71SHong Zhang   b    = (Mat_SeqBAIJ*)(B)->data;
13742efa7f71SHong Zhang   b->free_a       = PETSC_TRUE;
13752efa7f71SHong Zhang   b->free_ij      = PETSC_TRUE;
13762efa7f71SHong Zhang   b->singlemalloc = PETSC_FALSE;
13772efa7f71SHong Zhang   b->a          = ba;
13782efa7f71SHong Zhang   b->j          = bj;
13792efa7f71SHong Zhang   b->i          = bi;
13802efa7f71SHong Zhang   b->diag       = bdiag;
13812efa7f71SHong Zhang   b->ilen       = 0;
13822efa7f71SHong Zhang   b->imax       = 0;
13832efa7f71SHong Zhang   b->row        = isrow;
13842efa7f71SHong Zhang   b->col        = iscol;
13852efa7f71SHong Zhang   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
13862efa7f71SHong Zhang   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
13872efa7f71SHong Zhang   b->icol       = isicol;
13882efa7f71SHong Zhang   ierr = PetscMalloc((bs*(mbs+1))*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
13892efa7f71SHong Zhang 
13902efa7f71SHong Zhang   ierr = PetscLogObjectMemory(B,nnz_max*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
13912efa7f71SHong Zhang   b->maxnz = nnz_max/bs2;
13922efa7f71SHong Zhang 
13932efa7f71SHong Zhang   (B)->factor                = MAT_FACTOR_ILUDT;
13942efa7f71SHong Zhang   (B)->info.factor_mallocs   = 0;
13952efa7f71SHong Zhang   (B)->info.fill_ratio_given = ((PetscReal)nnz_max)/((PetscReal)(ai[mbs]*bs2));
13962efa7f71SHong Zhang   CHKMEMQ;
13972efa7f71SHong Zhang   /* ------- end of symbolic factorization ---------*/
13982efa7f71SHong Zhang   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
13992efa7f71SHong Zhang   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
14002efa7f71SHong Zhang   ics  = ic;
14012efa7f71SHong Zhang 
14022efa7f71SHong Zhang   /* linked list for storing column indices of the active row */
14032efa7f71SHong Zhang   nlnk = mbs + 1;
14042efa7f71SHong Zhang   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
14052efa7f71SHong Zhang 
14062efa7f71SHong Zhang   /* im: used by PetscLLAddSortedLU(); jtmp: working array for column indices of active row */
1407fca92195SBarry Smith   ierr = PetscMalloc2(mbs,PetscInt,&im,mbs,PetscInt,&jtmp);CHKERRQ(ierr);
14082efa7f71SHong Zhang   /* rtmp, vtmp: working arrays for sparse and contiguous row entries of active row */
1409fca92195SBarry Smith   ierr = PetscMalloc2(mbs*bs2,MatScalar,&rtmp,mbs*bs2,MatScalar,&vtmp);CHKERRQ(ierr);
1410021822bcSHong Zhang   ierr = PetscMalloc((mbs+1)*sizeof(PetscReal),&vtmp_abs);CHKERRQ(ierr);
1411fca92195SBarry Smith   ierr = PetscMalloc3(bs,MatScalar,&v_work,bs2,MatScalar,&multiplier,bs,PetscInt,&v_pivots);CHKERRQ(ierr);
14122efa7f71SHong Zhang 
1413e60cf9a0SBarry Smith   bi[0]    = 0;
14142efa7f71SHong Zhang   bdiag[0] = (nnz_max/bs2)-1; /* location of diagonal in factor B */
14156bce7ff8SHong Zhang   bi[2*mbs+1] = bdiag[0]+1; /* endof bj and ba array */
14162efa7f71SHong Zhang   for (i=0; i<mbs; i++) {
14172efa7f71SHong Zhang     /* copy initial fill into linked list */
14182efa7f71SHong Zhang     nzi = 0; /* nonzeros for active row i */
14192efa7f71SHong Zhang     nzi = ai[r[i]+1] - ai[r[i]];
14202efa7f71SHong 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);
14212efa7f71SHong Zhang     nzi_al = adiag[r[i]] - ai[r[i]];
14222efa7f71SHong Zhang     nzi_au = ai[r[i]+1] - adiag[r[i]] -1;
14236da515a1SHong Zhang     /* printf("row %d, nzi_al/au %d %d\n",i,nzi_al,nzi_au); */
14242efa7f71SHong Zhang 
14252efa7f71SHong Zhang     /* load in initial unfactored row */
14262efa7f71SHong Zhang     ajtmp = aj + ai[r[i]];
14272efa7f71SHong Zhang     ierr = PetscLLAddPerm(nzi,ajtmp,ic,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
1428fca92195SBarry Smith     ierr = PetscMemzero(rtmp,mbs*bs2*sizeof(PetscScalar));CHKERRQ(ierr);
14292efa7f71SHong Zhang     aatmp = a->a + bs2*ai[r[i]];
14302efa7f71SHong Zhang     for (j=0; j<nzi; j++) {
14312efa7f71SHong Zhang       ierr = PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],aatmp+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr);
14322efa7f71SHong Zhang     }
14332efa7f71SHong Zhang 
14342efa7f71SHong Zhang     /* add pivot rows into linked list */
14352efa7f71SHong Zhang     row = lnk[mbs];
14362efa7f71SHong Zhang     while (row < i) {
14372efa7f71SHong Zhang       nzi_bl = bi[row+1] - bi[row] + 1;
14382efa7f71SHong Zhang       bjtmp = bj + bdiag[row+1]+1; /* points to 1st column next to the diagonal in U */
14392efa7f71SHong Zhang       ierr  = PetscLLAddSortedLU(bjtmp,row,nlnk,lnk,lnkbt,i,nzi_bl,im);CHKERRQ(ierr);
14402efa7f71SHong Zhang       nzi  += nlnk;
14412efa7f71SHong Zhang       row   = lnk[row];
14422efa7f71SHong Zhang     }
14432efa7f71SHong Zhang 
14442efa7f71SHong Zhang     /* copy data from lnk into jtmp, then initialize lnk */
14452efa7f71SHong Zhang     ierr = PetscLLClean(mbs,mbs,nzi,lnk,jtmp,lnkbt);CHKERRQ(ierr);
14462efa7f71SHong Zhang 
14472efa7f71SHong Zhang     /* numerical factorization */
14482efa7f71SHong Zhang     bjtmp = jtmp;
14492efa7f71SHong Zhang     row   = *bjtmp++; /* 1st pivot row */
14502efa7f71SHong Zhang 
14512efa7f71SHong Zhang     while  (row < i) {
14522efa7f71SHong Zhang       pc = rtmp + bs2*row;
14532efa7f71SHong Zhang       pv = ba + bs2*bdiag[row]; /* inv(diag) of the pivot row */
14542efa7f71SHong Zhang       Kernel_A_gets_A_times_B(bs,pc,pv,multiplier); /* pc= multiplier = pc*inv(diag[row]) */
14552efa7f71SHong Zhang       ierr = BlockAbs_private(1,bs2,pc,vtmp_abs);CHKERRQ(ierr);
14562efa7f71SHong Zhang       if (vtmp_abs[0] > dt){ /* apply tolerance dropping rule */
14572efa7f71SHong Zhang         pj         = bj + bdiag[row+1] + 1; /* point to 1st entry of U(row,:) */
14582efa7f71SHong Zhang         pv         = ba + bs2*(bdiag[row+1] + 1);
14592efa7f71SHong Zhang         nz         = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:), excluding diagonal */
14602efa7f71SHong Zhang         for (j=0; j<nz; j++){
14612efa7f71SHong Zhang           Kernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j);
14622efa7f71SHong Zhang         }
14636da515a1SHong Zhang         /* ierr = PetscLogFlops(bslog*(nz+1.0)-bs);CHKERRQ(ierr); */
14642efa7f71SHong Zhang       }
14652efa7f71SHong Zhang       row = *bjtmp++;
14662efa7f71SHong Zhang     }
14672efa7f71SHong Zhang 
14682efa7f71SHong Zhang     /* copy sparse rtmp into contiguous vtmp; separate L and U part */
14692efa7f71SHong Zhang     nzi_bl = 0; j = 0;
14702efa7f71SHong Zhang     while (jtmp[j] < i){ /* L-part. Note: jtmp is sorted */
14712efa7f71SHong Zhang       ierr = PetscMemcpy(vtmp+bs2*j,rtmp+bs2*jtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
14722efa7f71SHong Zhang       nzi_bl++; j++;
14732efa7f71SHong Zhang     }
14742efa7f71SHong Zhang     nzi_bu = nzi - nzi_bl -1;
14756da515a1SHong Zhang     /* printf("nzi %d, nzi_bl %d, nzi_bu %d\n",nzi,nzi_bl,nzi_bu); */
14762efa7f71SHong Zhang 
14772efa7f71SHong Zhang     while (j < nzi){ /* U-part */
14782efa7f71SHong Zhang       ierr = PetscMemcpy(vtmp+bs2*j,rtmp+bs2*jtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
14792efa7f71SHong Zhang       /*
14802efa7f71SHong Zhang       printf(" col %d: ",jtmp[j]);
14812efa7f71SHong Zhang       for (j1=0; j1<bs2; j1++) printf(" %g",*(vtmp+bs2*j+j1));
14822efa7f71SHong Zhang       printf(" \n");
14832efa7f71SHong Zhang       */
14842efa7f71SHong Zhang       j++;
14852efa7f71SHong Zhang     }
14862efa7f71SHong Zhang 
14872efa7f71SHong Zhang     ierr = BlockAbs_private(nzi,bs2,vtmp,vtmp_abs);CHKERRQ(ierr);
14882efa7f71SHong Zhang     /*
14892efa7f71SHong Zhang     printf(" row %d, nzi %d, vtmp_abs\n",i,nzi);
14902efa7f71SHong Zhang     for (j1=0; j1<nzi; j1++) printf(" (%d %g),",jtmp[j1],vtmp_abs[j1]);
14912efa7f71SHong Zhang     printf(" \n");
14922efa7f71SHong Zhang     */
14932efa7f71SHong Zhang     bjtmp = bj + bi[i];
14942efa7f71SHong Zhang     batmp = ba + bs2*bi[i];
14952efa7f71SHong Zhang     /* apply level dropping rule to L part */
14962efa7f71SHong Zhang     ncut = nzi_al + dtcount;
14972efa7f71SHong Zhang     if (ncut < nzi_bl){
1498021822bcSHong Zhang       ierr = PetscSortSplitReal(ncut,nzi_bl,vtmp_abs,jtmp);CHKERRQ(ierr);
14992efa7f71SHong Zhang       ierr = PetscSortIntWithScalarArray(ncut,jtmp,vtmp);CHKERRQ(ierr);
15002efa7f71SHong Zhang     } else {
15012efa7f71SHong Zhang       ncut = nzi_bl;
15022efa7f71SHong Zhang     }
15032efa7f71SHong Zhang     for (j=0; j<ncut; j++){
15042efa7f71SHong Zhang       bjtmp[j] = jtmp[j];
15052efa7f71SHong Zhang       ierr = PetscMemcpy(batmp+bs2*j,rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
15062efa7f71SHong Zhang       /*
15072efa7f71SHong Zhang       printf(" col %d: ",bjtmp[j]);
15082efa7f71SHong Zhang       for (j1=0; j1<bs2; j1++) printf(" %g,",*(batmp+bs2*j+j1));
15092efa7f71SHong Zhang       printf("\n");
15102efa7f71SHong Zhang       */
15112efa7f71SHong Zhang     }
15122efa7f71SHong Zhang     bi[i+1] = bi[i] + ncut;
15132efa7f71SHong Zhang     nzi = ncut + 1;
15142efa7f71SHong Zhang 
15152efa7f71SHong Zhang     /* apply level dropping rule to U part */
15162efa7f71SHong Zhang     ncut = nzi_au + dtcount;
15172efa7f71SHong Zhang     if (ncut < nzi_bu){
1518021822bcSHong Zhang       ierr = PetscSortSplitReal(ncut,nzi_bu,vtmp_abs+nzi_bl+1,jtmp+nzi_bl+1);CHKERRQ(ierr);
15192efa7f71SHong Zhang       ierr = PetscSortIntWithScalarArray(ncut,jtmp+nzi_bl+1,vtmp+nzi_bl+1);CHKERRQ(ierr);
15202efa7f71SHong Zhang     } else {
15212efa7f71SHong Zhang       ncut = nzi_bu;
15222efa7f71SHong Zhang     }
15232efa7f71SHong Zhang     nzi += ncut;
15242efa7f71SHong Zhang 
15252efa7f71SHong Zhang     /* mark bdiagonal */
15262efa7f71SHong Zhang     bdiag[i+1]    = bdiag[i] - (ncut + 1);
15276bce7ff8SHong Zhang     bi[2*mbs - i] = bi[2*mbs - i +1] - (ncut + 1);
15286bce7ff8SHong Zhang 
15292efa7f71SHong Zhang     bjtmp = bj + bdiag[i];
15302efa7f71SHong Zhang     batmp = ba + bs2*bdiag[i];
15312efa7f71SHong Zhang     ierr = PetscMemcpy(batmp,rtmp+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
15322efa7f71SHong Zhang     *bjtmp = i;
15332efa7f71SHong Zhang     /*
15342efa7f71SHong Zhang     printf(" diag %d: ",*bjtmp);
15352efa7f71SHong Zhang     for (j=0; j<bs2; j++){
15362efa7f71SHong Zhang       printf(" %g,",batmp[j]);
15372efa7f71SHong Zhang     }
15382efa7f71SHong Zhang     printf("\n");
15392efa7f71SHong Zhang     */
15402efa7f71SHong Zhang     bjtmp = bj + bdiag[i+1]+1;
15412efa7f71SHong Zhang     batmp = ba + (bdiag[i+1]+1)*bs2;
15422efa7f71SHong Zhang 
15432efa7f71SHong Zhang     for (k=0; k<ncut; k++){
15442efa7f71SHong Zhang       bjtmp[k] = jtmp[nzi_bl+1+k];
15452efa7f71SHong Zhang       ierr = PetscMemcpy(batmp+bs2*k,rtmp+bs2*bjtmp[k],bs2*sizeof(MatScalar));CHKERRQ(ierr);
15462efa7f71SHong Zhang       /*
15472efa7f71SHong Zhang       printf(" col %d:",bjtmp[k]);
15482efa7f71SHong Zhang       for (j1=0; j1<bs2; j1++) printf(" %g,",*(batmp+bs2*k+j1));
15492efa7f71SHong Zhang       printf("\n");
15502efa7f71SHong Zhang       */
15512efa7f71SHong Zhang     }
15522efa7f71SHong Zhang 
15532efa7f71SHong Zhang     im[i] = nzi; /* used by PetscLLAddSortedLU() */
15542efa7f71SHong Zhang 
15552efa7f71SHong Zhang     /* invert diagonal block for simplier triangular solves - add shift??? */
15562efa7f71SHong Zhang     batmp = ba + bs2*bdiag[i];
15572efa7f71SHong Zhang     ierr = Kernel_A_gets_inverse_A(bs,batmp,v_pivots,v_work);CHKERRQ(ierr);
15582efa7f71SHong Zhang   } /* for (i=0; i<mbs; i++) */
1559fca92195SBarry Smith   ierr = PetscFree3(v_work,multiplier,v_pivots);CHKERRQ(ierr);
15602efa7f71SHong Zhang 
15616da515a1SHong Zhang   /* printf("end of L %d, beginning of U %d\n",bi[mbs],bdiag[mbs]); */
15622efa7f71SHong 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]);
15632efa7f71SHong Zhang 
15642efa7f71SHong Zhang   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
15652efa7f71SHong Zhang   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
15662efa7f71SHong Zhang 
15672efa7f71SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
15682efa7f71SHong Zhang 
1569fca92195SBarry Smith   ierr = PetscFree2(im,jtmp);CHKERRQ(ierr);
1570fca92195SBarry Smith   ierr = PetscFree2(rtmp,vtmp);CHKERRQ(ierr);
15712efa7f71SHong Zhang 
15722efa7f71SHong Zhang   ierr = PetscLogFlops(bs2*B->cmap->n);CHKERRQ(ierr);
15732efa7f71SHong Zhang   b->maxnz = b->nz = bi[mbs] + bdiag[0] - bdiag[mbs];
15742efa7f71SHong Zhang 
15752efa7f71SHong Zhang   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
15762efa7f71SHong Zhang   ierr = ISIdentity(isicol,&icol_identity);CHKERRQ(ierr);
15772efa7f71SHong Zhang   both_identity = (PetscTruth) (row_identity && icol_identity);
15782efa7f71SHong Zhang   if (row_identity && icol_identity) {
157984a281e5SHong Zhang     B->ops->solve = MatSolve_SeqBAIJ_N_NaturalOrdering_newdatastruct;
15802efa7f71SHong Zhang   } else {
158184a281e5SHong Zhang     B->ops->solve = MatSolve_SeqBAIJ_N_newdatastruct;
15822efa7f71SHong Zhang   }
15832efa7f71SHong Zhang 
15842efa7f71SHong Zhang   B->ops->solveadd          = 0;
15852efa7f71SHong Zhang   B->ops->solvetranspose    = 0;
15862efa7f71SHong Zhang   B->ops->solvetransposeadd = 0;
15872efa7f71SHong Zhang   B->ops->matsolve          = 0;
15882efa7f71SHong Zhang   B->assembled              = PETSC_TRUE;
15892efa7f71SHong Zhang   B->preallocated           = PETSC_TRUE;
1590c8342467SHong Zhang   PetscFunctionReturn(0);
1591c8342467SHong Zhang }
1592