xref: /petsc/src/mat/impls/baij/seq/baijfact.c (revision 2e92ee13a8395f820cc1e3fd74a7607ed52efa2a)
1be1d678aSKris Buschelman 
25d517e7eSBarry Smith /*
3ec3a8b7bSBarry Smith     Factorization code for BAIJ format.
45d517e7eSBarry Smith */
5c6db04a5SJed Brown #include <../src/mat/impls/baij/seq/baij.h>
6af0996ceSBarry Smith #include <petsc/private/kernels/blockinvert.h>
7ec3a8b7bSBarry Smith 
8b588c5a2SHong Zhang #undef __FUNCT__
94dd39f65SShri Abhyankar #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2"
104dd39f65SShri Abhyankar PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2(Mat B,Mat A,const MatFactorInfo *info)
11b588c5a2SHong Zhang {
12b588c5a2SHong Zhang   Mat            C     =B;
13b588c5a2SHong Zhang   Mat_SeqBAIJ    *a    =(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data;
14b588c5a2SHong Zhang   IS             isrow = b->row,isicol = b->icol;
15b588c5a2SHong Zhang   PetscErrorCode ierr;
165a586d82SBarry Smith   const PetscInt *r,*ic;
17bbd65245SShri Abhyankar   PetscInt       i,j,k,nz,nzL,row,*pj;
18bbd65245SShri Abhyankar   const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,bs2=a->bs2;
19bbd65245SShri Abhyankar   const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag;
20bbd65245SShri Abhyankar   MatScalar      *rtmp,*pc,*mwork,*pv;
21bbd65245SShri Abhyankar   MatScalar      *aa=a->a,*v;
22bbd65245SShri Abhyankar   PetscInt       flg;
23182b8fbaSHong Zhang   PetscReal      shift = info->shiftamount;
24*2e92ee13SHong Zhang   PetscBool      zeropivotdetected;
25b588c5a2SHong Zhang 
26b588c5a2SHong Zhang   PetscFunctionBegin;
27b588c5a2SHong Zhang   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
28b588c5a2SHong Zhang   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
29b588c5a2SHong Zhang 
304c000e2eSHong Zhang   /* generate work space needed by the factorization */
31dcca6d9dSJed Brown   ierr = PetscMalloc2(bs2*n,&rtmp,bs2,&mwork);CHKERRQ(ierr);
324c000e2eSHong Zhang   ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr);
33b588c5a2SHong Zhang 
34b588c5a2SHong Zhang   for (i=0; i<n; i++) {
35b588c5a2SHong Zhang     /* zero rtmp */
36b588c5a2SHong Zhang     /* L part */
37b588c5a2SHong Zhang     nz    = bi[i+1] - bi[i];
38b588c5a2SHong Zhang     bjtmp = bj + bi[i];
39b588c5a2SHong Zhang     for  (j=0; j<nz; j++) {
40b588c5a2SHong Zhang       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
41b588c5a2SHong Zhang     }
42b588c5a2SHong Zhang 
43b588c5a2SHong Zhang     /* U part */
440c4413a7SShri Abhyankar     nz    = bdiag[i] - bdiag[i+1];
450c4413a7SShri Abhyankar     bjtmp = bj + bdiag[i+1]+1;
460c4413a7SShri Abhyankar     for  (j=0; j<nz; j++) {
470c4413a7SShri Abhyankar       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
480c4413a7SShri Abhyankar     }
490c4413a7SShri Abhyankar 
500c4413a7SShri Abhyankar     /* load in initial (unfactored row) */
510c4413a7SShri Abhyankar     nz    = ai[r[i]+1] - ai[r[i]];
520c4413a7SShri Abhyankar     ajtmp = aj + ai[r[i]];
530c4413a7SShri Abhyankar     v     = aa + bs2*ai[r[i]];
540c4413a7SShri Abhyankar     for (j=0; j<nz; j++) {
550c4413a7SShri Abhyankar       ierr = PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr);
560c4413a7SShri Abhyankar     }
570c4413a7SShri Abhyankar 
580c4413a7SShri Abhyankar     /* elimination */
590c4413a7SShri Abhyankar     bjtmp = bj + bi[i];
600c4413a7SShri Abhyankar     nzL   = bi[i+1] - bi[i];
610c4413a7SShri Abhyankar     for (k=0; k < nzL; k++) {
620c4413a7SShri Abhyankar       row = bjtmp[k];
630c4413a7SShri Abhyankar       pc  = rtmp + bs2*row;
64c35f09e5SBarry Smith       for (flg=0,j=0; j<bs2; j++) {
65c35f09e5SBarry Smith         if (pc[j] != (PetscScalar)0.0) {
66c35f09e5SBarry Smith           flg = 1;
67c35f09e5SBarry Smith           break;
68c35f09e5SBarry Smith         }
69c35f09e5SBarry Smith       }
700c4413a7SShri Abhyankar       if (flg) {
710c4413a7SShri Abhyankar         pv = b->a + bs2*bdiag[row];
7296b95a6bSBarry Smith         /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
7396b95a6bSBarry Smith         ierr = PetscKernel_A_gets_A_times_B_2(pc,pv,mwork);CHKERRQ(ierr);
740c4413a7SShri Abhyankar 
750c4413a7SShri Abhyankar         pj = b->j + bdiag[row+1]+1; /* begining of U(row,:) */
760c4413a7SShri Abhyankar         pv = b->a + bs2*(bdiag[row+1]+1);
770c4413a7SShri Abhyankar         nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries inU(row,:), excluding diag */
780c4413a7SShri Abhyankar         for (j=0; j<nz; j++) {
7996b95a6bSBarry Smith           /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
800c4413a7SShri Abhyankar           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
810c4413a7SShri Abhyankar           v    = rtmp + 4*pj[j];
8296b95a6bSBarry Smith           ierr = PetscKernel_A_gets_A_minus_B_times_C_2(v,pc,pv);CHKERRQ(ierr);
830c4413a7SShri Abhyankar           pv  += 4;
840c4413a7SShri Abhyankar         }
850c4413a7SShri Abhyankar         ierr = PetscLogFlops(16*nz+12);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
860c4413a7SShri Abhyankar       }
870c4413a7SShri Abhyankar     }
880c4413a7SShri Abhyankar 
890c4413a7SShri Abhyankar     /* finished row so stick it into b->a */
900c4413a7SShri Abhyankar     /* L part */
910c4413a7SShri Abhyankar     pv = b->a + bs2*bi[i];
920c4413a7SShri Abhyankar     pj = b->j + bi[i];
930c4413a7SShri Abhyankar     nz = bi[i+1] - bi[i];
940c4413a7SShri Abhyankar     for (j=0; j<nz; j++) {
950c4413a7SShri Abhyankar       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
960c4413a7SShri Abhyankar     }
970c4413a7SShri Abhyankar 
980c4413a7SShri Abhyankar     /* Mark diagonal and invert diagonal for simplier triangular solves */
990c4413a7SShri Abhyankar     pv   = b->a + bs2*bdiag[i];
1000c4413a7SShri Abhyankar     pj   = b->j + bdiag[i];
1010c4413a7SShri Abhyankar     ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr);
10296b95a6bSBarry Smith     /* ierr = PetscKernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); */
103*2e92ee13SHong Zhang     ierr = PetscKernel_A_gets_inverse_A_2(pv,shift,!A->erroriffailure,&zeropivotdetected);CHKERRQ(ierr);
104*2e92ee13SHong Zhang     if (zeropivotdetected) B->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1050c4413a7SShri Abhyankar 
1060c4413a7SShri Abhyankar     /* U part */
1070c4413a7SShri Abhyankar     pv = b->a + bs2*(bdiag[i+1]+1);
1080c4413a7SShri Abhyankar     pj = b->j + bdiag[i+1]+1;
1090c4413a7SShri Abhyankar     nz = bdiag[i] - bdiag[i+1] - 1;
1100c4413a7SShri Abhyankar     for (j=0; j<nz; j++) {
1110c4413a7SShri Abhyankar       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
1120c4413a7SShri Abhyankar     }
1130c4413a7SShri Abhyankar   }
1140c4413a7SShri Abhyankar 
115fca92195SBarry Smith   ierr = PetscFree2(rtmp,mwork);CHKERRQ(ierr);
1160c4413a7SShri Abhyankar   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
1170c4413a7SShri Abhyankar   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
11826fbe8dcSKarl Rupp 
1194dd39f65SShri Abhyankar   C->ops->solve          = MatSolve_SeqBAIJ_2;
1204dd39f65SShri Abhyankar   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2;
1210c4413a7SShri Abhyankar   C->assembled           = PETSC_TRUE;
12226fbe8dcSKarl Rupp 
123766f9fbaSBarry Smith   ierr = PetscLogFlops(1.333333333333*2*2*2*n);CHKERRQ(ierr); /* from inverting diagonal blocks */
1240c4413a7SShri Abhyankar   PetscFunctionReturn(0);
1250c4413a7SShri Abhyankar }
1260c4413a7SShri Abhyankar 
1274c000e2eSHong Zhang #undef __FUNCT__
1284dd39f65SShri Abhyankar #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering"
1294dd39f65SShri Abhyankar PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info)
1304c000e2eSHong Zhang {
1314c000e2eSHong Zhang   Mat            C =B;
1324c000e2eSHong Zhang   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data;
1334c000e2eSHong Zhang   PetscErrorCode ierr;
134bbd65245SShri Abhyankar   PetscInt       i,j,k,nz,nzL,row,*pj;
135bbd65245SShri Abhyankar   const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,bs2=a->bs2;
136bbd65245SShri Abhyankar   const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag;
137bbd65245SShri Abhyankar   MatScalar      *rtmp,*pc,*mwork,*pv;
138bbd65245SShri Abhyankar   MatScalar      *aa=a->a,*v;
139bbd65245SShri Abhyankar   PetscInt       flg;
140182b8fbaSHong Zhang   PetscReal      shift = info->shiftamount;
141*2e92ee13SHong Zhang   PetscBool      zeropivotdetected;
1424c000e2eSHong Zhang 
1434c000e2eSHong Zhang   PetscFunctionBegin;
1444c000e2eSHong Zhang   /* generate work space needed by the factorization */
145dcca6d9dSJed Brown   ierr = PetscMalloc2(bs2*n,&rtmp,bs2,&mwork);CHKERRQ(ierr);
1464c000e2eSHong Zhang   ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr);
1474c000e2eSHong Zhang 
1484c000e2eSHong Zhang   for (i=0; i<n; i++) {
1494c000e2eSHong Zhang     /* zero rtmp */
1504c000e2eSHong Zhang     /* L part */
1514c000e2eSHong Zhang     nz    = bi[i+1] - bi[i];
1524c000e2eSHong Zhang     bjtmp = bj + bi[i];
1534c000e2eSHong Zhang     for  (j=0; j<nz; j++) {
1544c000e2eSHong Zhang       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
1554c000e2eSHong Zhang     }
1564c000e2eSHong Zhang 
1574c000e2eSHong Zhang     /* U part */
158b2b2dd24SShri Abhyankar     nz    = bdiag[i] - bdiag[i+1];
159b2b2dd24SShri Abhyankar     bjtmp = bj + bdiag[i+1]+1;
160b2b2dd24SShri Abhyankar     for  (j=0; j<nz; j++) {
161b2b2dd24SShri Abhyankar       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
162b2b2dd24SShri Abhyankar     }
163b2b2dd24SShri Abhyankar 
164b2b2dd24SShri Abhyankar     /* load in initial (unfactored row) */
165b2b2dd24SShri Abhyankar     nz    = ai[i+1] - ai[i];
166b2b2dd24SShri Abhyankar     ajtmp = aj + ai[i];
167b2b2dd24SShri Abhyankar     v     = aa + bs2*ai[i];
168b2b2dd24SShri Abhyankar     for (j=0; j<nz; j++) {
169b2b2dd24SShri Abhyankar       ierr = PetscMemcpy(rtmp+bs2*ajtmp[j],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr);
170b2b2dd24SShri Abhyankar     }
171b2b2dd24SShri Abhyankar 
172b2b2dd24SShri Abhyankar     /* elimination */
173b2b2dd24SShri Abhyankar     bjtmp = bj + bi[i];
174b2b2dd24SShri Abhyankar     nzL   = bi[i+1] - bi[i];
175b2b2dd24SShri Abhyankar     for (k=0; k < nzL; k++) {
176b2b2dd24SShri Abhyankar       row = bjtmp[k];
177b2b2dd24SShri Abhyankar       pc  = rtmp + bs2*row;
178c35f09e5SBarry Smith       for (flg=0,j=0; j<bs2; j++) {
179c35f09e5SBarry Smith         if (pc[j]!=(PetscScalar)0.0) {
180c35f09e5SBarry Smith           flg = 1;
181c35f09e5SBarry Smith           break;
182c35f09e5SBarry Smith         }
183c35f09e5SBarry Smith       }
184b2b2dd24SShri Abhyankar       if (flg) {
185b2b2dd24SShri Abhyankar         pv = b->a + bs2*bdiag[row];
18696b95a6bSBarry Smith         /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
18796b95a6bSBarry Smith         ierr = PetscKernel_A_gets_A_times_B_2(pc,pv,mwork);CHKERRQ(ierr);
188b2b2dd24SShri Abhyankar 
189b2b2dd24SShri Abhyankar         pj = b->j + bdiag[row+1]+1; /* beginning of U(row,:) */
190b2b2dd24SShri Abhyankar         pv = b->a + bs2*(bdiag[row+1]+1);
191b2b2dd24SShri Abhyankar         nz = bdiag[row]-bdiag[row+1] - 1; /* num of entries in U(row,:) excluding diag */
192b2b2dd24SShri Abhyankar         for (j=0; j<nz; j++) {
19396b95a6bSBarry Smith           /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
194b2b2dd24SShri Abhyankar           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
195b2b2dd24SShri Abhyankar           v    = rtmp + 4*pj[j];
19696b95a6bSBarry Smith           ierr = PetscKernel_A_gets_A_minus_B_times_C_2(v,pc,pv);CHKERRQ(ierr);
197b2b2dd24SShri Abhyankar           pv  += 4;
198b2b2dd24SShri Abhyankar         }
199b2b2dd24SShri Abhyankar         ierr = PetscLogFlops(16*nz+12);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
200b2b2dd24SShri Abhyankar       }
201b2b2dd24SShri Abhyankar     }
202b2b2dd24SShri Abhyankar 
203b2b2dd24SShri Abhyankar     /* finished row so stick it into b->a */
204b2b2dd24SShri Abhyankar     /* L part */
205b2b2dd24SShri Abhyankar     pv = b->a + bs2*bi[i];
206b2b2dd24SShri Abhyankar     pj = b->j + bi[i];
207b2b2dd24SShri Abhyankar     nz = bi[i+1] - bi[i];
208b2b2dd24SShri Abhyankar     for (j=0; j<nz; j++) {
209b2b2dd24SShri Abhyankar       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
210b2b2dd24SShri Abhyankar     }
211b2b2dd24SShri Abhyankar 
212b2b2dd24SShri Abhyankar     /* Mark diagonal and invert diagonal for simplier triangular solves */
213b2b2dd24SShri Abhyankar     pv   = b->a + bs2*bdiag[i];
214b2b2dd24SShri Abhyankar     pj   = b->j + bdiag[i];
215b2b2dd24SShri Abhyankar     ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr);
21696b95a6bSBarry Smith     /* ierr = PetscKernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); */
217*2e92ee13SHong Zhang     ierr = PetscKernel_A_gets_inverse_A_2(pv,shift,!A->erroriffailure,&zeropivotdetected);CHKERRQ(ierr);
218*2e92ee13SHong Zhang     if (zeropivotdetected) B->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
219b2b2dd24SShri Abhyankar 
220b2b2dd24SShri Abhyankar     /* U part */
221b2b2dd24SShri Abhyankar     /*
222b2b2dd24SShri Abhyankar     pv = b->a + bs2*bi[2*n-i];
223b2b2dd24SShri Abhyankar     pj = b->j + bi[2*n-i];
224b2b2dd24SShri Abhyankar     nz = bi[2*n-i+1] - bi[2*n-i] - 1;
225b2b2dd24SShri Abhyankar     */
226b2b2dd24SShri Abhyankar     pv = b->a + bs2*(bdiag[i+1]+1);
227b2b2dd24SShri Abhyankar     pj = b->j + bdiag[i+1]+1;
228b2b2dd24SShri Abhyankar     nz = bdiag[i] - bdiag[i+1] - 1;
229b2b2dd24SShri Abhyankar     for (j=0; j<nz; j++) {
230b2b2dd24SShri Abhyankar       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
231b2b2dd24SShri Abhyankar     }
232b2b2dd24SShri Abhyankar   }
233fca92195SBarry Smith   ierr = PetscFree2(rtmp,mwork);CHKERRQ(ierr);
234ae3d28f0SHong Zhang 
2354dd39f65SShri Abhyankar   C->ops->solve          = MatSolve_SeqBAIJ_2_NaturalOrdering;
2369f5c0bcdSBarry Smith   C->ops->forwardsolve   = MatForwardSolve_SeqBAIJ_2_NaturalOrdering;
2379f5c0bcdSBarry Smith   C->ops->backwardsolve  = MatBackwardSolve_SeqBAIJ_2_NaturalOrdering;
2384dd39f65SShri Abhyankar   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_NaturalOrdering;
239b2b2dd24SShri Abhyankar   C->assembled           = PETSC_TRUE;
24026fbe8dcSKarl Rupp 
241766f9fbaSBarry Smith   ierr = PetscLogFlops(1.333333333333*2*2*2*n);CHKERRQ(ierr); /* from inverting diagonal blocks */
242b2b2dd24SShri Abhyankar   PetscFunctionReturn(0);
243b2b2dd24SShri Abhyankar }
244b2b2dd24SShri Abhyankar 
245b2b2dd24SShri Abhyankar #undef __FUNCT__
24606e38f1dSHong Zhang #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2_inplace"
24706e38f1dSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_inplace(Mat B,Mat A,const MatFactorInfo *info)
2484eeb42bcSBarry Smith {
249719d5645SBarry Smith   Mat            C     = B;
2504eeb42bcSBarry Smith   Mat_SeqBAIJ    *a    = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)C->data;
2517cdf2dbbSSatish Balay   IS             isrow = b->row,isicol = b->icol;
2526849ba73SBarry Smith   PetscErrorCode ierr;
2535d0c19d7SBarry Smith   const PetscInt *r,*ic;
2545d0c19d7SBarry Smith   PetscInt       i,j,n = a->mbs,*bi = b->i,*bj = b->j;
255690b6cddSBarry Smith   PetscInt       *ajtmpold,*ajtmp,nz,row;
256690b6cddSBarry Smith   PetscInt       *diag_offset=b->diag,idx,*ai=a->i,*aj=a->j,*pj;
257329f5518SBarry Smith   MatScalar      *pv,*v,*rtmp,m1,m2,m3,m4,*pc,*w,*x,x1,x2,x3,x4;
2582515f8d2SSatish Balay   MatScalar      p1,p2,p3,p4;
2593f1db9ecSBarry Smith   MatScalar      *ba   = b->a,*aa = a->a;
260182b8fbaSHong Zhang   PetscReal      shift = info->shiftamount;
261*2e92ee13SHong Zhang   PetscBool      zeropivotdetected;
2624eeb42bcSBarry Smith 
2633a40ed3dSBarry Smith   PetscFunctionBegin;
2644eeb42bcSBarry Smith   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
2654eeb42bcSBarry Smith   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
266785e854fSJed Brown   ierr = PetscMalloc1(4*(n+1),&rtmp);CHKERRQ(ierr);
2674eeb42bcSBarry Smith 
2684eeb42bcSBarry Smith   for (i=0; i<n; i++) {
2694078e994SBarry Smith     nz    = bi[i+1] - bi[i];
2704078e994SBarry Smith     ajtmp = bj + bi[i];
2714eeb42bcSBarry Smith     for  (j=0; j<nz; j++) {
2724eeb42bcSBarry Smith       x = rtmp+4*ajtmp[j]; x[0] = x[1] = x[2] = x[3] = 0.0;
2734eeb42bcSBarry Smith     }
2744eeb42bcSBarry Smith     /* load in initial (unfactored row) */
2754eeb42bcSBarry Smith     idx      = r[i];
2764078e994SBarry Smith     nz       = ai[idx+1] - ai[idx];
2774078e994SBarry Smith     ajtmpold = aj + ai[idx];
2784078e994SBarry Smith     v        = aa + 4*ai[idx];
2794eeb42bcSBarry Smith     for (j=0; j<nz; j++) {
2804eeb42bcSBarry Smith       x    = rtmp+4*ic[ajtmpold[j]];
2814eeb42bcSBarry Smith       x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3];
2824eeb42bcSBarry Smith       v   += 4;
2834eeb42bcSBarry Smith     }
2844eeb42bcSBarry Smith     row = *ajtmp++;
2854eeb42bcSBarry Smith     while (row < i) {
2864eeb42bcSBarry Smith       pc = rtmp + 4*row;
2874eeb42bcSBarry Smith       p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3];
288d4a378daSJed Brown       if (p1 != (PetscScalar)0.0 || p2 != (PetscScalar)0.0 || p3 != (PetscScalar)0.0 || p4 != (PetscScalar)0.0) {
2894078e994SBarry Smith         pv    = ba + 4*diag_offset[row];
2904078e994SBarry Smith         pj    = bj + diag_offset[row] + 1;
2914eeb42bcSBarry Smith         x1    = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
2924eeb42bcSBarry Smith         pc[0] = m1 = p1*x1 + p3*x2;
2934eeb42bcSBarry Smith         pc[1] = m2 = p2*x1 + p4*x2;
2944eeb42bcSBarry Smith         pc[2] = m3 = p1*x3 + p3*x4;
2954eeb42bcSBarry Smith         pc[3] = m4 = p2*x3 + p4*x4;
2964078e994SBarry Smith         nz    = bi[row+1] - diag_offset[row] - 1;
2974eeb42bcSBarry Smith         pv   += 4;
2984eeb42bcSBarry Smith         for (j=0; j<nz; j++) {
2994eeb42bcSBarry Smith           x1    = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
3004eeb42bcSBarry Smith           x     = rtmp + 4*pj[j];
3014eeb42bcSBarry Smith           x[0] -= m1*x1 + m3*x2;
3024eeb42bcSBarry Smith           x[1] -= m2*x1 + m4*x2;
3034eeb42bcSBarry Smith           x[2] -= m1*x3 + m3*x4;
3044eeb42bcSBarry Smith           x[3] -= m2*x3 + m4*x4;
3054eeb42bcSBarry Smith           pv   += 4;
3064eeb42bcSBarry Smith         }
307dc0b31edSSatish Balay         ierr = PetscLogFlops(16.0*nz+12.0);CHKERRQ(ierr);
3084eeb42bcSBarry Smith       }
3094eeb42bcSBarry Smith       row = *ajtmp++;
3104eeb42bcSBarry Smith     }
3114eeb42bcSBarry Smith     /* finished row so stick it into b->a */
3124078e994SBarry Smith     pv = ba + 4*bi[i];
3134078e994SBarry Smith     pj = bj + bi[i];
3144078e994SBarry Smith     nz = bi[i+1] - bi[i];
3154eeb42bcSBarry Smith     for (j=0; j<nz; j++) {
3164eeb42bcSBarry Smith       x     = rtmp+4*pj[j];
3174eeb42bcSBarry Smith       pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3];
3184eeb42bcSBarry Smith       pv   += 4;
3194eeb42bcSBarry Smith     }
3204eeb42bcSBarry Smith     /* invert diagonal block */
3214078e994SBarry Smith     w    = ba + 4*diag_offset[i];
322*2e92ee13SHong Zhang     ierr = PetscKernel_A_gets_inverse_A_2(w,shift,!A->erroriffailure,&zeropivotdetected);CHKERRQ(ierr);
323*2e92ee13SHong Zhang     if (zeropivotdetected) C->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
3244eeb42bcSBarry Smith   }
3254eeb42bcSBarry Smith 
326606d414cSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
3274eeb42bcSBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
3284eeb42bcSBarry Smith   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
32926fbe8dcSKarl Rupp 
33006e38f1dSHong Zhang   C->ops->solve          = MatSolve_SeqBAIJ_2_inplace;
33106e38f1dSHong Zhang   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_inplace;
3324eeb42bcSBarry Smith   C->assembled           = PETSC_TRUE;
33326fbe8dcSKarl Rupp 
334766f9fbaSBarry Smith   ierr = PetscLogFlops(1.333333333333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
3353a40ed3dSBarry Smith   PetscFunctionReturn(0);
3364eeb42bcSBarry Smith }
3374cd38560SBarry Smith /*
3384cd38560SBarry Smith       Version for when blocks are 2 by 2 Using natural ordering
3394cd38560SBarry Smith */
3404a2ae208SSatish Balay #undef __FUNCT__
34106e38f1dSHong Zhang #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_inplace"
34206e38f1dSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_inplace(Mat C,Mat A,const MatFactorInfo *info)
3434cd38560SBarry Smith {
3444cd38560SBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)C->data;
345dfbe8321SBarry Smith   PetscErrorCode ierr;
346690b6cddSBarry Smith   PetscInt       i,j,n = a->mbs,*bi = b->i,*bj = b->j;
347690b6cddSBarry Smith   PetscInt       *ajtmpold,*ajtmp,nz,row;
348690b6cddSBarry Smith   PetscInt       *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj;
349329f5518SBarry Smith   MatScalar      *pv,*v,*rtmp,*pc,*w,*x;
3502515f8d2SSatish Balay   MatScalar      p1,p2,p3,p4,m1,m2,m3,m4,x1,x2,x3,x4;
3514cd38560SBarry Smith   MatScalar      *ba   = b->a,*aa = a->a;
352182b8fbaSHong Zhang   PetscReal      shift = info->shiftamount;
353*2e92ee13SHong Zhang   PetscBool      zeropivotdetected;
3544cd38560SBarry Smith 
3554cd38560SBarry Smith   PetscFunctionBegin;
356785e854fSJed Brown   ierr = PetscMalloc1(4*(n+1),&rtmp);CHKERRQ(ierr);
3574cd38560SBarry Smith   for (i=0; i<n; i++) {
3584cd38560SBarry Smith     nz    = bi[i+1] - bi[i];
3594cd38560SBarry Smith     ajtmp = bj + bi[i];
3604cd38560SBarry Smith     for  (j=0; j<nz; j++) {
3614cd38560SBarry Smith       x    = rtmp+4*ajtmp[j];
3624cd38560SBarry Smith       x[0] = x[1]  = x[2]  = x[3]  = 0.0;
3634cd38560SBarry Smith     }
3644cd38560SBarry Smith     /* load in initial (unfactored row) */
3654cd38560SBarry Smith     nz       = ai[i+1] - ai[i];
3664cd38560SBarry Smith     ajtmpold = aj + ai[i];
3674cd38560SBarry Smith     v        = aa + 4*ai[i];
3684cd38560SBarry Smith     for (j=0; j<nz; j++) {
3694cd38560SBarry Smith       x    = rtmp+4*ajtmpold[j];
3704cd38560SBarry Smith       x[0] = v[0];  x[1]  = v[1];  x[2]  = v[2];  x[3]  = v[3];
3714cd38560SBarry Smith       v   += 4;
3724cd38560SBarry Smith     }
3734cd38560SBarry Smith     row = *ajtmp++;
3744cd38560SBarry Smith     while (row < i) {
3754cd38560SBarry Smith       pc = rtmp + 4*row;
3764cd38560SBarry Smith       p1 = pc[0];  p2  = pc[1];  p3  = pc[2];  p4  = pc[3];
377d4a378daSJed Brown       if (p1 != (PetscScalar)0.0 || p2 != (PetscScalar)0.0 || p3 != (PetscScalar)0.0 || p4 != (PetscScalar)0.0) {
3784cd38560SBarry Smith         pv    = ba + 4*diag_offset[row];
3794cd38560SBarry Smith         pj    = bj + diag_offset[row] + 1;
3804cd38560SBarry Smith         x1    = pv[0];  x2  = pv[1];  x3  = pv[2];  x4  = pv[3];
3814cd38560SBarry Smith         pc[0] = m1 = p1*x1 + p3*x2;
3824cd38560SBarry Smith         pc[1] = m2 = p2*x1 + p4*x2;
3834cd38560SBarry Smith         pc[2] = m3 = p1*x3 + p3*x4;
3844cd38560SBarry Smith         pc[3] = m4 = p2*x3 + p4*x4;
3854cd38560SBarry Smith         nz    = bi[row+1] - diag_offset[row] - 1;
3864cd38560SBarry Smith         pv   += 4;
3874cd38560SBarry Smith         for (j=0; j<nz; j++) {
3884cd38560SBarry Smith           x1    = pv[0];  x2  = pv[1];   x3 = pv[2];  x4  = pv[3];
3894cd38560SBarry Smith           x     = rtmp + 4*pj[j];
3904cd38560SBarry Smith           x[0] -= m1*x1 + m3*x2;
3914cd38560SBarry Smith           x[1] -= m2*x1 + m4*x2;
3924cd38560SBarry Smith           x[2] -= m1*x3 + m3*x4;
3934cd38560SBarry Smith           x[3] -= m2*x3 + m4*x4;
3944cd38560SBarry Smith           pv   += 4;
3954cd38560SBarry Smith         }
396dc0b31edSSatish Balay         ierr = PetscLogFlops(16.0*nz+12.0);CHKERRQ(ierr);
3974cd38560SBarry Smith       }
3984cd38560SBarry Smith       row = *ajtmp++;
3994cd38560SBarry Smith     }
4004cd38560SBarry Smith     /* finished row so stick it into b->a */
4014cd38560SBarry Smith     pv = ba + 4*bi[i];
4024cd38560SBarry Smith     pj = bj + bi[i];
4034cd38560SBarry Smith     nz = bi[i+1] - bi[i];
4044cd38560SBarry Smith     for (j=0; j<nz; j++) {
4054cd38560SBarry Smith       x     = rtmp+4*pj[j];
4064cd38560SBarry Smith       pv[0] = x[0];  pv[1]  = x[1];  pv[2]  = x[2];  pv[3]  = x[3];
4072efa7f71SHong Zhang       /*
4082efa7f71SHong Zhang       printf(" col %d:",pj[j]);
4092efa7f71SHong Zhang       PetscInt j1;
4102efa7f71SHong Zhang       for (j1=0; j1<4; j1++) printf(" %g,",*(pv+j1));
4112efa7f71SHong Zhang       printf("\n");
4122efa7f71SHong Zhang       */
4134cd38560SBarry Smith       pv += 4;
4144cd38560SBarry Smith     }
4154cd38560SBarry Smith     /* invert diagonal block */
4164cd38560SBarry Smith     w = ba + 4*diag_offset[i];
4172efa7f71SHong Zhang     /*
4182efa7f71SHong Zhang     printf(" \n%d -th: diag: ",i);
4192efa7f71SHong Zhang     for (j=0; j<4; j++) {
4202efa7f71SHong Zhang       printf(" %g,",w[j]);
4212efa7f71SHong Zhang     }
4222efa7f71SHong Zhang     printf("\n----------------------------\n");
4232efa7f71SHong Zhang     */
424*2e92ee13SHong Zhang     ierr = PetscKernel_A_gets_inverse_A_2(w,shift,!A->erroriffailure,&zeropivotdetected);CHKERRQ(ierr);
425*2e92ee13SHong Zhang     if (zeropivotdetected) C->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
4264cd38560SBarry Smith   }
4274cd38560SBarry Smith 
428606d414cSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
42926fbe8dcSKarl Rupp 
43006e38f1dSHong Zhang   C->ops->solve          = MatSolve_SeqBAIJ_2_NaturalOrdering_inplace;
43106e38f1dSHong Zhang   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_NaturalOrdering_inplace;
4324cd38560SBarry Smith   C->assembled           = PETSC_TRUE;
43326fbe8dcSKarl Rupp 
434766f9fbaSBarry Smith   ierr = PetscLogFlops(1.333333333333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
4354cd38560SBarry Smith   PetscFunctionReturn(0);
4364cd38560SBarry Smith }
4377fc0212eSBarry Smith 
4387fc0212eSBarry Smith /* ----------------------------------------------------------- */
4397fc0212eSBarry Smith /*
4407fc0212eSBarry Smith      Version for when blocks are 1 by 1.
4417fc0212eSBarry Smith */
4424a2ae208SSatish Balay #undef __FUNCT__
443048b5e81SShri Abhyankar #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_1"
444048b5e81SShri Abhyankar PetscErrorCode MatLUFactorNumeric_SeqBAIJ_1(Mat B,Mat A,const MatFactorInfo *info)
445048b5e81SShri Abhyankar {
446048b5e81SShri Abhyankar   Mat             C     =B;
447048b5e81SShri Abhyankar   Mat_SeqBAIJ     *a    =(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data;
448048b5e81SShri Abhyankar   IS              isrow = b->row,isicol = b->icol;
449048b5e81SShri Abhyankar   PetscErrorCode  ierr;
450048b5e81SShri Abhyankar   const PetscInt  *r,*ic,*ics;
451048b5e81SShri Abhyankar   const PetscInt  n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bdiag=b->diag;
452048b5e81SShri Abhyankar   PetscInt        i,j,k,nz,nzL,row,*pj;
453048b5e81SShri Abhyankar   const PetscInt  *ajtmp,*bjtmp;
454048b5e81SShri Abhyankar   MatScalar       *rtmp,*pc,multiplier,*pv;
455048b5e81SShri Abhyankar   const MatScalar *aa=a->a,*v;
456ace3abfcSBarry Smith   PetscBool       row_identity,col_identity;
457048b5e81SShri Abhyankar   FactorShiftCtx  sctx;
458048b5e81SShri Abhyankar   const PetscInt  *ddiag;
459048b5e81SShri Abhyankar   PetscReal       rs;
460048b5e81SShri Abhyankar   MatScalar       d;
461048b5e81SShri Abhyankar 
462048b5e81SShri Abhyankar   PetscFunctionBegin;
463048b5e81SShri Abhyankar   /* MatPivotSetUp(): initialize shift context sctx */
464048b5e81SShri Abhyankar   ierr = PetscMemzero(&sctx,sizeof(FactorShiftCtx));CHKERRQ(ierr);
465048b5e81SShri Abhyankar 
466048b5e81SShri Abhyankar   if (info->shifttype == (PetscReal) MAT_SHIFT_POSITIVE_DEFINITE) { /* set sctx.shift_top=max{rs} */
467048b5e81SShri Abhyankar     ddiag          = a->diag;
468048b5e81SShri Abhyankar     sctx.shift_top = info->zeropivot;
469048b5e81SShri Abhyankar     for (i=0; i<n; i++) {
470048b5e81SShri Abhyankar       /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
471048b5e81SShri Abhyankar       d  = (aa)[ddiag[i]];
472048b5e81SShri Abhyankar       rs = -PetscAbsScalar(d) - PetscRealPart(d);
473048b5e81SShri Abhyankar       v  = aa+ai[i];
474048b5e81SShri Abhyankar       nz = ai[i+1] - ai[i];
47526fbe8dcSKarl Rupp       for (j=0; j<nz; j++) rs += PetscAbsScalar(v[j]);
476048b5e81SShri Abhyankar       if (rs>sctx.shift_top) sctx.shift_top = rs;
477048b5e81SShri Abhyankar     }
478048b5e81SShri Abhyankar     sctx.shift_top *= 1.1;
479048b5e81SShri Abhyankar     sctx.nshift_max = 5;
480048b5e81SShri Abhyankar     sctx.shift_lo   = 0.;
481048b5e81SShri Abhyankar     sctx.shift_hi   = 1.;
482048b5e81SShri Abhyankar   }
483048b5e81SShri Abhyankar 
484048b5e81SShri Abhyankar   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
485048b5e81SShri Abhyankar   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
486854ce69bSBarry Smith   ierr = PetscMalloc1(n+1,&rtmp);CHKERRQ(ierr);
487048b5e81SShri Abhyankar   ics  = ic;
488048b5e81SShri Abhyankar 
489048b5e81SShri Abhyankar   do {
490048b5e81SShri Abhyankar     sctx.newshift = PETSC_FALSE;
491048b5e81SShri Abhyankar     for (i=0; i<n; i++) {
492048b5e81SShri Abhyankar       /* zero rtmp */
493048b5e81SShri Abhyankar       /* L part */
494048b5e81SShri Abhyankar       nz    = bi[i+1] - bi[i];
495048b5e81SShri Abhyankar       bjtmp = bj + bi[i];
496048b5e81SShri Abhyankar       for  (j=0; j<nz; j++) rtmp[bjtmp[j]] = 0.0;
497048b5e81SShri Abhyankar 
498048b5e81SShri Abhyankar       /* U part */
499048b5e81SShri Abhyankar       nz    = bdiag[i]-bdiag[i+1];
500048b5e81SShri Abhyankar       bjtmp = bj + bdiag[i+1]+1;
501048b5e81SShri Abhyankar       for  (j=0; j<nz; j++) rtmp[bjtmp[j]] = 0.0;
502048b5e81SShri Abhyankar 
503048b5e81SShri Abhyankar       /* load in initial (unfactored row) */
504048b5e81SShri Abhyankar       nz    = ai[r[i]+1] - ai[r[i]];
505048b5e81SShri Abhyankar       ajtmp = aj + ai[r[i]];
506048b5e81SShri Abhyankar       v     = aa + ai[r[i]];
50726fbe8dcSKarl Rupp       for (j=0; j<nz; j++) rtmp[ics[ajtmp[j]]] = v[j];
50826fbe8dcSKarl Rupp 
509048b5e81SShri Abhyankar       /* ZeropivotApply() */
510048b5e81SShri Abhyankar       rtmp[i] += sctx.shift_amount;  /* shift the diagonal of the matrix */
511048b5e81SShri Abhyankar 
512048b5e81SShri Abhyankar       /* elimination */
513048b5e81SShri Abhyankar       bjtmp = bj + bi[i];
514048b5e81SShri Abhyankar       row   = *bjtmp++;
515048b5e81SShri Abhyankar       nzL   = bi[i+1] - bi[i];
516048b5e81SShri Abhyankar       for (k=0; k < nzL; k++) {
517048b5e81SShri Abhyankar         pc = rtmp + row;
518d4a378daSJed Brown         if (*pc != (PetscScalar)0.0) {
519048b5e81SShri Abhyankar           pv         = b->a + bdiag[row];
520048b5e81SShri Abhyankar           multiplier = *pc * (*pv);
521048b5e81SShri Abhyankar           *pc        = multiplier;
52226fbe8dcSKarl Rupp 
523048b5e81SShri Abhyankar           pj = b->j + bdiag[row+1]+1; /* beginning of U(row,:) */
524048b5e81SShri Abhyankar           pv = b->a + bdiag[row+1]+1;
525048b5e81SShri Abhyankar           nz = bdiag[row]-bdiag[row+1]-1; /* num of entries in U(row,:) excluding diag */
526048b5e81SShri Abhyankar           for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j];
527048b5e81SShri Abhyankar           ierr = PetscLogFlops(2.0*nz);CHKERRQ(ierr);
528048b5e81SShri Abhyankar         }
529048b5e81SShri Abhyankar         row = *bjtmp++;
530048b5e81SShri Abhyankar       }
531048b5e81SShri Abhyankar 
532048b5e81SShri Abhyankar       /* finished row so stick it into b->a */
533048b5e81SShri Abhyankar       rs = 0.0;
534048b5e81SShri Abhyankar       /* L part */
535048b5e81SShri Abhyankar       pv = b->a + bi[i];
536048b5e81SShri Abhyankar       pj = b->j + bi[i];
537048b5e81SShri Abhyankar       nz = bi[i+1] - bi[i];
538048b5e81SShri Abhyankar       for (j=0; j<nz; j++) {
539048b5e81SShri Abhyankar         pv[j] = rtmp[pj[j]]; rs += PetscAbsScalar(pv[j]);
540048b5e81SShri Abhyankar       }
541048b5e81SShri Abhyankar 
542048b5e81SShri Abhyankar       /* U part */
543048b5e81SShri Abhyankar       pv = b->a + bdiag[i+1]+1;
544048b5e81SShri Abhyankar       pj = b->j + bdiag[i+1]+1;
545048b5e81SShri Abhyankar       nz = bdiag[i] - bdiag[i+1]-1;
546048b5e81SShri Abhyankar       for (j=0; j<nz; j++) {
547048b5e81SShri Abhyankar         pv[j] = rtmp[pj[j]]; rs += PetscAbsScalar(pv[j]);
548048b5e81SShri Abhyankar       }
549048b5e81SShri Abhyankar 
550048b5e81SShri Abhyankar       sctx.rs = rs;
551048b5e81SShri Abhyankar       sctx.pv = rtmp[i];
5524cccfbddSHong Zhang       ierr    = MatPivotCheck(B,A,info,&sctx,i);CHKERRQ(ierr);
553048b5e81SShri Abhyankar       if (sctx.newshift) break; /* break for-loop */
554048b5e81SShri Abhyankar       rtmp[i] = sctx.pv; /* sctx.pv might be updated in the case of MAT_SHIFT_INBLOCKS */
555048b5e81SShri Abhyankar 
556048b5e81SShri Abhyankar       /* Mark diagonal and invert diagonal for simplier triangular solves */
557048b5e81SShri Abhyankar       pv  = b->a + bdiag[i];
558d4a378daSJed Brown       *pv = (PetscScalar)1.0/rtmp[i];
559048b5e81SShri Abhyankar 
560048b5e81SShri Abhyankar     } /* endof for (i=0; i<n; i++) { */
561048b5e81SShri Abhyankar 
562048b5e81SShri Abhyankar     /* MatPivotRefine() */
563048b5e81SShri Abhyankar     if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE && !sctx.newshift && sctx.shift_fraction>0 && sctx.nshift<sctx.nshift_max) {
564048b5e81SShri Abhyankar       /*
565048b5e81SShri Abhyankar        * if no shift in this attempt & shifting & started shifting & can refine,
566048b5e81SShri Abhyankar        * then try lower shift
567048b5e81SShri Abhyankar        */
568048b5e81SShri Abhyankar       sctx.shift_hi       = sctx.shift_fraction;
569048b5e81SShri Abhyankar       sctx.shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.;
570048b5e81SShri Abhyankar       sctx.shift_amount   = sctx.shift_fraction * sctx.shift_top;
571048b5e81SShri Abhyankar       sctx.newshift       = PETSC_TRUE;
572048b5e81SShri Abhyankar       sctx.nshift++;
573048b5e81SShri Abhyankar     }
574048b5e81SShri Abhyankar   } while (sctx.newshift);
575048b5e81SShri Abhyankar 
576048b5e81SShri Abhyankar   ierr = PetscFree(rtmp);CHKERRQ(ierr);
577048b5e81SShri Abhyankar   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
578048b5e81SShri Abhyankar   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
579048b5e81SShri Abhyankar 
580048b5e81SShri Abhyankar   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
581048b5e81SShri Abhyankar   ierr = ISIdentity(isicol,&col_identity);CHKERRQ(ierr);
582048b5e81SShri Abhyankar   if (row_identity && col_identity) {
583048b5e81SShri Abhyankar     C->ops->solve          = MatSolve_SeqBAIJ_1_NaturalOrdering;
5849f5c0bcdSBarry Smith     C->ops->forwardsolve   = MatForwardSolve_SeqBAIJ_1_NaturalOrdering;
5859f5c0bcdSBarry Smith     C->ops->backwardsolve  = MatBackwardSolve_SeqBAIJ_1_NaturalOrdering;
58693fd935bSShri Abhyankar     C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1_NaturalOrdering;
587048b5e81SShri Abhyankar   } else {
588048b5e81SShri Abhyankar     C->ops->solve          = MatSolve_SeqBAIJ_1;
58993fd935bSShri Abhyankar     C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1;
590048b5e81SShri Abhyankar   }
591048b5e81SShri Abhyankar   C->assembled = PETSC_TRUE;
592048b5e81SShri Abhyankar   ierr         = PetscLogFlops(C->cmap->n);CHKERRQ(ierr);
593048b5e81SShri Abhyankar 
594048b5e81SShri Abhyankar   /* MatShiftView(A,info,&sctx) */
595048b5e81SShri Abhyankar   if (sctx.nshift) {
596048b5e81SShri Abhyankar     if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
59757622a8eSBarry Smith       ierr = PetscInfo4(A,"number of shift_pd tries %D, shift_amount %g, diagonal shifted up by %e fraction top_value %e\n",sctx.nshift,(double)sctx.shift_amount,(double)sctx.shift_fraction,(double)sctx.shift_top);CHKERRQ(ierr);
598048b5e81SShri Abhyankar     } else if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
59957622a8eSBarry Smith       ierr = PetscInfo2(A,"number of shift_nz tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr);
600048b5e81SShri Abhyankar     } else if (info->shifttype == (PetscReal)MAT_SHIFT_INBLOCKS) {
60157622a8eSBarry Smith       ierr = PetscInfo2(A,"number of shift_inblocks applied %D, each shift_amount %g\n",sctx.nshift,(double)info->shiftamount);CHKERRQ(ierr);
602048b5e81SShri Abhyankar     }
603048b5e81SShri Abhyankar   }
604048b5e81SShri Abhyankar   PetscFunctionReturn(0);
605048b5e81SShri Abhyankar }
606048b5e81SShri Abhyankar 
607048b5e81SShri Abhyankar #undef __FUNCT__
60806e38f1dSHong Zhang #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_1_inplace"
60906e38f1dSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqBAIJ_1_inplace(Mat C,Mat A,const MatFactorInfo *info)
6107fc0212eSBarry Smith {
6117fc0212eSBarry Smith   Mat_SeqBAIJ    *a    = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)C->data;
6127cdf2dbbSSatish Balay   IS             isrow = b->row,isicol = b->icol;
6136849ba73SBarry Smith   PetscErrorCode ierr;
6145d0c19d7SBarry Smith   const PetscInt *r,*ic;
6155d0c19d7SBarry Smith   PetscInt       i,j,n = a->mbs,*bi = b->i,*bj = b->j;
616690b6cddSBarry Smith   PetscInt       *ajtmpold,*ajtmp,nz,row,*ai = a->i,*aj = a->j;
617690b6cddSBarry Smith   PetscInt       *diag_offset = b->diag,diag,*pj;
618329f5518SBarry Smith   MatScalar      *pv,*v,*rtmp,multiplier,*pc;
6193f1db9ecSBarry Smith   MatScalar      *ba = b->a,*aa = a->a;
620ace3abfcSBarry Smith   PetscBool      row_identity, col_identity;
6217fc0212eSBarry Smith 
6223a40ed3dSBarry Smith   PetscFunctionBegin;
6237fc0212eSBarry Smith   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
6247fc0212eSBarry Smith   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
625854ce69bSBarry Smith   ierr = PetscMalloc1(n+1,&rtmp);CHKERRQ(ierr);
6267fc0212eSBarry Smith 
6277fc0212eSBarry Smith   for (i=0; i<n; i++) {
6284078e994SBarry Smith     nz    = bi[i+1] - bi[i];
6294078e994SBarry Smith     ajtmp = bj + bi[i];
6307fc0212eSBarry Smith     for  (j=0; j<nz; j++) rtmp[ajtmp[j]] = 0.0;
6317fc0212eSBarry Smith 
6327fc0212eSBarry Smith     /* load in initial (unfactored row) */
6334078e994SBarry Smith     nz       = ai[r[i]+1] - ai[r[i]];
6344078e994SBarry Smith     ajtmpold = aj + ai[r[i]];
6354078e994SBarry Smith     v        = aa + ai[r[i]];
6367fc0212eSBarry Smith     for (j=0; j<nz; j++) rtmp[ic[ajtmpold[j]]] =  v[j];
6377fc0212eSBarry Smith 
6387fc0212eSBarry Smith     row = *ajtmp++;
6397fc0212eSBarry Smith     while (row < i) {
6407fc0212eSBarry Smith       pc = rtmp + row;
6417fc0212eSBarry Smith       if (*pc != 0.0) {
6424078e994SBarry Smith         pv         = ba + diag_offset[row];
6434078e994SBarry Smith         pj         = bj + diag_offset[row] + 1;
6447fc0212eSBarry Smith         multiplier = *pc * *pv++;
6457fc0212eSBarry Smith         *pc        = multiplier;
6464078e994SBarry Smith         nz         = bi[row+1] - diag_offset[row] - 1;
6477fc0212eSBarry Smith         for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j];
648dc0b31edSSatish Balay         ierr = PetscLogFlops(1.0+2.0*nz);CHKERRQ(ierr);
6497fc0212eSBarry Smith       }
6507fc0212eSBarry Smith       row = *ajtmp++;
6517fc0212eSBarry Smith     }
6527fc0212eSBarry Smith     /* finished row so stick it into b->a */
6534078e994SBarry Smith     pv = ba + bi[i];
6544078e994SBarry Smith     pj = bj + bi[i];
6554078e994SBarry Smith     nz = bi[i+1] - bi[i];
65626fbe8dcSKarl Rupp     for (j=0; j<nz; j++) pv[j] = rtmp[pj[j]];
6574078e994SBarry Smith     diag = diag_offset[i] - bi[i];
6587fc0212eSBarry Smith     /* check pivot entry for current row */
65965e19b50SBarry Smith     if (pv[diag] == 0.0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot: row in original ordering %D in permuted ordering %D",r[i],i);
6607fc0212eSBarry Smith     pv[diag] = 1.0/pv[diag];
6617fc0212eSBarry Smith   }
6627fc0212eSBarry Smith 
663606d414cSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
6647fc0212eSBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
6657fc0212eSBarry Smith   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
666f58c8c74SBarry Smith   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
667f58c8c74SBarry Smith   ierr = ISIdentity(isicol,&col_identity);CHKERRQ(ierr);
668f58c8c74SBarry Smith   if (row_identity && col_identity) {
66906e38f1dSHong Zhang     C->ops->solve          = MatSolve_SeqBAIJ_1_NaturalOrdering_inplace;
67006e38f1dSHong Zhang     C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1_NaturalOrdering_inplace;
671f58c8c74SBarry Smith   } else {
67206e38f1dSHong Zhang     C->ops->solve          = MatSolve_SeqBAIJ_1_inplace;
67306e38f1dSHong Zhang     C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1_inplace;
674f58c8c74SBarry Smith   }
6757fc0212eSBarry Smith   C->assembled = PETSC_TRUE;
676d0f46423SBarry Smith   ierr         = PetscLogFlops(C->cmap->n);CHKERRQ(ierr);
6773a40ed3dSBarry Smith   PetscFunctionReturn(0);
6787fc0212eSBarry Smith }
6797fc0212eSBarry Smith 
680b24902e0SBarry Smith #undef __FUNCT__
681b24902e0SBarry Smith #define __FUNCT__ "MatGetFactor_seqbaij_petsc"
6829a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatGetFactor_seqbaij_petsc(Mat A,const MatFactorType ftype,Mat *B)
683b24902e0SBarry Smith {
684d0f46423SBarry Smith   PetscInt       n = A->rmap->n;
685b24902e0SBarry Smith   PetscErrorCode ierr;
686b24902e0SBarry Smith 
687b24902e0SBarry Smith   PetscFunctionBegin;
688b499a2aeSHong Zhang #if defined(PETSC_USE_COMPLEX)
689b499a2aeSHong Zhang   if (A->hermitian && (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Hermitian Factor is not supported");
690b499a2aeSHong Zhang #endif
691ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),B);CHKERRQ(ierr);
692b24902e0SBarry Smith   ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr);
693c8342467SHong Zhang   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU || ftype == MAT_FACTOR_ILUDT) {
694b24902e0SBarry Smith     ierr = MatSetType(*B,MATSEQBAIJ);CHKERRQ(ierr);
69526fbe8dcSKarl Rupp 
6964dd39f65SShri Abhyankar     (*B)->ops->lufactorsymbolic  = MatLUFactorSymbolic_SeqBAIJ;
6974dd39f65SShri Abhyankar     (*B)->ops->ilufactorsymbolic = MatILUFactorSymbolic_SeqBAIJ;
698b24902e0SBarry Smith   } else if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) {
699b24902e0SBarry Smith     ierr = MatSetType(*B,MATSEQSBAIJ);CHKERRQ(ierr);
7000298fd71SBarry Smith     ierr = MatSeqSBAIJSetPreallocation(*B,A->rmap->bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr);
70126fbe8dcSKarl Rupp 
7025c9eb25fSBarry Smith     (*B)->ops->iccfactorsymbolic      = MatICCFactorSymbolic_SeqBAIJ;
7035c9eb25fSBarry Smith     (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqBAIJ;
704e32f2f54SBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
705d5f3da31SBarry Smith   (*B)->factortype = ftype;
706b24902e0SBarry Smith   PetscFunctionReturn(0);
707b24902e0SBarry Smith }
708273d9f13SBarry Smith 
7095d517e7eSBarry Smith /* ----------------------------------------------------------- */
7104a2ae208SSatish Balay #undef __FUNCT__
7114a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqBAIJ"
7120481f469SBarry Smith PetscErrorCode MatLUFactor_SeqBAIJ(Mat A,IS row,IS col,const MatFactorInfo *info)
7135d517e7eSBarry Smith {
714dfbe8321SBarry Smith   PetscErrorCode ierr;
7155d517e7eSBarry Smith   Mat            C;
7165d517e7eSBarry Smith 
7173a40ed3dSBarry Smith   PetscFunctionBegin;
7182692d6eeSBarry Smith   ierr = MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_LU,&C);CHKERRQ(ierr);
719719d5645SBarry Smith   ierr = MatLUFactorSymbolic(C,A,row,col,info);CHKERRQ(ierr);
720719d5645SBarry Smith   ierr = MatLUFactorNumeric(C,A,info);CHKERRQ(ierr);
72126fbe8dcSKarl Rupp 
722db4efbfdSBarry Smith   A->ops->solve          = C->ops->solve;
723db4efbfdSBarry Smith   A->ops->solvetranspose = C->ops->solvetranspose;
72426fbe8dcSKarl Rupp 
72528be2f97SBarry Smith   ierr = MatHeaderMerge(A,&C);CHKERRQ(ierr);
7263bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)A,(PetscObject)((Mat_SeqBAIJ*)(A->data))->icol);CHKERRQ(ierr);
7273a40ed3dSBarry Smith   PetscFunctionReturn(0);
7285d517e7eSBarry Smith }
7294cd38560SBarry Smith 
730c6db04a5SJed Brown #include <../src/mat/impls/sbaij/seq/sbaij.h>
731c05c3958SHong Zhang #undef __FUNCT__
732c05c3958SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N"
7330481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N(Mat C,Mat A,const MatFactorInfo *info)
734c05c3958SHong Zhang {
735f3086b4bSHong Zhang   PetscErrorCode ierr;
73678910aadSHong Zhang   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data;
73778910aadSHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
73878910aadSHong Zhang   IS             ip=b->row;
7395d0c19d7SBarry Smith   const PetscInt *rip;
7405d0c19d7SBarry Smith   PetscInt       i,j,mbs=a->mbs,bs=A->rmap->bs,*bi=b->i,*bj=b->j,*bcol;
74178910aadSHong Zhang   PetscInt       *ai=a->i,*aj=a->j;
74278910aadSHong Zhang   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
74378910aadSHong Zhang   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
74407b50cabSHong Zhang   PetscReal      rs;
7450e95ead3SHong Zhang   FactorShiftCtx sctx;
74678910aadSHong Zhang 
747c05c3958SHong Zhang   PetscFunctionBegin;
74807b50cabSHong Zhang   if (bs > 1) { /* convert A to a SBAIJ matrix and apply Cholesky factorization from it */
7496ad2eaddSHong Zhang     if (!a->sbaijMat) {
750ceb03754SKris Buschelman       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
7516ad2eaddSHong Zhang     }
752719d5645SBarry Smith     ierr = (a->sbaijMat)->ops->choleskyfactornumeric(C,a->sbaijMat,info);CHKERRQ(ierr);
7536bf464f9SBarry Smith     ierr = MatDestroy(&a->sbaijMat);CHKERRQ(ierr);
7546ad2eaddSHong Zhang     PetscFunctionReturn(0);
7556ad2eaddSHong Zhang   }
75678910aadSHong Zhang 
75707b50cabSHong Zhang   /* MatPivotSetUp(): initialize shift context sctx */
75807b50cabSHong Zhang   ierr = PetscMemzero(&sctx,sizeof(FactorShiftCtx));CHKERRQ(ierr);
7593cea4cbeSHong Zhang 
7606ad2eaddSHong Zhang   ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr);
761dcca6d9dSJed Brown   ierr = PetscMalloc3(mbs,&rtmp,mbs,&il,mbs,&jl);CHKERRQ(ierr);
76278910aadSHong Zhang 
76375567043SBarry Smith   sctx.shift_amount = 0.;
7643cea4cbeSHong Zhang   sctx.nshift       = 0;
76578910aadSHong Zhang   do {
76607b50cabSHong Zhang     sctx.newshift = PETSC_FALSE;
76778910aadSHong Zhang     for (i=0; i<mbs; i++) {
768e60cf9a0SBarry Smith       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
76978910aadSHong Zhang     }
77078910aadSHong Zhang 
77178910aadSHong Zhang     for (k = 0; k<mbs; k++) {
77278910aadSHong Zhang       bval = ba + bi[k];
77378910aadSHong Zhang       /* initialize k-th row by the perm[k]-th row of A */
77478910aadSHong Zhang       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
77578910aadSHong Zhang       for (j = jmin; j < jmax; j++) {
77678910aadSHong Zhang         col = rip[aj[j]];
77778910aadSHong Zhang         if (col >= k) { /* only take upper triangular entry */
77878910aadSHong Zhang           rtmp[col] = aa[j];
77978910aadSHong Zhang           *bval++   = 0.0; /* for in-place factorization */
78078910aadSHong Zhang         }
78178910aadSHong Zhang       }
7823cea4cbeSHong Zhang 
7833cea4cbeSHong Zhang       /* shift the diagonal of the matrix */
7843cea4cbeSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
78578910aadSHong Zhang 
78678910aadSHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
78778910aadSHong Zhang       dk = rtmp[k];
78878910aadSHong Zhang       i  = jl[k]; /* first row to be added to k_th row  */
78978910aadSHong Zhang 
79078910aadSHong Zhang       while (i < k) {
79178910aadSHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
79278910aadSHong Zhang 
79378910aadSHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
79478910aadSHong Zhang         ili     = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
79578910aadSHong Zhang         uikdi   = -ba[ili]*ba[bi[i]]; /* diagonal(k) */
79678910aadSHong Zhang         dk     += uikdi*ba[ili];
79778910aadSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
79878910aadSHong Zhang 
79978910aadSHong Zhang         /* add multiple of row i to k-th row */
80078910aadSHong Zhang         jmin = ili + 1; jmax = bi[i+1];
80178910aadSHong Zhang         if (jmin < jmax) {
80278910aadSHong Zhang           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
80378910aadSHong Zhang           /* update il and jl for row i */
80478910aadSHong Zhang           il[i] = jmin;
80578910aadSHong Zhang           j     = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
80678910aadSHong Zhang         }
80778910aadSHong Zhang         i = nexti;
80878910aadSHong Zhang       }
80978910aadSHong Zhang 
8103cea4cbeSHong Zhang       /* shift the diagonals when zero pivot is detected */
8113cea4cbeSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
8123cea4cbeSHong Zhang       rs   = 0.0;
81378910aadSHong Zhang       jmin = bi[k]+1;
81478910aadSHong Zhang       nz   = bi[k+1] - jmin;
81578910aadSHong Zhang       if (nz) {
81678910aadSHong Zhang         bcol = bj + jmin;
81778910aadSHong Zhang         while (nz--) {
81889f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
81989f845c8SHong Zhang           bcol++;
82078910aadSHong Zhang         }
82178910aadSHong Zhang       }
8223cea4cbeSHong Zhang 
8233cea4cbeSHong Zhang       sctx.rs = rs;
8243cea4cbeSHong Zhang       sctx.pv = dk;
8254cccfbddSHong Zhang       ierr    = MatPivotCheck(C,A,info,&sctx,k);CHKERRQ(ierr);
82607b50cabSHong Zhang       if (sctx.newshift) break;
8270e95ead3SHong Zhang       dk = sctx.pv;
82878910aadSHong Zhang 
82978910aadSHong Zhang       /* copy data into U(k,:) */
83078910aadSHong Zhang       ba[bi[k]] = 1.0/dk; /* U(k,k) */
83178910aadSHong Zhang       jmin      = bi[k]+1; jmax = bi[k+1];
83278910aadSHong Zhang       if (jmin < jmax) {
83378910aadSHong Zhang         for (j=jmin; j<jmax; j++) {
83478910aadSHong Zhang           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
83578910aadSHong Zhang         }
83678910aadSHong Zhang         /* add the k-th row into il and jl */
83778910aadSHong Zhang         il[k] = jmin;
83878910aadSHong Zhang         i     = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
83978910aadSHong Zhang       }
84078910aadSHong Zhang     }
84107b50cabSHong Zhang   } while (sctx.newshift);
842fca92195SBarry Smith   ierr = PetscFree3(rtmp,il,jl);CHKERRQ(ierr);
84378910aadSHong Zhang 
84478910aadSHong Zhang   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
84526fbe8dcSKarl Rupp 
84678910aadSHong Zhang   C->assembled    = PETSC_TRUE;
84778910aadSHong Zhang   C->preallocated = PETSC_TRUE;
84826fbe8dcSKarl Rupp 
849d0f46423SBarry Smith   ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr);
8503cea4cbeSHong Zhang   if (sctx.nshift) {
851f4db908eSBarry Smith     if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
85257622a8eSBarry Smith       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr);
853f4db908eSBarry Smith     } else if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
85457622a8eSBarry Smith       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr);
85578910aadSHong Zhang     }
85678910aadSHong Zhang   }
857c05c3958SHong Zhang   PetscFunctionReturn(0);
858c05c3958SHong Zhang }
8594cd38560SBarry Smith 
860c05c3958SHong Zhang #undef __FUNCT__
861c05c3958SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering"
8620481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info)
863c05c3958SHong Zhang {
86478910aadSHong Zhang   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data;
86578910aadSHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
86678910aadSHong Zhang   PetscErrorCode ierr;
8673cea4cbeSHong Zhang   PetscInt       i,j,am=a->mbs;
86878910aadSHong Zhang   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
8693cea4cbeSHong Zhang   PetscInt       k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz;
87078910aadSHong Zhang   MatScalar      *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval;
8710e95ead3SHong Zhang   PetscReal      rs;
8720e95ead3SHong Zhang   FactorShiftCtx sctx;
87378910aadSHong Zhang 
874c05c3958SHong Zhang   PetscFunctionBegin;
8750e95ead3SHong Zhang   /* MatPivotSetUp(): initialize shift context sctx */
8760e95ead3SHong Zhang   ierr = PetscMemzero(&sctx,sizeof(FactorShiftCtx));CHKERRQ(ierr);
8773cea4cbeSHong Zhang 
878dcca6d9dSJed Brown   ierr = PetscMalloc3(am,&rtmp,am,&il,am,&jl);CHKERRQ(ierr);
87978910aadSHong Zhang 
88078910aadSHong Zhang   do {
88107b50cabSHong Zhang     sctx.newshift = PETSC_FALSE;
88278910aadSHong Zhang     for (i=0; i<am; i++) {
883e60cf9a0SBarry Smith       rtmp[i] = 0.0; jl[i] = am; il[0] = 0;
88478910aadSHong Zhang     }
88578910aadSHong Zhang 
88678910aadSHong Zhang     for (k = 0; k<am; k++) {
88778910aadSHong Zhang       /* initialize k-th row with elements nonzero in row perm(k) of A */
88878910aadSHong Zhang       nz   = ai[k+1] - ai[k];
88978910aadSHong Zhang       acol = aj + ai[k];
89078910aadSHong Zhang       aval = aa + ai[k];
89178910aadSHong Zhang       bval = ba + bi[k];
89278910aadSHong Zhang       while (nz--) {
89378910aadSHong Zhang         if (*acol < k) { /* skip lower triangular entries */
89478910aadSHong Zhang           acol++; aval++;
89578910aadSHong Zhang         } else {
89678910aadSHong Zhang           rtmp[*acol++] = *aval++;
89778910aadSHong Zhang           *bval++       = 0.0; /* for in-place factorization */
89878910aadSHong Zhang         }
89978910aadSHong Zhang       }
9003cea4cbeSHong Zhang 
9013cea4cbeSHong Zhang       /* shift the diagonal of the matrix */
9023cea4cbeSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
90378910aadSHong Zhang 
90478910aadSHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
90578910aadSHong Zhang       dk = rtmp[k];
90678910aadSHong Zhang       i  = jl[k]; /* first row to be added to k_th row  */
90778910aadSHong Zhang 
90878910aadSHong Zhang       while (i < k) {
90978910aadSHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
91078910aadSHong Zhang         /* compute multiplier, update D(k) and U(i,k) */
91178910aadSHong Zhang         ili     = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
91278910aadSHong Zhang         uikdi   = -ba[ili]*ba[bi[i]];
91378910aadSHong Zhang         dk     += uikdi*ba[ili];
91478910aadSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
91578910aadSHong Zhang 
91678910aadSHong Zhang         /* add multiple of row i to k-th row ... */
91778910aadSHong Zhang         jmin = ili + 1;
91878910aadSHong Zhang         nz   = bi[i+1] - jmin;
91978910aadSHong Zhang         if (nz > 0) {
92078910aadSHong Zhang           bcol = bj + jmin;
92178910aadSHong Zhang           bval = ba + jmin;
92278910aadSHong Zhang           while (nz--) rtmp[*bcol++] += uikdi*(*bval++);
92378910aadSHong Zhang           /* update il and jl for i-th row */
92478910aadSHong Zhang           il[i] = jmin;
92578910aadSHong Zhang           j     = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
92678910aadSHong Zhang         }
92778910aadSHong Zhang         i = nexti;
92878910aadSHong Zhang       }
92978910aadSHong Zhang 
9303cea4cbeSHong Zhang       /* shift the diagonals when zero pivot is detected */
9313cea4cbeSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
9323cea4cbeSHong Zhang       rs   = 0.0;
93378910aadSHong Zhang       jmin = bi[k]+1;
93478910aadSHong Zhang       nz   = bi[k+1] - jmin;
93578910aadSHong Zhang       if (nz) {
93678910aadSHong Zhang         bcol = bj + jmin;
93778910aadSHong Zhang         while (nz--) {
93889f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
93989f845c8SHong Zhang           bcol++;
94078910aadSHong Zhang         }
94178910aadSHong Zhang       }
9423cea4cbeSHong Zhang 
9433cea4cbeSHong Zhang       sctx.rs = rs;
9443cea4cbeSHong Zhang       sctx.pv = dk;
9454cccfbddSHong Zhang       ierr    = MatPivotCheck(C,A,info,&sctx,k);CHKERRQ(ierr);
94607b50cabSHong Zhang       if (sctx.newshift) break;    /* sctx.shift_amount is updated */
9470e95ead3SHong Zhang       dk = sctx.pv;
94878910aadSHong Zhang 
94978910aadSHong Zhang       /* copy data into U(k,:) */
95078910aadSHong Zhang       ba[bi[k]] = 1.0/dk;
95178910aadSHong Zhang       jmin      = bi[k]+1;
95278910aadSHong Zhang       nz        = bi[k+1] - jmin;
95378910aadSHong Zhang       if (nz) {
95478910aadSHong Zhang         bcol = bj + jmin;
95578910aadSHong Zhang         bval = ba + jmin;
95678910aadSHong Zhang         while (nz--) {
95778910aadSHong Zhang           *bval++       = rtmp[*bcol];
95878910aadSHong Zhang           rtmp[*bcol++] = 0.0;
95978910aadSHong Zhang         }
96078910aadSHong Zhang         /* add k-th row into il and jl */
96178910aadSHong Zhang         il[k] = jmin;
96278910aadSHong Zhang         i     = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
96378910aadSHong Zhang       }
96478910aadSHong Zhang     }
96507b50cabSHong Zhang   } while (sctx.newshift);
966fca92195SBarry Smith   ierr = PetscFree3(rtmp,il,jl);CHKERRQ(ierr);
96778910aadSHong Zhang 
9680a3c351aSHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
9690a3c351aSHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
97078910aadSHong Zhang   C->assembled           = PETSC_TRUE;
97178910aadSHong Zhang   C->preallocated        = PETSC_TRUE;
97226fbe8dcSKarl Rupp 
973d0f46423SBarry Smith   ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr);
9743cea4cbeSHong Zhang   if (sctx.nshift) {
975f4db908eSBarry Smith     if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
97657622a8eSBarry Smith       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr);
977f4db908eSBarry Smith     } else if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
97857622a8eSBarry Smith       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr);
97978910aadSHong Zhang     }
98078910aadSHong Zhang   }
981c05c3958SHong Zhang   PetscFunctionReturn(0);
982c05c3958SHong Zhang }
983c05c3958SHong Zhang 
984c6db04a5SJed Brown #include <petscbt.h>
985c6db04a5SJed Brown #include <../src/mat/utils/freespace.h>
986c05c3958SHong Zhang #undef __FUNCT__
987c05c3958SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqBAIJ"
9880481f469SBarry Smith PetscErrorCode MatICCFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
989c05c3958SHong Zhang {
99078910aadSHong Zhang   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data;
99178910aadSHong Zhang   Mat_SeqSBAIJ       *b;
99278910aadSHong Zhang   Mat                B;
99378910aadSHong Zhang   PetscErrorCode     ierr;
99448dd3d27SHong Zhang   PetscBool          perm_identity,missing;
9955d0c19d7SBarry Smith   PetscInt           reallocs=0,i,*ai=a->i,*aj=a->j,am=a->mbs,bs=A->rmap->bs,*ui;
9965d0c19d7SBarry Smith   const PetscInt     *rip;
99778910aadSHong Zhang   PetscInt           jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow;
9980298fd71SBarry Smith   PetscInt           nlnk,*lnk,*lnk_lvl=NULL,ncols,ncols_upper,*cols,*cols_lvl,*uj,**uj_ptr,**uj_lvl_ptr;
99978910aadSHong Zhang   PetscReal          fill          =info->fill,levels=info->levels;
10000298fd71SBarry Smith   PetscFreeSpaceList free_space    =NULL,current_space=NULL;
10010298fd71SBarry Smith   PetscFreeSpaceList free_space_lvl=NULL,current_space_lvl=NULL;
100278910aadSHong Zhang   PetscBT            lnkbt;
100378910aadSHong Zhang 
1004c05c3958SHong Zhang   PetscFunctionBegin;
100548dd3d27SHong Zhang   ierr = MatMissingDiagonal(A,&missing,&i);CHKERRQ(ierr);
100648dd3d27SHong Zhang   if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i);
100748dd3d27SHong Zhang 
10086ad2eaddSHong Zhang   if (bs > 1) {
10096ad2eaddSHong Zhang     if (!a->sbaijMat) {
1010ceb03754SKris Buschelman       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
10116ad2eaddSHong Zhang     }
1012719d5645SBarry Smith     (fact)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqSBAIJ;  /* undue the change made in MatGetFactor_seqbaij_petsc */
101326fbe8dcSKarl Rupp 
1014719d5645SBarry Smith     ierr = MatICCFactorSymbolic(fact,a->sbaijMat,perm,info);CHKERRQ(ierr);
10156ad2eaddSHong Zhang     PetscFunctionReturn(0);
10166ad2eaddSHong Zhang   }
10176ad2eaddSHong Zhang 
101878910aadSHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
101978910aadSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
102078910aadSHong Zhang 
102178910aadSHong Zhang   /* special case that simply copies fill pattern */
102278910aadSHong Zhang   if (!levels && perm_identity) {
1023854ce69bSBarry Smith     ierr = PetscMalloc1(am+1,&ui);CHKERRQ(ierr);
102426fbe8dcSKarl Rupp     for (i=0; i<am; i++) ui[i] = ai[i+1] - a->diag[i]; /* ui: rowlengths - changes when !perm_identity */
1025719d5645SBarry Smith     B    = fact;
102678910aadSHong Zhang     ierr = MatSeqSBAIJSetPreallocation(B,1,0,ui);CHKERRQ(ierr);
102778910aadSHong Zhang 
1028b24902e0SBarry Smith 
102978910aadSHong Zhang     b  = (Mat_SeqSBAIJ*)B->data;
103078910aadSHong Zhang     uj = b->j;
103178910aadSHong Zhang     for (i=0; i<am; i++) {
103278910aadSHong Zhang       aj = a->j + a->diag[i];
103326fbe8dcSKarl Rupp       for (j=0; j<ui[i]; j++) *uj++ = *aj++;
103478910aadSHong Zhang       b->ilen[i] = ui[i];
103578910aadSHong Zhang     }
103678910aadSHong Zhang     ierr = PetscFree(ui);CHKERRQ(ierr);
103726fbe8dcSKarl Rupp 
1038d5f3da31SBarry Smith     B->factortype = MAT_FACTOR_NONE;
103926fbe8dcSKarl Rupp 
104078910aadSHong Zhang     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
104178910aadSHong Zhang     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1042d5f3da31SBarry Smith     B->factortype = MAT_FACTOR_ICC;
104378910aadSHong Zhang 
104478910aadSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
104578910aadSHong Zhang     PetscFunctionReturn(0);
104678910aadSHong Zhang   }
104778910aadSHong Zhang 
104878910aadSHong Zhang   /* initialization */
1049854ce69bSBarry Smith   ierr  = PetscMalloc1(am+1,&ui);CHKERRQ(ierr);
1050e60cf9a0SBarry Smith   ui[0] = 0;
1051854ce69bSBarry Smith   ierr  = PetscMalloc1(2*am+1,&cols_lvl);CHKERRQ(ierr);
105278910aadSHong Zhang 
105378910aadSHong Zhang   /* jl: linked list for storing indices of the pivot rows
105478910aadSHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
1055dcca6d9dSJed Brown   ierr = PetscMalloc4(am,&uj_ptr,am,&uj_lvl_ptr,am,&il,am,&jl);CHKERRQ(ierr);
105678910aadSHong Zhang   for (i=0; i<am; i++) {
1057e60cf9a0SBarry Smith     jl[i] = am; il[i] = 0;
105878910aadSHong Zhang   }
105978910aadSHong Zhang 
106078910aadSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
106178910aadSHong Zhang   nlnk = am + 1;
106278910aadSHong Zhang   ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
106378910aadSHong Zhang 
106495b5ac22SHong Zhang   /* initial FreeSpace size is fill*(ai[am]+am)/2 */
1065f91af8c7SBarry Smith   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am]/2,am/2)),&free_space);CHKERRQ(ierr);
106626fbe8dcSKarl Rupp 
106778910aadSHong Zhang   current_space = free_space;
106826fbe8dcSKarl Rupp 
1069f91af8c7SBarry Smith   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am]/2,am/2)),&free_space_lvl);CHKERRQ(ierr);
107078910aadSHong Zhang   current_space_lvl = free_space_lvl;
107178910aadSHong Zhang 
107278910aadSHong Zhang   for (k=0; k<am; k++) {  /* for each active row k */
107378910aadSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
107478910aadSHong Zhang     nzk         = 0;
107578910aadSHong Zhang     ncols       = ai[rip[k]+1] - ai[rip[k]];
107678910aadSHong Zhang     ncols_upper = 0;
107778910aadSHong Zhang     cols        = cols_lvl + am;
107878910aadSHong Zhang     for (j=0; j<ncols; j++) {
107978910aadSHong Zhang       i = rip[*(aj + ai[rip[k]] + j)];
108078910aadSHong Zhang       if (i >= k) { /* only take upper triangular entry */
108178910aadSHong Zhang         cols[ncols_upper]     = i;
108278910aadSHong Zhang         cols_lvl[ncols_upper] = -1;  /* initialize level for nonzero entries */
108378910aadSHong Zhang         ncols_upper++;
108478910aadSHong Zhang       }
108578910aadSHong Zhang     }
108678910aadSHong Zhang     ierr = PetscIncompleteLLAdd(ncols_upper,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
108778910aadSHong Zhang     nzk += nlnk;
108878910aadSHong Zhang 
108978910aadSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
109078910aadSHong Zhang     prow = jl[k]; /* 1st pivot row */
109178910aadSHong Zhang 
109278910aadSHong Zhang     while (prow < k) {
109378910aadSHong Zhang       nextprow = jl[prow];
109478910aadSHong Zhang 
109578910aadSHong Zhang       /* merge prow into k-th row */
109678910aadSHong Zhang       jmin  = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:am-1) */
109778910aadSHong Zhang       jmax  = ui[prow+1];
109878910aadSHong Zhang       ncols = jmax-jmin;
109978910aadSHong Zhang       i     = jmin - ui[prow];
110078910aadSHong Zhang       cols  = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
110178910aadSHong Zhang       for (j=0; j<ncols; j++) cols_lvl[j] = *(uj_lvl_ptr[prow] + i + j);
11025a8e39fbSHong Zhang       ierr = PetscIncompleteLLAddSorted(ncols,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
110378910aadSHong Zhang       nzk += nlnk;
110478910aadSHong Zhang 
110578910aadSHong Zhang       /* update il and jl for prow */
110678910aadSHong Zhang       if (jmin < jmax) {
110778910aadSHong Zhang         il[prow] = jmin;
110826fbe8dcSKarl Rupp 
110978910aadSHong Zhang         j = *cols; jl[prow] = jl[j]; jl[j] = prow;
111078910aadSHong Zhang       }
111178910aadSHong Zhang       prow = nextprow;
111278910aadSHong Zhang     }
111378910aadSHong Zhang 
111478910aadSHong Zhang     /* if free space is not available, make more free space */
111578910aadSHong Zhang     if (current_space->local_remaining<nzk) {
111678910aadSHong Zhang       i    = am - k + 1; /* num of unfactored rows */
1117f91af8c7SBarry Smith       i    = PetscMin(PetscIntMultTruncate(i,nzk), PetscIntMultTruncate(i,i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
1118a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
1119a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(i,&current_space_lvl);CHKERRQ(ierr);
112078910aadSHong Zhang       reallocs++;
112178910aadSHong Zhang     }
112278910aadSHong Zhang 
112378910aadSHong Zhang     /* copy data into free_space and free_space_lvl, then initialize lnk */
112478910aadSHong Zhang     ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
112578910aadSHong Zhang 
112678910aadSHong Zhang     /* add the k-th row into il and jl */
112778910aadSHong Zhang     if (nzk-1 > 0) {
112878910aadSHong Zhang       i     = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
112978910aadSHong Zhang       jl[k] = jl[i]; jl[i] = k;
113078910aadSHong Zhang       il[k] = ui[k] + 1;
113178910aadSHong Zhang     }
113278910aadSHong Zhang     uj_ptr[k]     = current_space->array;
113378910aadSHong Zhang     uj_lvl_ptr[k] = current_space_lvl->array;
113478910aadSHong Zhang 
113578910aadSHong Zhang     current_space->array           += nzk;
113678910aadSHong Zhang     current_space->local_used      += nzk;
113778910aadSHong Zhang     current_space->local_remaining -= nzk;
113878910aadSHong Zhang 
113978910aadSHong Zhang     current_space_lvl->array           += nzk;
114078910aadSHong Zhang     current_space_lvl->local_used      += nzk;
114178910aadSHong Zhang     current_space_lvl->local_remaining -= nzk;
114278910aadSHong Zhang 
114378910aadSHong Zhang     ui[k+1] = ui[k] + nzk;
114478910aadSHong Zhang   }
114578910aadSHong Zhang 
114678910aadSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
1147fca92195SBarry Smith   ierr = PetscFree4(uj_ptr,uj_lvl_ptr,il,jl);CHKERRQ(ierr);
114878910aadSHong Zhang   ierr = PetscFree(cols_lvl);CHKERRQ(ierr);
114978910aadSHong Zhang 
11509263d837SHong Zhang   /* copy free_space into uj and free free_space; set uj in new datastructure; */
1151854ce69bSBarry Smith   ierr = PetscMalloc1(ui[am]+1,&uj);CHKERRQ(ierr);
1152a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
115378910aadSHong Zhang   ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1154a1a86e44SBarry Smith   ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
115578910aadSHong Zhang 
115678910aadSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
1157719d5645SBarry Smith   B    = fact;
11580298fd71SBarry Smith   ierr = MatSeqSBAIJSetPreallocation(B,1,MAT_SKIP_ALLOCATION,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;
116426fbe8dcSKarl Rupp 
1165854ce69bSBarry Smith   ierr = PetscMalloc1(ui[am]+1,&b->a);CHKERRQ(ierr);
116626fbe8dcSKarl Rupp 
116778910aadSHong Zhang   b->j             = uj;
116878910aadSHong Zhang   b->i             = ui;
116978910aadSHong Zhang   b->diag          = 0;
117078910aadSHong Zhang   b->ilen          = 0;
117178910aadSHong Zhang   b->imax          = 0;
117278910aadSHong Zhang   b->row           = perm;
117378910aadSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
117426fbe8dcSKarl Rupp 
117578910aadSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
117626fbe8dcSKarl Rupp 
117778910aadSHong Zhang   b->icol = perm;
117826fbe8dcSKarl Rupp 
117978910aadSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
1180854ce69bSBarry Smith   ierr    = PetscMalloc1(am+1,&b->solve_work);CHKERRQ(ierr);
11813bb1ff40SBarry Smith   ierr    = PetscLogObjectMemory((PetscObject)B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
118226fbe8dcSKarl Rupp 
118378910aadSHong Zhang   b->maxnz = b->nz = ui[am];
118478910aadSHong Zhang 
118578910aadSHong Zhang   B->info.factor_mallocs   = reallocs;
118678910aadSHong Zhang   B->info.fill_ratio_given = fill;
118775567043SBarry Smith   if (ai[am] != 0.) {
118895b5ac22SHong Zhang     /* nonzeros in lower triangular part of A (includign diagonals)= (ai[am]+am)/2 */
118995b5ac22SHong Zhang     B->info.fill_ratio_needed = ((PetscReal)2*ui[am])/(ai[am]+am);
119078910aadSHong Zhang   } else {
119178910aadSHong Zhang     B->info.fill_ratio_needed = 0.0;
119278910aadSHong Zhang   }
11939263d837SHong Zhang #if defined(PETSC_USE_INFO)
11949263d837SHong Zhang   if (ai[am] != 0) {
119595b5ac22SHong Zhang     PetscReal af = B->info.fill_ratio_needed;
119657622a8eSBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)fill,(double)af);CHKERRQ(ierr);
119757622a8eSBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);CHKERRQ(ierr);
119857622a8eSBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%g) for best performance.\n",(double)af);CHKERRQ(ierr);
11999263d837SHong Zhang   } else {
12009263d837SHong Zhang     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
12019263d837SHong Zhang   }
12029263d837SHong Zhang #endif
120378910aadSHong Zhang   if (perm_identity) {
12040a3c351aSHong Zhang     B->ops->solve                 = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
12050a3c351aSHong Zhang     B->ops->solvetranspose        = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
120678910aadSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
120778910aadSHong Zhang   } else {
1208719d5645SBarry Smith     (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
120978910aadSHong Zhang   }
1210c05c3958SHong Zhang   PetscFunctionReturn(0);
1211c05c3958SHong Zhang }
1212c05c3958SHong Zhang 
1213c05c3958SHong Zhang #undef __FUNCT__
1214c05c3958SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqBAIJ"
12150481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
1216c05c3958SHong Zhang {
121778910aadSHong Zhang   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data;
121878910aadSHong Zhang   Mat_SeqSBAIJ       *b;
121978910aadSHong Zhang   Mat                B;
122078910aadSHong Zhang   PetscErrorCode     ierr;
12219186b105SHong Zhang   PetscBool          perm_identity,missing;
122278910aadSHong Zhang   PetscReal          fill = info->fill;
12235d0c19d7SBarry Smith   const PetscInt     *rip;
12245d0c19d7SBarry Smith   PetscInt           i,mbs=a->mbs,bs=A->rmap->bs,*ai=a->i,*aj=a->j,reallocs=0,prow;
122578910aadSHong Zhang   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
122678910aadSHong Zhang   PetscInt           nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr;
12270298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
122878910aadSHong Zhang   PetscBT            lnkbt;
122978910aadSHong Zhang 
1230c05c3958SHong Zhang   PetscFunctionBegin;
12316ad2eaddSHong Zhang   if (bs > 1) { /* convert to seqsbaij */
12326ad2eaddSHong Zhang     if (!a->sbaijMat) {
1233ceb03754SKris Buschelman       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
12346ad2eaddSHong Zhang     }
1235719d5645SBarry Smith     (fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqSBAIJ; /* undue the change made in MatGetFactor_seqbaij_petsc */
123626fbe8dcSKarl Rupp 
1237719d5645SBarry Smith     ierr = MatCholeskyFactorSymbolic(fact,a->sbaijMat,perm,info);CHKERRQ(ierr);
12386ad2eaddSHong Zhang     PetscFunctionReturn(0);
12396ad2eaddSHong Zhang   }
12406ad2eaddSHong Zhang 
12419186b105SHong Zhang   ierr = MatMissingDiagonal(A,&missing,&i);CHKERRQ(ierr);
12429186b105SHong Zhang   if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i);
12439186b105SHong Zhang 
124478910aadSHong Zhang   /* check whether perm is the identity mapping */
124578910aadSHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
1246e32f2f54SBarry Smith   if (!perm_identity) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix reordering is not supported");
124778910aadSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
124878910aadSHong Zhang 
124978910aadSHong Zhang   /* initialization */
1250854ce69bSBarry Smith   ierr  = PetscMalloc1(mbs+1,&ui);CHKERRQ(ierr);
1251e60cf9a0SBarry Smith   ui[0] = 0;
125278910aadSHong Zhang 
125378910aadSHong Zhang   /* jl: linked list for storing indices of the pivot rows
125478910aadSHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
1255dcca6d9dSJed Brown   ierr = PetscMalloc4(mbs,&ui_ptr,mbs,&il,mbs,&jl,mbs,&cols);CHKERRQ(ierr);
125678910aadSHong Zhang   for (i=0; i<mbs; i++) {
1257e60cf9a0SBarry Smith     jl[i] = mbs; il[i] = 0;
125878910aadSHong Zhang   }
125978910aadSHong Zhang 
126078910aadSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
126178910aadSHong Zhang   nlnk = mbs + 1;
126278910aadSHong Zhang   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
126378910aadSHong Zhang 
12646a69fef8SHong Zhang   /* initial FreeSpace size is fill* (ai[mbs]+mbs)/2 */
1265f91af8c7SBarry Smith   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[mbs]/2,mbs/2)),&free_space);CHKERRQ(ierr);
126626fbe8dcSKarl Rupp 
126778910aadSHong Zhang   current_space = free_space;
126878910aadSHong Zhang 
126978910aadSHong Zhang   for (k=0; k<mbs; k++) {  /* for each active row k */
127078910aadSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
127178910aadSHong Zhang     nzk         = 0;
127278910aadSHong Zhang     ncols       = ai[rip[k]+1] - ai[rip[k]];
1273e60cf9a0SBarry Smith     ncols_upper = 0;
127478910aadSHong Zhang     for (j=0; j<ncols; j++) {
127578910aadSHong Zhang       i = rip[*(aj + ai[rip[k]] + j)];
127678910aadSHong Zhang       if (i >= k) { /* only take upper triangular entry */
127778910aadSHong Zhang         cols[ncols_upper] = i;
127878910aadSHong Zhang         ncols_upper++;
127978910aadSHong Zhang       }
128078910aadSHong Zhang     }
128178910aadSHong Zhang     ierr = PetscLLAdd(ncols_upper,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
128278910aadSHong Zhang     nzk += nlnk;
128378910aadSHong Zhang 
128478910aadSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
128578910aadSHong Zhang     prow = jl[k]; /* 1st pivot row */
128678910aadSHong Zhang 
128778910aadSHong Zhang     while (prow < k) {
128878910aadSHong Zhang       nextprow = jl[prow];
128978910aadSHong Zhang       /* merge prow into k-th row */
129078910aadSHong Zhang       jmin   = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
129178910aadSHong Zhang       jmax   = ui[prow+1];
129278910aadSHong Zhang       ncols  = jmax-jmin;
129378910aadSHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
12945a8e39fbSHong Zhang       ierr   = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
129578910aadSHong Zhang       nzk   += nlnk;
129678910aadSHong Zhang 
129778910aadSHong Zhang       /* update il and jl for prow */
129878910aadSHong Zhang       if (jmin < jmax) {
129978910aadSHong Zhang         il[prow] = jmin;
130026fbe8dcSKarl Rupp         j        = *uj_ptr;
130126fbe8dcSKarl Rupp         jl[prow] = jl[j];
130226fbe8dcSKarl Rupp         jl[j]    = prow;
130378910aadSHong Zhang       }
130478910aadSHong Zhang       prow = nextprow;
130578910aadSHong Zhang     }
130678910aadSHong Zhang 
130778910aadSHong Zhang     /* if free space is not available, make more free space */
130878910aadSHong Zhang     if (current_space->local_remaining<nzk) {
130978910aadSHong Zhang       i    = mbs - k + 1; /* num of unfactored rows */
1310f91af8c7SBarry Smith       i    = PetscMin(PetscIntMultTruncate(i,nzk), PetscIntMultTruncate(i,i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
1311a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
131278910aadSHong Zhang       reallocs++;
131378910aadSHong Zhang     }
131478910aadSHong Zhang 
131578910aadSHong Zhang     /* copy data into free space, then initialize lnk */
131678910aadSHong Zhang     ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
131778910aadSHong Zhang 
131878910aadSHong Zhang     /* add the k-th row into il and jl */
131978910aadSHong Zhang     if (nzk-1 > 0) {
132078910aadSHong Zhang       i     = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
132178910aadSHong Zhang       jl[k] = jl[i]; jl[i] = k;
132278910aadSHong Zhang       il[k] = ui[k] + 1;
132378910aadSHong Zhang     }
132478910aadSHong Zhang     ui_ptr[k]                       = current_space->array;
132578910aadSHong Zhang     current_space->array           += nzk;
132678910aadSHong Zhang     current_space->local_used      += nzk;
132778910aadSHong Zhang     current_space->local_remaining -= nzk;
132878910aadSHong Zhang 
132978910aadSHong Zhang     ui[k+1] = ui[k] + nzk;
133078910aadSHong Zhang   }
133178910aadSHong Zhang 
133278910aadSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
1333fca92195SBarry Smith   ierr = PetscFree4(ui_ptr,il,jl,cols);CHKERRQ(ierr);
133478910aadSHong Zhang 
13359263d837SHong Zhang   /* copy free_space into uj and free free_space; set uj in new datastructure; */
1336854ce69bSBarry Smith   ierr = PetscMalloc1(ui[mbs]+1,&uj);CHKERRQ(ierr);
1337a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
133878910aadSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
133978910aadSHong Zhang 
134078910aadSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
1341719d5645SBarry Smith   B    = fact;
13420298fd71SBarry Smith   ierr = MatSeqSBAIJSetPreallocation(B,bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr);
134378910aadSHong Zhang 
134478910aadSHong Zhang   b               = (Mat_SeqSBAIJ*)B->data;
134578910aadSHong Zhang   b->singlemalloc = PETSC_FALSE;
1346e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
1347e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
134826fbe8dcSKarl Rupp 
1349854ce69bSBarry Smith   ierr = PetscMalloc1(ui[mbs]+1,&b->a);CHKERRQ(ierr);
135026fbe8dcSKarl Rupp 
135178910aadSHong Zhang   b->j             = uj;
135278910aadSHong Zhang   b->i             = ui;
135378910aadSHong Zhang   b->diag          = 0;
135478910aadSHong Zhang   b->ilen          = 0;
135578910aadSHong Zhang   b->imax          = 0;
135678910aadSHong Zhang   b->row           = perm;
135778910aadSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
135826fbe8dcSKarl Rupp 
135978910aadSHong Zhang   ierr     = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
136078910aadSHong Zhang   b->icol  = perm;
136178910aadSHong Zhang   ierr     = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
1362854ce69bSBarry Smith   ierr     = PetscMalloc1(mbs+1,&b->solve_work);CHKERRQ(ierr);
13633bb1ff40SBarry Smith   ierr     = PetscLogObjectMemory((PetscObject)B,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
136478910aadSHong Zhang   b->maxnz = b->nz = ui[mbs];
136578910aadSHong Zhang 
136678910aadSHong Zhang   B->info.factor_mallocs   = reallocs;
136778910aadSHong Zhang   B->info.fill_ratio_given = fill;
136875567043SBarry Smith   if (ai[mbs] != 0.) {
136995b5ac22SHong Zhang     /* nonzeros in lower triangular part of A = (ai[mbs]+mbs)/2 */
137095b5ac22SHong Zhang     B->info.fill_ratio_needed = ((PetscReal)2*ui[mbs])/(ai[mbs]+mbs);
137178910aadSHong Zhang   } else {
137278910aadSHong Zhang     B->info.fill_ratio_needed = 0.0;
137378910aadSHong Zhang   }
13749263d837SHong Zhang #if defined(PETSC_USE_INFO)
13759263d837SHong Zhang   if (ai[mbs] != 0.) {
13769263d837SHong Zhang     PetscReal af = B->info.fill_ratio_needed;
137757622a8eSBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)fill,(double)af);CHKERRQ(ierr);
137857622a8eSBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);CHKERRQ(ierr);
137957622a8eSBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%g) for best performance.\n",(double)af);CHKERRQ(ierr);
13809263d837SHong Zhang   } else {
13819263d837SHong Zhang     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
13829263d837SHong Zhang   }
13839263d837SHong Zhang #endif
138478910aadSHong Zhang   if (perm_identity) {
13856ad2eaddSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
138678910aadSHong Zhang   } else {
13876ad2eaddSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
138878910aadSHong Zhang   }
1389c05c3958SHong Zhang   PetscFunctionReturn(0);
1390c05c3958SHong Zhang }
1391c8342467SHong Zhang 
13922efa7f71SHong Zhang #undef __FUNCT__
13934dd39f65SShri Abhyankar #define __FUNCT__ "MatSolve_SeqBAIJ_N_NaturalOrdering"
13944dd39f65SShri Abhyankar PetscErrorCode MatSolve_SeqBAIJ_N_NaturalOrdering(Mat A,Vec bb,Vec xx)
13951a83e813SShri Abhyankar {
13961a83e813SShri Abhyankar   Mat_SeqBAIJ       *a=(Mat_SeqBAIJ*)A->data;
13971a83e813SShri Abhyankar   PetscErrorCode    ierr;
13981a83e813SShri Abhyankar   const PetscInt    *ai=a->i,*aj=a->j,*adiag=a->diag,*vi;
13991a83e813SShri Abhyankar   PetscInt          i,k,n=a->mbs;
14001a83e813SShri Abhyankar   PetscInt          nz,bs=A->rmap->bs,bs2=a->bs2;
1401d9ca1df4SBarry Smith   const MatScalar   *aa=a->a,*v;
1402d9ca1df4SBarry Smith   PetscScalar       *x,*s,*t,*ls;
1403d9ca1df4SBarry Smith   const PetscScalar *b;
14041a83e813SShri Abhyankar 
14051a83e813SShri Abhyankar   PetscFunctionBegin;
1406d9ca1df4SBarry Smith   ierr = VecGetArrayRead(bb,&b);CHKERRQ(ierr);
14071a83e813SShri Abhyankar   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
14081a83e813SShri Abhyankar   t    = a->solve_work;
14091a83e813SShri Abhyankar 
14101a83e813SShri Abhyankar   /* forward solve the lower triangular */
14111a83e813SShri Abhyankar   ierr = PetscMemcpy(t,b,bs*sizeof(PetscScalar));CHKERRQ(ierr); /* copy 1st block of b to t */
14121a83e813SShri Abhyankar 
14131a83e813SShri Abhyankar   for (i=1; i<n; i++) {
14141a83e813SShri Abhyankar     v    = aa + bs2*ai[i];
14151a83e813SShri Abhyankar     vi   = aj + ai[i];
14161a83e813SShri Abhyankar     nz   = ai[i+1] - ai[i];
14171a83e813SShri Abhyankar     s    = t + bs*i;
14181a83e813SShri Abhyankar     ierr = PetscMemcpy(s,b+bs*i,bs*sizeof(PetscScalar));CHKERRQ(ierr); /* copy i_th block of b to t */
14191a83e813SShri Abhyankar     for (k=0;k<nz;k++) {
142096b95a6bSBarry Smith       PetscKernel_v_gets_v_minus_A_times_w(bs,s,v,t+bs*vi[k]);
14211a83e813SShri Abhyankar       v += bs2;
14221a83e813SShri Abhyankar     }
14231a83e813SShri Abhyankar   }
14241a83e813SShri Abhyankar 
14251a83e813SShri Abhyankar   /* backward solve the upper triangular */
14261a83e813SShri Abhyankar   ls = a->solve_work + A->cmap->n;
14271a83e813SShri Abhyankar   for (i=n-1; i>=0; i--) {
14281a83e813SShri Abhyankar     v    = aa + bs2*(adiag[i+1]+1);
14291a83e813SShri Abhyankar     vi   = aj + adiag[i+1]+1;
14301a83e813SShri Abhyankar     nz   = adiag[i] - adiag[i+1]-1;
14311a83e813SShri Abhyankar     ierr = PetscMemcpy(ls,t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr);
14321a83e813SShri Abhyankar     for (k=0; k<nz; k++) {
143396b95a6bSBarry Smith       PetscKernel_v_gets_v_minus_A_times_w(bs,ls,v,t+bs*vi[k]);
14341a83e813SShri Abhyankar       v += bs2;
14351a83e813SShri Abhyankar     }
143696b95a6bSBarry Smith     PetscKernel_w_gets_A_times_v(bs,ls,aa+bs2*adiag[i],t+i*bs); /* *inv(diagonal[i]) */
14371a83e813SShri Abhyankar     ierr = PetscMemcpy(x+i*bs,t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr);
14381a83e813SShri Abhyankar   }
14391a83e813SShri Abhyankar 
1440d9ca1df4SBarry Smith   ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr);
14411a83e813SShri Abhyankar   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
14421a83e813SShri Abhyankar   ierr = PetscLogFlops(2.0*(a->bs2)*(a->nz) - A->rmap->bs*A->cmap->n);CHKERRQ(ierr);
14431a83e813SShri Abhyankar   PetscFunctionReturn(0);
14441a83e813SShri Abhyankar }
14451a83e813SShri Abhyankar 
14462efa7f71SHong Zhang #undef __FUNCT__
14474dd39f65SShri Abhyankar #define __FUNCT__ "MatSolve_SeqBAIJ_N"
14484dd39f65SShri Abhyankar PetscErrorCode MatSolve_SeqBAIJ_N(Mat A,Vec bb,Vec xx)
144935aa4fcfSShri Abhyankar {
145035aa4fcfSShri Abhyankar   Mat_SeqBAIJ        *a   =(Mat_SeqBAIJ*)A->data;
145135aa4fcfSShri Abhyankar   IS                 iscol=a->col,isrow=a->row;
145235aa4fcfSShri Abhyankar   PetscErrorCode     ierr;
145335aa4fcfSShri Abhyankar   const PetscInt     *r,*c,*rout,*cout,*ai=a->i,*aj=a->j,*adiag=a->diag,*vi;
145435aa4fcfSShri Abhyankar   PetscInt           i,m,n=a->mbs;
145535aa4fcfSShri Abhyankar   PetscInt           nz,bs=A->rmap->bs,bs2=a->bs2;
1456d9ca1df4SBarry Smith   const MatScalar    *aa=a->a,*v;
1457d9ca1df4SBarry Smith   PetscScalar        *x,*s,*t,*ls;
1458d9ca1df4SBarry Smith   const PetscScalar  *b;
145935aa4fcfSShri Abhyankar 
146035aa4fcfSShri Abhyankar   PetscFunctionBegin;
1461d9ca1df4SBarry Smith   ierr = VecGetArrayRead(bb,&b);CHKERRQ(ierr);
146235aa4fcfSShri Abhyankar   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
146335aa4fcfSShri Abhyankar   t    = a->solve_work;
146435aa4fcfSShri Abhyankar 
146535aa4fcfSShri Abhyankar   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
146635aa4fcfSShri Abhyankar   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
146735aa4fcfSShri Abhyankar 
146835aa4fcfSShri Abhyankar   /* forward solve the lower triangular */
146935aa4fcfSShri Abhyankar   ierr = PetscMemcpy(t,b+bs*r[0],bs*sizeof(PetscScalar));CHKERRQ(ierr);
147035aa4fcfSShri Abhyankar   for (i=1; i<n; i++) {
147135aa4fcfSShri Abhyankar     v    = aa + bs2*ai[i];
147235aa4fcfSShri Abhyankar     vi   = aj + ai[i];
147335aa4fcfSShri Abhyankar     nz   = ai[i+1] - ai[i];
147435aa4fcfSShri Abhyankar     s    = t + bs*i;
147535aa4fcfSShri Abhyankar     ierr = PetscMemcpy(s,b+bs*r[i],bs*sizeof(PetscScalar));CHKERRQ(ierr);
147635aa4fcfSShri Abhyankar     for (m=0; m<nz; m++) {
147796b95a6bSBarry Smith       PetscKernel_v_gets_v_minus_A_times_w(bs,s,v,t+bs*vi[m]);
147835aa4fcfSShri Abhyankar       v += bs2;
147935aa4fcfSShri Abhyankar     }
148035aa4fcfSShri Abhyankar   }
148135aa4fcfSShri Abhyankar 
148235aa4fcfSShri Abhyankar   /* backward solve the upper triangular */
148335aa4fcfSShri Abhyankar   ls = a->solve_work + A->cmap->n;
148435aa4fcfSShri Abhyankar   for (i=n-1; i>=0; i--) {
148535aa4fcfSShri Abhyankar     v    = aa + bs2*(adiag[i+1]+1);
148635aa4fcfSShri Abhyankar     vi   = aj + adiag[i+1]+1;
148735aa4fcfSShri Abhyankar     nz   = adiag[i] - adiag[i+1] - 1;
148835aa4fcfSShri Abhyankar     ierr = PetscMemcpy(ls,t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr);
148935aa4fcfSShri Abhyankar     for (m=0; m<nz; m++) {
149096b95a6bSBarry Smith       PetscKernel_v_gets_v_minus_A_times_w(bs,ls,v,t+bs*vi[m]);
149135aa4fcfSShri Abhyankar       v += bs2;
149235aa4fcfSShri Abhyankar     }
149396b95a6bSBarry Smith     PetscKernel_w_gets_A_times_v(bs,ls,v,t+i*bs); /* *inv(diagonal[i]) */
149435aa4fcfSShri Abhyankar     ierr = PetscMemcpy(x + bs*c[i],t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr);
149535aa4fcfSShri Abhyankar   }
149635aa4fcfSShri Abhyankar   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
149735aa4fcfSShri Abhyankar   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
1498d9ca1df4SBarry Smith   ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr);
149935aa4fcfSShri Abhyankar   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
150035aa4fcfSShri Abhyankar   ierr = PetscLogFlops(2.0*(a->bs2)*(a->nz) - A->rmap->bs*A->cmap->n);CHKERRQ(ierr);
150135aa4fcfSShri Abhyankar   PetscFunctionReturn(0);
150235aa4fcfSShri Abhyankar }
150335aa4fcfSShri Abhyankar 
150435aa4fcfSShri Abhyankar #undef __FUNCT__
1505a25532f0SBarry Smith #define __FUNCT__ "MatBlockAbs_privat"
1506a25532f0SBarry Smith /*
1507a25532f0SBarry Smith     For each block in an block array saves the largest absolute value in the block into another array
1508a25532f0SBarry Smith */
1509a25532f0SBarry Smith static PetscErrorCode MatBlockAbs_private(PetscInt nbs,PetscInt bs2,PetscScalar *blockarray,PetscReal *absarray)
15102efa7f71SHong Zhang {
15112efa7f71SHong Zhang   PetscErrorCode ierr;
15122efa7f71SHong Zhang   PetscInt       i,j;
15135fd66863SKarl Rupp 
15142efa7f71SHong Zhang   PetscFunctionBegin;
15152efa7f71SHong Zhang   ierr = PetscMemzero(absarray,(nbs+1)*sizeof(PetscReal));CHKERRQ(ierr);
15162efa7f71SHong Zhang   for (i=0; i<nbs; i++) {
15172efa7f71SHong Zhang     for (j=0; j<bs2; j++) {
15182efa7f71SHong Zhang       if (absarray[i] < PetscAbsScalar(blockarray[i*nbs+j])) absarray[i] = PetscAbsScalar(blockarray[i*nbs+j]);
15192efa7f71SHong Zhang     }
15202efa7f71SHong Zhang   }
15212efa7f71SHong Zhang   PetscFunctionReturn(0);
15222efa7f71SHong Zhang }
15232efa7f71SHong Zhang 
1524c8342467SHong Zhang #undef __FUNCT__
1525c8342467SHong Zhang #define __FUNCT__ "MatILUDTFactor_SeqBAIJ"
1526fe97e370SBarry Smith /*
1527fe97e370SBarry Smith      This needs to be renamed and called by the regular MatILUFactor_SeqBAIJ when drop tolerance is used
1528fe97e370SBarry Smith */
1529c8342467SHong Zhang PetscErrorCode MatILUDTFactor_SeqBAIJ(Mat A,IS isrow,IS iscol,const MatFactorInfo *info,Mat *fact)
1530c8342467SHong Zhang {
15312efa7f71SHong Zhang   Mat            B = *fact;
15322efa7f71SHong Zhang   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data,*b;
15332efa7f71SHong Zhang   IS             isicol;
15342efa7f71SHong Zhang   PetscErrorCode ierr;
15352efa7f71SHong Zhang   const PetscInt *r,*ic;
15362efa7f71SHong Zhang   PetscInt       i,mbs=a->mbs,bs=A->rmap->bs,bs2=a->bs2,*ai=a->i,*aj=a->j,*ajtmp,*adiag;
15372efa7f71SHong Zhang   PetscInt       *bi,*bj,*bdiag;
15382efa7f71SHong Zhang 
15392efa7f71SHong Zhang   PetscInt  row,nzi,nzi_bl,nzi_bu,*im,dtcount,nzi_al,nzi_au;
15402efa7f71SHong Zhang   PetscInt  nlnk,*lnk;
15412efa7f71SHong Zhang   PetscBT   lnkbt;
1542ace3abfcSBarry Smith   PetscBool row_identity,icol_identity;
15432efa7f71SHong Zhang   MatScalar *aatmp,*pv,*batmp,*ba,*rtmp,*pc,*multiplier,*vtmp;
15442efa7f71SHong Zhang   PetscInt  j,nz,*pj,*bjtmp,k,ncut,*jtmp;
15452efa7f71SHong Zhang 
1546182b8fbaSHong Zhang   PetscReal dt=info->dt;          /* shift=info->shiftamount; */
15472efa7f71SHong Zhang   PetscInt  nnz_max;
1548ace3abfcSBarry Smith   PetscBool missing;
1549021822bcSHong Zhang   PetscReal *vtmp_abs;
155097ef94ebSSatish Balay   MatScalar *v_work;
155197ef94ebSSatish Balay   PetscInt  *v_pivots;
15522efa7f71SHong Zhang 
1553c8342467SHong Zhang   PetscFunctionBegin;
15542efa7f71SHong Zhang   /* ------- symbolic factorization, can be reused ---------*/
15552efa7f71SHong Zhang   ierr = MatMissingDiagonal(A,&missing,&i);CHKERRQ(ierr);
1556e32f2f54SBarry Smith   if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i);
15572efa7f71SHong Zhang   adiag=a->diag;
15582efa7f71SHong Zhang 
15592efa7f71SHong Zhang   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
15602efa7f71SHong Zhang 
15612efa7f71SHong Zhang   /* bdiag is location of diagonal in factor */
1562854ce69bSBarry Smith   ierr = PetscMalloc1(mbs+1,&bdiag);CHKERRQ(ierr);
15632efa7f71SHong Zhang 
15642efa7f71SHong Zhang   /* allocate row pointers bi */
1565854ce69bSBarry Smith   ierr = PetscMalloc1(2*mbs+2,&bi);CHKERRQ(ierr);
15662efa7f71SHong Zhang 
15672efa7f71SHong Zhang   /* allocate bj and ba; max num of nonzero entries is (ai[n]+2*n*dtcount+2) */
15682efa7f71SHong Zhang   dtcount = (PetscInt)info->dtcount;
15696bce7ff8SHong Zhang   if (dtcount > mbs-1) dtcount = mbs-1;
15706bce7ff8SHong Zhang   nnz_max = ai[mbs]+2*mbs*dtcount +2;
15716da515a1SHong Zhang   /* printf("MatILUDTFactor_SeqBAIJ, bs %d, ai[mbs] %d, nnz_max  %d, dtcount %d\n",bs,ai[mbs],nnz_max,dtcount); */
1572785e854fSJed Brown   ierr    = PetscMalloc1(nnz_max,&bj);CHKERRQ(ierr);
15736bce7ff8SHong Zhang   nnz_max = nnz_max*bs2;
1574785e854fSJed Brown   ierr    = PetscMalloc1(nnz_max,&ba);CHKERRQ(ierr);
15752efa7f71SHong Zhang 
15762efa7f71SHong Zhang   /* put together the new matrix */
15770298fd71SBarry Smith   ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(B,bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr);
15783bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)B,(PetscObject)isicol);CHKERRQ(ierr);
157926fbe8dcSKarl Rupp 
15802efa7f71SHong Zhang   b               = (Mat_SeqBAIJ*)(B)->data;
15812efa7f71SHong Zhang   b->free_a       = PETSC_TRUE;
15822efa7f71SHong Zhang   b->free_ij      = PETSC_TRUE;
15832efa7f71SHong Zhang   b->singlemalloc = PETSC_FALSE;
158426fbe8dcSKarl Rupp 
15852efa7f71SHong Zhang   b->a    = ba;
15862efa7f71SHong Zhang   b->j    = bj;
15872efa7f71SHong Zhang   b->i    = bi;
15882efa7f71SHong Zhang   b->diag = bdiag;
15892efa7f71SHong Zhang   b->ilen = 0;
15902efa7f71SHong Zhang   b->imax = 0;
15912efa7f71SHong Zhang   b->row  = isrow;
15922efa7f71SHong Zhang   b->col  = iscol;
159326fbe8dcSKarl Rupp 
15942efa7f71SHong Zhang   ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
15952efa7f71SHong Zhang   ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
159626fbe8dcSKarl Rupp 
15972efa7f71SHong Zhang   b->icol  = isicol;
1598854ce69bSBarry Smith   ierr     = PetscMalloc1(bs*(mbs+1),&b->solve_work);CHKERRQ(ierr);
15993bb1ff40SBarry Smith   ierr     = PetscLogObjectMemory((PetscObject)B,nnz_max*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
16002efa7f71SHong Zhang   b->maxnz = nnz_max/bs2;
16012efa7f71SHong Zhang 
1602d5f3da31SBarry Smith   (B)->factortype            = MAT_FACTOR_ILUDT;
16032efa7f71SHong Zhang   (B)->info.factor_mallocs   = 0;
16042efa7f71SHong Zhang   (B)->info.fill_ratio_given = ((PetscReal)nnz_max)/((PetscReal)(ai[mbs]*bs2));
16052efa7f71SHong Zhang   /* ------- end of symbolic factorization ---------*/
16062efa7f71SHong Zhang   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
16072efa7f71SHong Zhang   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
16082efa7f71SHong Zhang 
16092efa7f71SHong Zhang   /* linked list for storing column indices of the active row */
16102efa7f71SHong Zhang   nlnk = mbs + 1;
16112efa7f71SHong Zhang   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
16122efa7f71SHong Zhang 
16132efa7f71SHong Zhang   /* im: used by PetscLLAddSortedLU(); jtmp: working array for column indices of active row */
1614dcca6d9dSJed Brown   ierr = PetscMalloc2(mbs,&im,mbs,&jtmp);CHKERRQ(ierr);
16152efa7f71SHong Zhang   /* rtmp, vtmp: working arrays for sparse and contiguous row entries of active row */
1616dcca6d9dSJed Brown   ierr = PetscMalloc2(mbs*bs2,&rtmp,mbs*bs2,&vtmp);CHKERRQ(ierr);
1617854ce69bSBarry Smith   ierr = PetscMalloc1(mbs+1,&vtmp_abs);CHKERRQ(ierr);
1618dcca6d9dSJed Brown   ierr = PetscMalloc3(bs,&v_work,bs2,&multiplier,bs,&v_pivots);CHKERRQ(ierr);
16192efa7f71SHong Zhang 
1620e60cf9a0SBarry Smith   bi[0]       = 0;
16212efa7f71SHong Zhang   bdiag[0]    = (nnz_max/bs2)-1; /* location of diagonal in factor B */
16226bce7ff8SHong Zhang   bi[2*mbs+1] = bdiag[0]+1; /* endof bj and ba array */
16232efa7f71SHong Zhang   for (i=0; i<mbs; i++) {
16242efa7f71SHong Zhang     /* copy initial fill into linked list */
16252efa7f71SHong Zhang     nzi = 0; /* nonzeros for active row i */
16262efa7f71SHong Zhang     nzi = ai[r[i]+1] - ai[r[i]];
1627e32f2f54SBarry Smith     if (!nzi) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i);
16282efa7f71SHong Zhang     nzi_al = adiag[r[i]] - ai[r[i]];
16292efa7f71SHong Zhang     nzi_au = ai[r[i]+1] - adiag[r[i]] -1;
16306da515a1SHong Zhang     /* printf("row %d, nzi_al/au %d %d\n",i,nzi_al,nzi_au); */
16312efa7f71SHong Zhang 
16322efa7f71SHong Zhang     /* load in initial unfactored row */
16332efa7f71SHong Zhang     ajtmp = aj + ai[r[i]];
16342efa7f71SHong Zhang     ierr  = PetscLLAddPerm(nzi,ajtmp,ic,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
1635fca92195SBarry Smith     ierr  = PetscMemzero(rtmp,mbs*bs2*sizeof(PetscScalar));CHKERRQ(ierr);
16362efa7f71SHong Zhang     aatmp = a->a + bs2*ai[r[i]];
16372efa7f71SHong Zhang     for (j=0; j<nzi; j++) {
16382efa7f71SHong Zhang       ierr = PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],aatmp+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr);
16392efa7f71SHong Zhang     }
16402efa7f71SHong Zhang 
16412efa7f71SHong Zhang     /* add pivot rows into linked list */
16422efa7f71SHong Zhang     row = lnk[mbs];
16432efa7f71SHong Zhang     while (row < i) {
16442efa7f71SHong Zhang       nzi_bl = bi[row+1] - bi[row] + 1;
16452efa7f71SHong Zhang       bjtmp  = bj + bdiag[row+1]+1; /* points to 1st column next to the diagonal in U */
16462efa7f71SHong Zhang       ierr   = PetscLLAddSortedLU(bjtmp,row,nlnk,lnk,lnkbt,i,nzi_bl,im);CHKERRQ(ierr);
16472efa7f71SHong Zhang       nzi   += nlnk;
16482efa7f71SHong Zhang       row    = lnk[row];
16492efa7f71SHong Zhang     }
16502efa7f71SHong Zhang 
16512efa7f71SHong Zhang     /* copy data from lnk into jtmp, then initialize lnk */
16522efa7f71SHong Zhang     ierr = PetscLLClean(mbs,mbs,nzi,lnk,jtmp,lnkbt);CHKERRQ(ierr);
16532efa7f71SHong Zhang 
16542efa7f71SHong Zhang     /* numerical factorization */
16552efa7f71SHong Zhang     bjtmp = jtmp;
16562efa7f71SHong Zhang     row   = *bjtmp++; /* 1st pivot row */
16572efa7f71SHong Zhang 
16582efa7f71SHong Zhang     while  (row < i) {
16592efa7f71SHong Zhang       pc = rtmp + bs2*row;
16602efa7f71SHong Zhang       pv = ba + bs2*bdiag[row]; /* inv(diag) of the pivot row */
166196b95a6bSBarry Smith       PetscKernel_A_gets_A_times_B(bs,pc,pv,multiplier); /* pc= multiplier = pc*inv(diag[row]) */
1662a25532f0SBarry Smith       ierr = MatBlockAbs_private(1,bs2,pc,vtmp_abs);CHKERRQ(ierr);
16632efa7f71SHong Zhang       if (vtmp_abs[0] > dt) { /* apply tolerance dropping rule */
16642efa7f71SHong Zhang         pj = bj + bdiag[row+1] + 1;         /* point to 1st entry of U(row,:) */
16652efa7f71SHong Zhang         pv = ba + bs2*(bdiag[row+1] + 1);
16662efa7f71SHong Zhang         nz = bdiag[row] - bdiag[row+1] - 1;         /* num of entries in U(row,:), excluding diagonal */
16672efa7f71SHong Zhang         for (j=0; j<nz; j++) {
166896b95a6bSBarry Smith           PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j);
16692efa7f71SHong Zhang         }
16706da515a1SHong Zhang         /* ierr = PetscLogFlops(bslog*(nz+1.0)-bs);CHKERRQ(ierr); */
16712efa7f71SHong Zhang       }
16722efa7f71SHong Zhang       row = *bjtmp++;
16732efa7f71SHong Zhang     }
16742efa7f71SHong Zhang 
16752efa7f71SHong Zhang     /* copy sparse rtmp into contiguous vtmp; separate L and U part */
16762efa7f71SHong Zhang     nzi_bl = 0; j = 0;
16772efa7f71SHong Zhang     while (jtmp[j] < i) { /* L-part. Note: jtmp is sorted */
16782efa7f71SHong Zhang       ierr = PetscMemcpy(vtmp+bs2*j,rtmp+bs2*jtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
16792efa7f71SHong Zhang       nzi_bl++; j++;
16802efa7f71SHong Zhang     }
16812efa7f71SHong Zhang     nzi_bu = nzi - nzi_bl -1;
16826da515a1SHong Zhang     /* printf("nzi %d, nzi_bl %d, nzi_bu %d\n",nzi,nzi_bl,nzi_bu); */
16832efa7f71SHong Zhang 
16842efa7f71SHong Zhang     while (j < nzi) { /* U-part */
16852efa7f71SHong Zhang       ierr = PetscMemcpy(vtmp+bs2*j,rtmp+bs2*jtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
16862efa7f71SHong Zhang       /*
16872efa7f71SHong Zhang       printf(" col %d: ",jtmp[j]);
16882efa7f71SHong Zhang       for (j1=0; j1<bs2; j1++) printf(" %g",*(vtmp+bs2*j+j1));
16892efa7f71SHong Zhang       printf(" \n");
16902efa7f71SHong Zhang       */
16912efa7f71SHong Zhang       j++;
16922efa7f71SHong Zhang     }
16932efa7f71SHong Zhang 
1694a25532f0SBarry Smith     ierr = MatBlockAbs_private(nzi,bs2,vtmp,vtmp_abs);CHKERRQ(ierr);
16952efa7f71SHong Zhang     /*
16962efa7f71SHong Zhang     printf(" row %d, nzi %d, vtmp_abs\n",i,nzi);
16972efa7f71SHong Zhang     for (j1=0; j1<nzi; j1++) printf(" (%d %g),",jtmp[j1],vtmp_abs[j1]);
16982efa7f71SHong Zhang     printf(" \n");
16992efa7f71SHong Zhang     */
17002efa7f71SHong Zhang     bjtmp = bj + bi[i];
17012efa7f71SHong Zhang     batmp = ba + bs2*bi[i];
17022efa7f71SHong Zhang     /* apply level dropping rule to L part */
17032efa7f71SHong Zhang     ncut = nzi_al + dtcount;
17042efa7f71SHong Zhang     if (ncut < nzi_bl) {
1705021822bcSHong Zhang       ierr = PetscSortSplitReal(ncut,nzi_bl,vtmp_abs,jtmp);CHKERRQ(ierr);
17062efa7f71SHong Zhang       ierr = PetscSortIntWithScalarArray(ncut,jtmp,vtmp);CHKERRQ(ierr);
17072efa7f71SHong Zhang     } else {
17082efa7f71SHong Zhang       ncut = nzi_bl;
17092efa7f71SHong Zhang     }
17102efa7f71SHong Zhang     for (j=0; j<ncut; j++) {
17112efa7f71SHong Zhang       bjtmp[j] = jtmp[j];
17122efa7f71SHong Zhang       ierr     = PetscMemcpy(batmp+bs2*j,rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
17132efa7f71SHong Zhang       /*
17142efa7f71SHong Zhang       printf(" col %d: ",bjtmp[j]);
17152efa7f71SHong Zhang       for (j1=0; j1<bs2; j1++) printf(" %g,",*(batmp+bs2*j+j1));
17162efa7f71SHong Zhang       printf("\n");
17172efa7f71SHong Zhang       */
17182efa7f71SHong Zhang     }
17192efa7f71SHong Zhang     bi[i+1] = bi[i] + ncut;
17202efa7f71SHong Zhang     nzi     = ncut + 1;
17212efa7f71SHong Zhang 
17222efa7f71SHong Zhang     /* apply level dropping rule to U part */
17232efa7f71SHong Zhang     ncut = nzi_au + dtcount;
17242efa7f71SHong Zhang     if (ncut < nzi_bu) {
1725021822bcSHong Zhang       ierr = PetscSortSplitReal(ncut,nzi_bu,vtmp_abs+nzi_bl+1,jtmp+nzi_bl+1);CHKERRQ(ierr);
17262efa7f71SHong Zhang       ierr = PetscSortIntWithScalarArray(ncut,jtmp+nzi_bl+1,vtmp+nzi_bl+1);CHKERRQ(ierr);
17272efa7f71SHong Zhang     } else {
17282efa7f71SHong Zhang       ncut = nzi_bu;
17292efa7f71SHong Zhang     }
17302efa7f71SHong Zhang     nzi += ncut;
17312efa7f71SHong Zhang 
17322efa7f71SHong Zhang     /* mark bdiagonal */
17332efa7f71SHong Zhang     bdiag[i+1]    = bdiag[i] - (ncut + 1);
17346bce7ff8SHong Zhang     bi[2*mbs - i] = bi[2*mbs - i +1] - (ncut + 1);
17356bce7ff8SHong Zhang 
17362efa7f71SHong Zhang     bjtmp  = bj + bdiag[i];
17372efa7f71SHong Zhang     batmp  = ba + bs2*bdiag[i];
17382efa7f71SHong Zhang     ierr   = PetscMemcpy(batmp,rtmp+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
17392efa7f71SHong Zhang     *bjtmp = i;
17402efa7f71SHong Zhang     /*
17412efa7f71SHong Zhang     printf(" diag %d: ",*bjtmp);
17422efa7f71SHong Zhang     for (j=0; j<bs2; j++) {
17432efa7f71SHong Zhang       printf(" %g,",batmp[j]);
17442efa7f71SHong Zhang     }
17452efa7f71SHong Zhang     printf("\n");
17462efa7f71SHong Zhang     */
17472efa7f71SHong Zhang     bjtmp = bj + bdiag[i+1]+1;
17482efa7f71SHong Zhang     batmp = ba + (bdiag[i+1]+1)*bs2;
17492efa7f71SHong Zhang 
17502efa7f71SHong Zhang     for (k=0; k<ncut; k++) {
17512efa7f71SHong Zhang       bjtmp[k] = jtmp[nzi_bl+1+k];
17522efa7f71SHong Zhang       ierr     = PetscMemcpy(batmp+bs2*k,rtmp+bs2*bjtmp[k],bs2*sizeof(MatScalar));CHKERRQ(ierr);
17532efa7f71SHong Zhang       /*
17542efa7f71SHong Zhang       printf(" col %d:",bjtmp[k]);
17552efa7f71SHong Zhang       for (j1=0; j1<bs2; j1++) printf(" %g,",*(batmp+bs2*k+j1));
17562efa7f71SHong Zhang       printf("\n");
17572efa7f71SHong Zhang       */
17582efa7f71SHong Zhang     }
17592efa7f71SHong Zhang 
17602efa7f71SHong Zhang     im[i] = nzi; /* used by PetscLLAddSortedLU() */
17612efa7f71SHong Zhang 
17622efa7f71SHong Zhang     /* invert diagonal block for simplier triangular solves - add shift??? */
17632efa7f71SHong Zhang     batmp = ba + bs2*bdiag[i];
176496b95a6bSBarry Smith     ierr  = PetscKernel_A_gets_inverse_A(bs,batmp,v_pivots,v_work);CHKERRQ(ierr);
17652efa7f71SHong Zhang   } /* for (i=0; i<mbs; i++) */
1766fca92195SBarry Smith   ierr = PetscFree3(v_work,multiplier,v_pivots);CHKERRQ(ierr);
17672efa7f71SHong Zhang 
17686da515a1SHong Zhang   /* printf("end of L %d, beginning of U %d\n",bi[mbs],bdiag[mbs]); */
1769e32f2f54SBarry Smith   if (bi[mbs] >= bdiag[mbs]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"end of L array %d cannot >= the beginning of U array %d",bi[mbs],bdiag[mbs]);
17702efa7f71SHong Zhang 
17712efa7f71SHong Zhang   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
17722efa7f71SHong Zhang   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
17732efa7f71SHong Zhang 
17742efa7f71SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
17752efa7f71SHong Zhang 
1776fca92195SBarry Smith   ierr = PetscFree2(im,jtmp);CHKERRQ(ierr);
1777fca92195SBarry Smith   ierr = PetscFree2(rtmp,vtmp);CHKERRQ(ierr);
17782efa7f71SHong Zhang 
17792efa7f71SHong Zhang   ierr     = PetscLogFlops(bs2*B->cmap->n);CHKERRQ(ierr);
17802efa7f71SHong Zhang   b->maxnz = b->nz = bi[mbs] + bdiag[0] - bdiag[mbs];
17812efa7f71SHong Zhang 
17822efa7f71SHong Zhang   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
17832efa7f71SHong Zhang   ierr = ISIdentity(isicol,&icol_identity);CHKERRQ(ierr);
17842efa7f71SHong Zhang   if (row_identity && icol_identity) {
17854dd39f65SShri Abhyankar     B->ops->solve = MatSolve_SeqBAIJ_N_NaturalOrdering;
17862efa7f71SHong Zhang   } else {
17874dd39f65SShri Abhyankar     B->ops->solve = MatSolve_SeqBAIJ_N;
17882efa7f71SHong Zhang   }
17892efa7f71SHong Zhang 
17902efa7f71SHong Zhang   B->ops->solveadd          = 0;
17912efa7f71SHong Zhang   B->ops->solvetranspose    = 0;
17922efa7f71SHong Zhang   B->ops->solvetransposeadd = 0;
17932efa7f71SHong Zhang   B->ops->matsolve          = 0;
17942efa7f71SHong Zhang   B->assembled              = PETSC_TRUE;
17952efa7f71SHong Zhang   B->preallocated           = PETSC_TRUE;
1796c8342467SHong Zhang   PetscFunctionReturn(0);
1797c8342467SHong Zhang }
1798