xref: /petsc/src/mat/impls/baij/seq/baijfact2.c (revision 6c4ed00291fe44f94936dd2f04c02ab3c442e77c)
1be1d678aSKris Buschelman 
24e2b4712SSatish Balay /*
34e2b4712SSatish Balay     Factorization code for BAIJ format.
44e2b4712SSatish Balay */
54e2b4712SSatish Balay 
6c6db04a5SJed Brown #include <../src/mat/impls/baij/seq/baij.h>
7af0996ceSBarry Smith #include <petsc/private/kernels/blockinvert.h>
8c6db04a5SJed Brown #include <petscbt.h>
9c6db04a5SJed Brown #include <../src/mat/utils/freespace.h>
104e2b4712SSatish Balay 
114e2b4712SSatish Balay /* ----------------------------------------------------------------*/
1209573ac7SBarry Smith extern PetscErrorCode MatDuplicateNoCreate_SeqBAIJ(Mat,Mat,MatDuplicateOption,PetscBool);
136bce7ff8SHong Zhang 
142b0b2ea7SShri Abhyankar #undef __FUNCT__
1529a97285SShri Abhyankar #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_15_NaturalOrdering"
16766f9fbaSBarry Smith /*
17766f9fbaSBarry Smith    This is not much faster than MatLUFactorNumeric_SeqBAIJ_N() but the solve is faster at least sometimes
18766f9fbaSBarry Smith */
1929a97285SShri Abhyankar PetscErrorCode MatLUFactorNumeric_SeqBAIJ_15_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info)
202b0b2ea7SShri Abhyankar {
212b0b2ea7SShri Abhyankar   Mat             C =B;
222b0b2ea7SShri Abhyankar   Mat_SeqBAIJ     *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data;
232b0b2ea7SShri Abhyankar   PetscErrorCode  ierr;
24766f9fbaSBarry Smith   PetscInt        i,j,k,ipvt[15];
25766f9fbaSBarry Smith   const PetscInt  n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*ajtmp,*bjtmp,*bdiag=b->diag,*pj;
26766f9fbaSBarry Smith   PetscInt        nz,nzL,row;
27766f9fbaSBarry Smith   MatScalar       *rtmp,*pc,*mwork,*pv,*vv,work[225];
28766f9fbaSBarry Smith   const MatScalar *v,*aa=a->a;
292b0b2ea7SShri Abhyankar   PetscInt        bs2 = a->bs2,bs=A->rmap->bs,flg;
300fa040f9SShri Abhyankar   PetscInt        sol_ver;
312b0b2ea7SShri Abhyankar 
322b0b2ea7SShri Abhyankar   PetscFunctionBegin;
330298fd71SBarry Smith   ierr = PetscOptionsGetInt(((PetscObject)A)->prefix,"-sol_ver",&sol_ver,NULL);CHKERRQ(ierr);
340fa040f9SShri Abhyankar 
352b0b2ea7SShri Abhyankar   /* generate work space needed by the factorization */
36dcca6d9dSJed Brown   ierr = PetscMalloc2(bs2*n,&rtmp,bs2,&mwork);CHKERRQ(ierr);
372b0b2ea7SShri Abhyankar   ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr);
382b0b2ea7SShri Abhyankar 
392b0b2ea7SShri Abhyankar   for (i=0; i<n; i++) {
402b0b2ea7SShri Abhyankar     /* zero rtmp */
412b0b2ea7SShri Abhyankar     /* L part */
422b0b2ea7SShri Abhyankar     nz    = bi[i+1] - bi[i];
432b0b2ea7SShri Abhyankar     bjtmp = bj + bi[i];
442b0b2ea7SShri Abhyankar     for  (j=0; j<nz; j++) {
452b0b2ea7SShri Abhyankar       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
462b0b2ea7SShri Abhyankar     }
472b0b2ea7SShri Abhyankar 
482b0b2ea7SShri Abhyankar     /* U part */
492b0b2ea7SShri Abhyankar     nz    = bdiag[i] - bdiag[i+1];
502b0b2ea7SShri Abhyankar     bjtmp = bj + bdiag[i+1]+1;
512b0b2ea7SShri Abhyankar     for  (j=0; j<nz; j++) {
522b0b2ea7SShri Abhyankar       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
532b0b2ea7SShri Abhyankar     }
542b0b2ea7SShri Abhyankar 
552b0b2ea7SShri Abhyankar     /* load in initial (unfactored row) */
5629a97285SShri Abhyankar     nz    = ai[i+1] - ai[i];
5729a97285SShri Abhyankar     ajtmp = aj + ai[i];
5829a97285SShri Abhyankar     v     = aa + bs2*ai[i];
592b0b2ea7SShri Abhyankar     for (j=0; j<nz; j++) {
6029a97285SShri Abhyankar       ierr = PetscMemcpy(rtmp+bs2*ajtmp[j],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr);
612b0b2ea7SShri Abhyankar     }
622b0b2ea7SShri Abhyankar 
632b0b2ea7SShri Abhyankar     /* elimination */
642b0b2ea7SShri Abhyankar     bjtmp = bj + bi[i];
652b0b2ea7SShri Abhyankar     nzL   = bi[i+1] - bi[i];
662b0b2ea7SShri Abhyankar     for (k=0; k < nzL; k++) {
672b0b2ea7SShri Abhyankar       row = bjtmp[k];
682b0b2ea7SShri Abhyankar       pc  = rtmp + bs2*row;
69c35f09e5SBarry Smith       for (flg=0,j=0; j<bs2; j++) {
70c35f09e5SBarry Smith         if (pc[j]!=0.0) {
71c35f09e5SBarry Smith           flg = 1;
72c35f09e5SBarry Smith           break;
73c35f09e5SBarry Smith         }
74c35f09e5SBarry Smith       }
752b0b2ea7SShri Abhyankar       if (flg) {
762b0b2ea7SShri Abhyankar         pv = b->a + bs2*bdiag[row];
7796b95a6bSBarry Smith         PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork);
7896b95a6bSBarry Smith         /*ierr = PetscKernel_A_gets_A_times_B_15(pc,pv,mwork);CHKERRQ(ierr);*/
792b0b2ea7SShri Abhyankar         pj = b->j + bdiag[row+1]+1; /* begining of U(row,:) */
802b0b2ea7SShri Abhyankar         pv = b->a + bs2*(bdiag[row+1]+1);
812b0b2ea7SShri Abhyankar         nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries inU(row,:), excluding diag */
822b0b2ea7SShri Abhyankar         for (j=0; j<nz; j++) {
83766f9fbaSBarry Smith           vv = rtmp + bs2*pj[j];
8496b95a6bSBarry Smith           PetscKernel_A_gets_A_minus_B_times_C(bs,vv,pc,pv);
8596b95a6bSBarry Smith           /* ierr = PetscKernel_A_gets_A_minus_B_times_C_15(vv,pc,pv);CHKERRQ(ierr); */
862b0b2ea7SShri Abhyankar           pv += bs2;
872b0b2ea7SShri Abhyankar         }
88766f9fbaSBarry Smith         ierr = PetscLogFlops(2*bs2*bs*(nz+1)-bs2);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
892b0b2ea7SShri Abhyankar       }
902b0b2ea7SShri Abhyankar     }
912b0b2ea7SShri Abhyankar 
922b0b2ea7SShri Abhyankar     /* finished row so stick it into b->a */
932b0b2ea7SShri Abhyankar     /* L part */
942b0b2ea7SShri Abhyankar     pv = b->a + bs2*bi[i];
952b0b2ea7SShri Abhyankar     pj = b->j + bi[i];
962b0b2ea7SShri Abhyankar     nz = bi[i+1] - bi[i];
972b0b2ea7SShri Abhyankar     for (j=0; j<nz; j++) {
982b0b2ea7SShri Abhyankar       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
992b0b2ea7SShri Abhyankar     }
1002b0b2ea7SShri Abhyankar 
1012b0b2ea7SShri Abhyankar     /* Mark diagonal and invert diagonal for simplier triangular solves */
1022b0b2ea7SShri Abhyankar     pv   = b->a + bs2*bdiag[i];
1032b0b2ea7SShri Abhyankar     pj   = b->j + bdiag[i];
1042b0b2ea7SShri Abhyankar     ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr);
10596b95a6bSBarry Smith     /* PetscKernel_A_gets_inverse_A(bs,pv,pivots,work); */
10696b95a6bSBarry Smith     ierr = PetscKernel_A_gets_inverse_A_15(pv,ipvt,work,info->shiftamount);CHKERRQ(ierr);
1072b0b2ea7SShri Abhyankar 
1082b0b2ea7SShri Abhyankar     /* U part */
1092b0b2ea7SShri Abhyankar     pv = b->a + bs2*(bdiag[i+1]+1);
1102b0b2ea7SShri Abhyankar     pj = b->j + bdiag[i+1]+1;
1112b0b2ea7SShri Abhyankar     nz = bdiag[i] - bdiag[i+1] - 1;
1122b0b2ea7SShri Abhyankar     for (j=0; j<nz; j++) {
1132b0b2ea7SShri Abhyankar       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
1142b0b2ea7SShri Abhyankar     }
1152b0b2ea7SShri Abhyankar   }
1162b0b2ea7SShri Abhyankar 
1172b0b2ea7SShri Abhyankar   ierr = PetscFree2(rtmp,mwork);CHKERRQ(ierr);
11826fbe8dcSKarl Rupp 
119832cc040SShri Abhyankar   C->ops->solve          = MatSolve_SeqBAIJ_15_NaturalOrdering_ver1;
120766f9fbaSBarry Smith   C->ops->solvetranspose = MatSolve_SeqBAIJ_N_NaturalOrdering;
1212b0b2ea7SShri Abhyankar   C->assembled           = PETSC_TRUE;
12226fbe8dcSKarl Rupp 
123766f9fbaSBarry Smith   ierr = PetscLogFlops(1.333333333333*bs*bs2*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
1242b0b2ea7SShri Abhyankar   PetscFunctionReturn(0);
1252b0b2ea7SShri Abhyankar }
1262b0b2ea7SShri Abhyankar 
1276bce7ff8SHong Zhang #undef __FUNCT__
1284dd39f65SShri Abhyankar #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_N"
1294dd39f65SShri Abhyankar PetscErrorCode MatLUFactorNumeric_SeqBAIJ_N(Mat B,Mat A,const MatFactorInfo *info)
1306bce7ff8SHong Zhang {
1316bce7ff8SHong Zhang   Mat            C     =B;
1326bce7ff8SHong Zhang   Mat_SeqBAIJ    *a    =(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data;
1336bce7ff8SHong Zhang   IS             isrow = b->row,isicol = b->icol;
1346bce7ff8SHong Zhang   PetscErrorCode ierr;
1355a586d82SBarry Smith   const PetscInt *r,*ic;
1366bce7ff8SHong Zhang   PetscInt       i,j,k,n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
1376bce7ff8SHong Zhang   PetscInt       *ajtmp,*bjtmp,nz,nzL,row,*bdiag=b->diag,*pj;
138b588c5a2SHong Zhang   MatScalar      *rtmp,*pc,*mwork,*v,*pv,*aa=a->a;
139914a18a2SHong Zhang   PetscInt       bs=A->rmap->bs,bs2 = a->bs2,*v_pivots,flg;
140914a18a2SHong Zhang   MatScalar      *v_work;
141ace3abfcSBarry Smith   PetscBool      col_identity,row_identity,both_identity;
1426bce7ff8SHong Zhang 
1436bce7ff8SHong Zhang   PetscFunctionBegin;
1446bce7ff8SHong Zhang   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
1456bce7ff8SHong Zhang   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
146ae3d28f0SHong Zhang 
147785e854fSJed Brown   ierr = PetscMalloc1(bs2*n,&rtmp);CHKERRQ(ierr);
148fca92195SBarry Smith   ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr);
1496bce7ff8SHong Zhang 
150914a18a2SHong Zhang   /* generate work space needed by dense LU factorization */
151dcca6d9dSJed Brown   ierr = PetscMalloc3(bs,&v_work,bs2,&mwork,bs,&v_pivots);CHKERRQ(ierr);
152914a18a2SHong Zhang 
1536bce7ff8SHong Zhang   for (i=0; i<n; i++) {
1546bce7ff8SHong Zhang     /* zero rtmp */
1556bce7ff8SHong Zhang     /* L part */
1566bce7ff8SHong Zhang     nz    = bi[i+1] - bi[i];
1576bce7ff8SHong Zhang     bjtmp = bj + bi[i];
158914a18a2SHong Zhang     for  (j=0; j<nz; j++) {
159914a18a2SHong Zhang       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
160914a18a2SHong Zhang     }
1616bce7ff8SHong Zhang 
1626bce7ff8SHong Zhang     /* U part */
1631a83e813SShri Abhyankar     nz    = bdiag[i] - bdiag[i+1];
1641a83e813SShri Abhyankar     bjtmp = bj + bdiag[i+1]+1;
1651a83e813SShri Abhyankar     for  (j=0; j<nz; j++) {
1661a83e813SShri Abhyankar       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
1671a83e813SShri Abhyankar     }
1681a83e813SShri Abhyankar 
1691a83e813SShri Abhyankar     /* load in initial (unfactored row) */
1701a83e813SShri Abhyankar     nz    = ai[r[i]+1] - ai[r[i]];
1711a83e813SShri Abhyankar     ajtmp = aj + ai[r[i]];
1721a83e813SShri Abhyankar     v     = aa + bs2*ai[r[i]];
1731a83e813SShri Abhyankar     for (j=0; j<nz; j++) {
1741a83e813SShri Abhyankar       ierr = PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1751a83e813SShri Abhyankar     }
1761a83e813SShri Abhyankar 
1771a83e813SShri Abhyankar     /* elimination */
1781a83e813SShri Abhyankar     bjtmp = bj + bi[i];
1791a83e813SShri Abhyankar     nzL   = bi[i+1] - bi[i];
1801a83e813SShri Abhyankar     for (k=0; k < nzL; k++) {
1811a83e813SShri Abhyankar       row = bjtmp[k];
1821a83e813SShri Abhyankar       pc  = rtmp + bs2*row;
183c35f09e5SBarry Smith       for (flg=0,j=0; j<bs2; j++) {
184c35f09e5SBarry Smith         if (pc[j]!=0.0) {
185c35f09e5SBarry Smith           flg = 1;
186c35f09e5SBarry Smith           break;
187c35f09e5SBarry Smith         }
188c35f09e5SBarry Smith       }
1891a83e813SShri Abhyankar       if (flg) {
1901a83e813SShri Abhyankar         pv = b->a + bs2*bdiag[row];
19196b95a6bSBarry Smith         PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); /* *pc = *pc * (*pv); */
1921a83e813SShri Abhyankar         pj = b->j + bdiag[row+1]+1;         /* begining of U(row,:) */
1931a83e813SShri Abhyankar         pv = b->a + bs2*(bdiag[row+1]+1);
1941a83e813SShri Abhyankar         nz = bdiag[row] - bdiag[row+1] - 1;         /* num of entries inU(row,:), excluding diag */
1951a83e813SShri Abhyankar         for (j=0; j<nz; j++) {
19696b95a6bSBarry Smith           PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j);
1971a83e813SShri Abhyankar         }
1981a83e813SShri Abhyankar         ierr = PetscLogFlops(2*bs2*bs*(nz+1)-bs2);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
1991a83e813SShri Abhyankar       }
2001a83e813SShri Abhyankar     }
2011a83e813SShri Abhyankar 
2021a83e813SShri Abhyankar     /* finished row so stick it into b->a */
2031a83e813SShri Abhyankar     /* L part */
2041a83e813SShri Abhyankar     pv = b->a + bs2*bi[i];
2051a83e813SShri Abhyankar     pj = b->j + bi[i];
2061a83e813SShri Abhyankar     nz = bi[i+1] - bi[i];
2071a83e813SShri Abhyankar     for (j=0; j<nz; j++) {
2081a83e813SShri Abhyankar       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
2091a83e813SShri Abhyankar     }
2101a83e813SShri Abhyankar 
2111a83e813SShri Abhyankar     /* Mark diagonal and invert diagonal for simplier triangular solves */
2121a83e813SShri Abhyankar     pv = b->a + bs2*bdiag[i];
2131a83e813SShri Abhyankar     pj = b->j + bdiag[i];
214e32f2f54SBarry Smith     /* if (*pj != i)SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"row %d != *pj %d",i,*pj); */
2151a83e813SShri 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);
2171a83e813SShri Abhyankar 
2181a83e813SShri Abhyankar     /* U part */
2191a83e813SShri Abhyankar     pv = b->a + bs2*(bdiag[i+1]+1);
2201a83e813SShri Abhyankar     pj = b->j + bdiag[i+1]+1;
2211a83e813SShri Abhyankar     nz = bdiag[i] - bdiag[i+1] - 1;
2221a83e813SShri Abhyankar     for (j=0; j<nz; j++) {
2231a83e813SShri Abhyankar       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
2241a83e813SShri Abhyankar     }
2251a83e813SShri Abhyankar   }
2261a83e813SShri Abhyankar 
2271a83e813SShri Abhyankar   ierr = PetscFree(rtmp);CHKERRQ(ierr);
228fca92195SBarry Smith   ierr = PetscFree3(v_work,mwork,v_pivots);CHKERRQ(ierr);
2291a83e813SShri Abhyankar   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
2301a83e813SShri Abhyankar   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
2311a83e813SShri Abhyankar 
232ae3d28f0SHong Zhang   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
233ae3d28f0SHong Zhang   ierr = ISIdentity(isicol,&col_identity);CHKERRQ(ierr);
23426fbe8dcSKarl Rupp 
235ace3abfcSBarry Smith   both_identity = (PetscBool) (row_identity && col_identity);
236ae3d28f0SHong Zhang   if (both_identity) {
2374dd39f65SShri Abhyankar     C->ops->solve = MatSolve_SeqBAIJ_N_NaturalOrdering;
238ae3d28f0SHong Zhang   } else {
2394dd39f65SShri Abhyankar     C->ops->solve = MatSolve_SeqBAIJ_N;
240ae3d28f0SHong Zhang   }
2414dd39f65SShri Abhyankar   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_N;
242ae3d28f0SHong Zhang 
2431a83e813SShri Abhyankar   C->assembled = PETSC_TRUE;
24426fbe8dcSKarl Rupp 
245766f9fbaSBarry Smith   ierr = PetscLogFlops(1.333333333333*bs*bs2*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
2461a83e813SShri Abhyankar   PetscFunctionReturn(0);
2471a83e813SShri Abhyankar }
2481a83e813SShri Abhyankar 
2496bce7ff8SHong Zhang /*
2506bce7ff8SHong Zhang    ilu(0) with natural ordering under new data structure.
2514dd39f65SShri Abhyankar    See MatILUFactorSymbolic_SeqAIJ_ilu0() for detailed description
2524dd39f65SShri Abhyankar    because this code is almost identical to MatILUFactorSymbolic_SeqAIJ_ilu0_inplace().
2536bce7ff8SHong Zhang */
254c0c7eb62SShri Abhyankar 
2556bce7ff8SHong Zhang #undef __FUNCT__
2564dd39f65SShri Abhyankar #define __FUNCT__ "MatILUFactorSymbolic_SeqBAIJ_ilu0"
2574dd39f65SShri Abhyankar PetscErrorCode MatILUFactorSymbolic_SeqBAIJ_ilu0(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
2586bce7ff8SHong Zhang {
2596bce7ff8SHong Zhang 
2606bce7ff8SHong Zhang   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b;
2616bce7ff8SHong Zhang   PetscErrorCode ierr;
26216a2bf60SHong Zhang   PetscInt       n=a->mbs,*ai=a->i,*aj,*adiag=a->diag,bs2 = a->bs2;
26335aa4fcfSShri Abhyankar   PetscInt       i,j,nz,*bi,*bj,*bdiag,bi_temp;
26435aa4fcfSShri Abhyankar 
26535aa4fcfSShri Abhyankar   PetscFunctionBegin;
26635aa4fcfSShri Abhyankar   ierr = MatDuplicateNoCreate_SeqBAIJ(fact,A,MAT_DO_NOT_COPY_VALUES,PETSC_FALSE);CHKERRQ(ierr);
26735aa4fcfSShri Abhyankar   b    = (Mat_SeqBAIJ*)(fact)->data;
26835aa4fcfSShri Abhyankar 
26935aa4fcfSShri Abhyankar   /* allocate matrix arrays for new data structure */
270dcca6d9dSJed Brown   ierr = PetscMalloc3(bs2*ai[n]+1,&b->a,ai[n]+1,&b->j,n+1,&b->i);CHKERRQ(ierr);
2713bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)fact,ai[n]*(bs2*sizeof(PetscScalar)+sizeof(PetscInt))+(n+1)*sizeof(PetscInt));CHKERRQ(ierr);
27226fbe8dcSKarl Rupp 
27335aa4fcfSShri Abhyankar   b->singlemalloc    = PETSC_TRUE;
274379be0ddSLisandro Dalcin   b->free_a          = PETSC_TRUE;
275379be0ddSLisandro Dalcin   b->free_ij         = PETSC_TRUE;
2761e40a84eSLisandro Dalcin   fact->preallocated = PETSC_TRUE;
2771e40a84eSLisandro Dalcin   fact->assembled    = PETSC_TRUE;
27835aa4fcfSShri Abhyankar   if (!b->diag) {
279854ce69bSBarry Smith     ierr = PetscMalloc1(n+1,&b->diag);CHKERRQ(ierr);
2803bb1ff40SBarry Smith     ierr = PetscLogObjectMemory((PetscObject)fact,(n+1)*sizeof(PetscInt));CHKERRQ(ierr);
28135aa4fcfSShri Abhyankar   }
28235aa4fcfSShri Abhyankar   bdiag = b->diag;
28335aa4fcfSShri Abhyankar 
28435aa4fcfSShri Abhyankar   if (n > 0) {
28535aa4fcfSShri Abhyankar     ierr = PetscMemzero(b->a,bs2*ai[n]*sizeof(MatScalar));CHKERRQ(ierr);
28635aa4fcfSShri Abhyankar   }
28735aa4fcfSShri Abhyankar 
28835aa4fcfSShri Abhyankar   /* set bi and bj with new data structure */
28935aa4fcfSShri Abhyankar   bi = b->i;
29035aa4fcfSShri Abhyankar   bj = b->j;
29135aa4fcfSShri Abhyankar 
29235aa4fcfSShri Abhyankar   /* L part */
29335aa4fcfSShri Abhyankar   bi[0] = 0;
29435aa4fcfSShri Abhyankar   for (i=0; i<n; i++) {
29535aa4fcfSShri Abhyankar     nz      = adiag[i] - ai[i];
29635aa4fcfSShri Abhyankar     bi[i+1] = bi[i] + nz;
29735aa4fcfSShri Abhyankar     aj      = a->j + ai[i];
29835aa4fcfSShri Abhyankar     for (j=0; j<nz; j++) {
29935aa4fcfSShri Abhyankar       *bj = aj[j]; bj++;
30035aa4fcfSShri Abhyankar     }
30135aa4fcfSShri Abhyankar   }
30235aa4fcfSShri Abhyankar 
30335aa4fcfSShri Abhyankar   /* U part */
30435aa4fcfSShri Abhyankar   bi_temp  = bi[n];
30535aa4fcfSShri Abhyankar   bdiag[n] = bi[n]-1;
30635aa4fcfSShri Abhyankar   for (i=n-1; i>=0; i--) {
30735aa4fcfSShri Abhyankar     nz      = ai[i+1] - adiag[i] - 1;
30835aa4fcfSShri Abhyankar     bi_temp = bi_temp + nz + 1;
30935aa4fcfSShri Abhyankar     aj      = a->j + adiag[i] + 1;
31035aa4fcfSShri Abhyankar     for (j=0; j<nz; j++) {
31135aa4fcfSShri Abhyankar       *bj = aj[j]; bj++;
31235aa4fcfSShri Abhyankar     }
31335aa4fcfSShri Abhyankar     /* diag[i] */
31435aa4fcfSShri Abhyankar     *bj      = i; bj++;
31535aa4fcfSShri Abhyankar     bdiag[i] = bi_temp - 1;
31635aa4fcfSShri Abhyankar   }
31735aa4fcfSShri Abhyankar   PetscFunctionReturn(0);
31835aa4fcfSShri Abhyankar }
31935aa4fcfSShri Abhyankar 
32035aa4fcfSShri Abhyankar #undef __FUNCT__
3214dd39f65SShri Abhyankar #define __FUNCT__ "MatILUFactorSymbolic_SeqBAIJ"
3224dd39f65SShri Abhyankar PetscErrorCode MatILUFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
32316a2bf60SHong Zhang {
32416a2bf60SHong Zhang   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data,*b;
32516a2bf60SHong Zhang   IS                 isicol;
32616a2bf60SHong Zhang   PetscErrorCode     ierr;
32716a2bf60SHong Zhang   const PetscInt     *r,*ic;
3287fa3a6a0SHong Zhang   PetscInt           n=a->mbs,*ai=a->i,*aj=a->j,d;
32916a2bf60SHong Zhang   PetscInt           *bi,*cols,nnz,*cols_lvl;
33016a2bf60SHong Zhang   PetscInt           *bdiag,prow,fm,nzbd,reallocs=0,dcount=0;
33116a2bf60SHong Zhang   PetscInt           i,levels,diagonal_fill;
332ace3abfcSBarry Smith   PetscBool          col_identity,row_identity,both_identity;
33316a2bf60SHong Zhang   PetscReal          f;
3340298fd71SBarry Smith   PetscInt           nlnk,*lnk,*lnk_lvl=NULL;
33516a2bf60SHong Zhang   PetscBT            lnkbt;
33616a2bf60SHong Zhang   PetscInt           nzi,*bj,**bj_ptr,**bjlvl_ptr;
3370298fd71SBarry Smith   PetscFreeSpaceList free_space    =NULL,current_space=NULL;
3380298fd71SBarry Smith   PetscFreeSpaceList free_space_lvl=NULL,current_space_lvl=NULL;
339ace3abfcSBarry Smith   PetscBool          missing;
3407fa3a6a0SHong Zhang   PetscInt           bs=A->rmap->bs,bs2=a->bs2;
34116a2bf60SHong Zhang 
34216a2bf60SHong Zhang   PetscFunctionBegin;
343e32f2f54SBarry Smith   if (A->rmap->n != A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Must be square matrix, rows %D columns %D",A->rmap->n,A->cmap->n);
3446ba06ab7SHong Zhang   if (bs>1) {  /* check shifttype */
345*6c4ed002SBarry Smith     if (info->shifttype == MAT_SHIFT_NONZERO || info->shifttype == MAT_SHIFT_POSITIVE_DEFINITE) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only MAT_SHIFT_NONE and MAT_SHIFT_INBLOCKS are supported for BAIJ matrix");
3466ba06ab7SHong Zhang   }
3476ba06ab7SHong Zhang 
34816a2bf60SHong Zhang   ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr);
349e32f2f54SBarry Smith   if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
35016a2bf60SHong Zhang 
35116a2bf60SHong Zhang   f             = info->fill;
35216a2bf60SHong Zhang   levels        = (PetscInt)info->levels;
35316a2bf60SHong Zhang   diagonal_fill = (PetscInt)info->diagonal_fill;
35426fbe8dcSKarl Rupp 
35516a2bf60SHong Zhang   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
35616a2bf60SHong Zhang 
35716a2bf60SHong Zhang   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
35816a2bf60SHong Zhang   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
35926fbe8dcSKarl Rupp 
360ace3abfcSBarry Smith   both_identity = (PetscBool) (row_identity && col_identity);
36116a2bf60SHong Zhang 
3627fa3a6a0SHong Zhang   if (!levels && both_identity) {
36316a2bf60SHong Zhang     /* special case: ilu(0) with natural ordering */
3644dd39f65SShri Abhyankar     ierr = MatILUFactorSymbolic_SeqBAIJ_ilu0(fact,A,isrow,iscol,info);CHKERRQ(ierr);
3654dd39f65SShri Abhyankar     ierr = MatSeqBAIJSetNumericFactorization(fact,both_identity);CHKERRQ(ierr);
36635aa4fcfSShri Abhyankar 
367d5f3da31SBarry Smith     fact->factortype               = MAT_FACTOR_ILU;
36835aa4fcfSShri Abhyankar     (fact)->info.factor_mallocs    = 0;
36935aa4fcfSShri Abhyankar     (fact)->info.fill_ratio_given  = info->fill;
37035aa4fcfSShri Abhyankar     (fact)->info.fill_ratio_needed = 1.0;
37126fbe8dcSKarl Rupp 
37235aa4fcfSShri Abhyankar     b                = (Mat_SeqBAIJ*)(fact)->data;
37335aa4fcfSShri Abhyankar     b->row           = isrow;
37435aa4fcfSShri Abhyankar     b->col           = iscol;
37535aa4fcfSShri Abhyankar     b->icol          = isicol;
37635aa4fcfSShri Abhyankar     ierr             = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
37735aa4fcfSShri Abhyankar     ierr             = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
37835aa4fcfSShri Abhyankar     b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE;
37926fbe8dcSKarl Rupp 
380785e854fSJed Brown     ierr = PetscMalloc1((n+1)*bs,&b->solve_work);CHKERRQ(ierr);
38135aa4fcfSShri Abhyankar     PetscFunctionReturn(0);
38235aa4fcfSShri Abhyankar   }
38335aa4fcfSShri Abhyankar 
38435aa4fcfSShri Abhyankar   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
38535aa4fcfSShri Abhyankar   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
38635aa4fcfSShri Abhyankar 
38735aa4fcfSShri Abhyankar   /* get new row pointers */
388854ce69bSBarry Smith   ierr  = PetscMalloc1(n+1,&bi);CHKERRQ(ierr);
38935aa4fcfSShri Abhyankar   bi[0] = 0;
39035aa4fcfSShri Abhyankar   /* bdiag is location of diagonal in factor */
391854ce69bSBarry Smith   ierr     = PetscMalloc1(n+1,&bdiag);CHKERRQ(ierr);
39235aa4fcfSShri Abhyankar   bdiag[0] = 0;
39335aa4fcfSShri Abhyankar 
394dcca6d9dSJed Brown   ierr = PetscMalloc2(n,&bj_ptr,n,&bjlvl_ptr);CHKERRQ(ierr);
39535aa4fcfSShri Abhyankar 
39635aa4fcfSShri Abhyankar   /* create a linked list for storing column indices of the active row */
39735aa4fcfSShri Abhyankar   nlnk = n + 1;
39835aa4fcfSShri Abhyankar   ierr = PetscIncompleteLLCreate(n,n,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
39935aa4fcfSShri Abhyankar 
40035aa4fcfSShri Abhyankar   /* initial FreeSpace size is f*(ai[n]+1) */
40135aa4fcfSShri Abhyankar   ierr              = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr);
40235aa4fcfSShri Abhyankar   current_space     = free_space;
40335aa4fcfSShri Abhyankar   ierr              = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space_lvl);CHKERRQ(ierr);
40435aa4fcfSShri Abhyankar   current_space_lvl = free_space_lvl;
40535aa4fcfSShri Abhyankar 
40635aa4fcfSShri Abhyankar   for (i=0; i<n; i++) {
40735aa4fcfSShri Abhyankar     nzi = 0;
40835aa4fcfSShri Abhyankar     /* copy current row into linked list */
40935aa4fcfSShri Abhyankar     nnz = ai[r[i]+1] - ai[r[i]];
410e32f2f54SBarry Smith     if (!nnz) 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);
41135aa4fcfSShri Abhyankar     cols   = aj + ai[r[i]];
41235aa4fcfSShri Abhyankar     lnk[i] = -1; /* marker to indicate if diagonal exists */
41335aa4fcfSShri Abhyankar     ierr   = PetscIncompleteLLInit(nnz,cols,n,ic,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
41435aa4fcfSShri Abhyankar     nzi   += nlnk;
41535aa4fcfSShri Abhyankar 
41635aa4fcfSShri Abhyankar     /* make sure diagonal entry is included */
41735aa4fcfSShri Abhyankar     if (diagonal_fill && lnk[i] == -1) {
41835aa4fcfSShri Abhyankar       fm = n;
41935aa4fcfSShri Abhyankar       while (lnk[fm] < i) fm = lnk[fm];
42035aa4fcfSShri Abhyankar       lnk[i]     = lnk[fm]; /* insert diagonal into linked list */
42135aa4fcfSShri Abhyankar       lnk[fm]    = i;
42235aa4fcfSShri Abhyankar       lnk_lvl[i] = 0;
42335aa4fcfSShri Abhyankar       nzi++; dcount++;
42435aa4fcfSShri Abhyankar     }
42535aa4fcfSShri Abhyankar 
42635aa4fcfSShri Abhyankar     /* add pivot rows into the active row */
42735aa4fcfSShri Abhyankar     nzbd = 0;
42835aa4fcfSShri Abhyankar     prow = lnk[n];
42935aa4fcfSShri Abhyankar     while (prow < i) {
43035aa4fcfSShri Abhyankar       nnz      = bdiag[prow];
43135aa4fcfSShri Abhyankar       cols     = bj_ptr[prow] + nnz + 1;
43235aa4fcfSShri Abhyankar       cols_lvl = bjlvl_ptr[prow] + nnz + 1;
43335aa4fcfSShri Abhyankar       nnz      = bi[prow+1] - bi[prow] - nnz - 1;
43426fbe8dcSKarl Rupp 
43535aa4fcfSShri Abhyankar       ierr = PetscILULLAddSorted(nnz,cols,levels,cols_lvl,prow,nlnk,lnk,lnk_lvl,lnkbt,prow);CHKERRQ(ierr);
43635aa4fcfSShri Abhyankar       nzi += nlnk;
43735aa4fcfSShri Abhyankar       prow = lnk[prow];
43835aa4fcfSShri Abhyankar       nzbd++;
43935aa4fcfSShri Abhyankar     }
44035aa4fcfSShri Abhyankar     bdiag[i] = nzbd;
44135aa4fcfSShri Abhyankar     bi[i+1]  = bi[i] + nzi;
44235aa4fcfSShri Abhyankar 
44335aa4fcfSShri Abhyankar     /* if free space is not available, make more free space */
44435aa4fcfSShri Abhyankar     if (current_space->local_remaining<nzi) {
44535aa4fcfSShri Abhyankar       nnz  = 2*nzi*(n - i); /* estimated and max additional space needed */
44635aa4fcfSShri Abhyankar       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
44735aa4fcfSShri Abhyankar       ierr = PetscFreeSpaceGet(nnz,&current_space_lvl);CHKERRQ(ierr);
44835aa4fcfSShri Abhyankar       reallocs++;
44935aa4fcfSShri Abhyankar     }
45035aa4fcfSShri Abhyankar 
45135aa4fcfSShri Abhyankar     /* copy data into free_space and free_space_lvl, then initialize lnk */
45235aa4fcfSShri Abhyankar     ierr = PetscIncompleteLLClean(n,n,nzi,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
45326fbe8dcSKarl Rupp 
45435aa4fcfSShri Abhyankar     bj_ptr[i]    = current_space->array;
45535aa4fcfSShri Abhyankar     bjlvl_ptr[i] = current_space_lvl->array;
45635aa4fcfSShri Abhyankar 
45735aa4fcfSShri Abhyankar     /* make sure the active row i has diagonal entry */
45865e19b50SBarry Smith     if (*(bj_ptr[i]+bdiag[i]) != i) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Row %D has missing diagonal in factored matrix\ntry running with -pc_factor_nonzeros_along_diagonal or -pc_factor_diagonal_fill",i);
45935aa4fcfSShri Abhyankar 
46035aa4fcfSShri Abhyankar     current_space->array           += nzi;
46135aa4fcfSShri Abhyankar     current_space->local_used      += nzi;
46235aa4fcfSShri Abhyankar     current_space->local_remaining -= nzi;
46326fbe8dcSKarl Rupp 
46435aa4fcfSShri Abhyankar     current_space_lvl->array           += nzi;
46535aa4fcfSShri Abhyankar     current_space_lvl->local_used      += nzi;
46635aa4fcfSShri Abhyankar     current_space_lvl->local_remaining -= nzi;
46735aa4fcfSShri Abhyankar   }
46835aa4fcfSShri Abhyankar 
46935aa4fcfSShri Abhyankar   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
47035aa4fcfSShri Abhyankar   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
47135aa4fcfSShri Abhyankar 
47235aa4fcfSShri Abhyankar   /* copy free_space into bj and free free_space; set bi, bj, bdiag in new datastructure; */
473854ce69bSBarry Smith   ierr = PetscMalloc1(bi[n]+1,&bj);CHKERRQ(ierr);
4742ce24eb6SHong Zhang   ierr = PetscFreeSpaceContiguous_LU(&free_space,bj,n,bi,bdiag);CHKERRQ(ierr);
47535aa4fcfSShri Abhyankar 
47635aa4fcfSShri Abhyankar   ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
47735aa4fcfSShri Abhyankar   ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
478fca92195SBarry Smith   ierr = PetscFree2(bj_ptr,bjlvl_ptr);CHKERRQ(ierr);
47935aa4fcfSShri Abhyankar 
48035aa4fcfSShri Abhyankar #if defined(PETSC_USE_INFO)
48135aa4fcfSShri Abhyankar   {
482aef85c9fSShri Abhyankar     PetscReal af = ((PetscReal)(bdiag[0]+1))/((PetscReal)ai[n]);
48357622a8eSBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)f,(double)af);CHKERRQ(ierr);
48457622a8eSBarry Smith     ierr = PetscInfo1(A,"Run with -[sub_]pc_factor_fill %g or use \n",(double)af);CHKERRQ(ierr);
48557622a8eSBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill([sub]pc,%g);\n",(double)af);CHKERRQ(ierr);
48635aa4fcfSShri Abhyankar     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
48735aa4fcfSShri Abhyankar     if (diagonal_fill) {
488955c1f14SBarry Smith       ierr = PetscInfo1(A,"Detected and replaced %D missing diagonals\n",dcount);CHKERRQ(ierr);
48935aa4fcfSShri Abhyankar     }
49035aa4fcfSShri Abhyankar   }
49135aa4fcfSShri Abhyankar #endif
49235aa4fcfSShri Abhyankar 
49335aa4fcfSShri Abhyankar   /* put together the new matrix */
4940298fd71SBarry Smith   ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(fact,bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr);
4953bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)fact,(PetscObject)isicol);CHKERRQ(ierr);
49626fbe8dcSKarl Rupp 
49735aa4fcfSShri Abhyankar   b               = (Mat_SeqBAIJ*)(fact)->data;
49835aa4fcfSShri Abhyankar   b->free_a       = PETSC_TRUE;
49935aa4fcfSShri Abhyankar   b->free_ij      = PETSC_TRUE;
50035aa4fcfSShri Abhyankar   b->singlemalloc = PETSC_FALSE;
50126fbe8dcSKarl Rupp 
502854ce69bSBarry Smith   ierr = PetscMalloc1(bs2*(bdiag[0]+1),&b->a);CHKERRQ(ierr);
50326fbe8dcSKarl Rupp 
50435aa4fcfSShri Abhyankar   b->j          = bj;
50535aa4fcfSShri Abhyankar   b->i          = bi;
50635aa4fcfSShri Abhyankar   b->diag       = bdiag;
50735aa4fcfSShri Abhyankar   b->free_diag  = PETSC_TRUE;
50835aa4fcfSShri Abhyankar   b->ilen       = 0;
50935aa4fcfSShri Abhyankar   b->imax       = 0;
51035aa4fcfSShri Abhyankar   b->row        = isrow;
51135aa4fcfSShri Abhyankar   b->col        = iscol;
51235aa4fcfSShri Abhyankar   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
51335aa4fcfSShri Abhyankar   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
51435aa4fcfSShri Abhyankar   b->icol       = isicol;
51526fbe8dcSKarl Rupp 
516854ce69bSBarry Smith   ierr = PetscMalloc1(bs*n+bs,&b->solve_work);CHKERRQ(ierr);
51735aa4fcfSShri Abhyankar   /* In b structure:  Free imax, ilen, old a, old j.
51835aa4fcfSShri Abhyankar      Allocate bdiag, solve_work, new a, new j */
5193bb1ff40SBarry Smith   ierr     = PetscLogObjectMemory((PetscObject)fact,(bdiag[0]+1) * (sizeof(PetscInt)+bs2*sizeof(PetscScalar)));CHKERRQ(ierr);
52035aa4fcfSShri Abhyankar   b->maxnz = b->nz = bdiag[0]+1;
52126fbe8dcSKarl Rupp 
522ae3d28f0SHong Zhang   fact->info.factor_mallocs    = reallocs;
523ae3d28f0SHong Zhang   fact->info.fill_ratio_given  = f;
524ae3d28f0SHong Zhang   fact->info.fill_ratio_needed = ((PetscReal)(bdiag[0]+1))/((PetscReal)ai[n]);
52526fbe8dcSKarl Rupp 
5264dd39f65SShri Abhyankar   ierr = MatSeqBAIJSetNumericFactorization(fact,both_identity);CHKERRQ(ierr);
52735aa4fcfSShri Abhyankar   PetscFunctionReturn(0);
52835aa4fcfSShri Abhyankar }
52935aa4fcfSShri Abhyankar 
5304e2b4712SSatish Balay /*
5314e2b4712SSatish Balay      This code is virtually identical to MatILUFactorSymbolic_SeqAIJ
5324e2b4712SSatish Balay    except that the data structure of Mat_SeqAIJ is slightly different.
5334e2b4712SSatish Balay    Not a good example of code reuse.
5344e2b4712SSatish Balay */
5354a2ae208SSatish Balay #undef __FUNCT__
53606e38f1dSHong Zhang #define __FUNCT__ "MatILUFactorSymbolic_SeqBAIJ_inplace"
53706e38f1dSHong Zhang PetscErrorCode MatILUFactorSymbolic_SeqBAIJ_inplace(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
5384e2b4712SSatish Balay {
5394e2b4712SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b;
5404e2b4712SSatish Balay   IS             isicol;
5416849ba73SBarry Smith   PetscErrorCode ierr;
5425d0c19d7SBarry Smith   const PetscInt *r,*ic,*ai = a->i,*aj = a->j,*xi;
5435d0c19d7SBarry Smith   PetscInt       prow,n = a->mbs,*ainew,*ajnew,jmax,*fill,nz,*im,*ajfill,*flev,*xitmp;
544a96a251dSBarry Smith   PetscInt       *dloc,idx,row,m,fm,nzf,nzi,reallocate = 0,dcount = 0;
545d0f46423SBarry Smith   PetscInt       incrlev,nnz,i,bs = A->rmap->bs,bs2 = a->bs2,levels,diagonal_fill,dd;
546ace3abfcSBarry Smith   PetscBool      col_identity,row_identity,both_identity,flg;
547329f5518SBarry Smith   PetscReal      f;
5484e2b4712SSatish Balay 
5494e2b4712SSatish Balay   PetscFunctionBegin;
5506bce7ff8SHong Zhang   ierr = MatMissingDiagonal_SeqBAIJ(A,&flg,&dd);CHKERRQ(ierr);
551e32f2f54SBarry Smith   if (flg) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix A is missing diagonal entry in row %D",dd);
5526bce7ff8SHong Zhang 
553435faa5fSBarry Smith   f             = info->fill;
554690b6cddSBarry Smith   levels        = (PetscInt)info->levels;
555690b6cddSBarry Smith   diagonal_fill = (PetscInt)info->diagonal_fill;
55626fbe8dcSKarl Rupp 
5574c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
55816a2bf60SHong Zhang 
559667159a5SBarry Smith   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
560667159a5SBarry Smith   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
561ace3abfcSBarry Smith   both_identity = (PetscBool) (row_identity && col_identity);
562309c388cSBarry Smith 
56341df41f0SMatthew Knepley   if (!levels && both_identity) {  /* special case copy the nonzero structure */
56416a2bf60SHong Zhang     ierr = MatDuplicateNoCreate_SeqBAIJ(fact,A,MAT_DO_NOT_COPY_VALUES,PETSC_TRUE);CHKERRQ(ierr);
5658b1456e3SHong Zhang     ierr = MatSeqBAIJSetNumericFactorization_inplace(fact,both_identity);CHKERRQ(ierr);
5666bce7ff8SHong Zhang 
567d5f3da31SBarry Smith     fact->factortype = MAT_FACTOR_ILU;
568ae3d28f0SHong Zhang     b                = (Mat_SeqBAIJ*)fact->data;
569bb3d539aSBarry Smith     b->row           = isrow;
570bb3d539aSBarry Smith     b->col           = iscol;
571bb3d539aSBarry Smith     ierr             = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
572bb3d539aSBarry Smith     ierr             = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
573bb3d539aSBarry Smith     b->icol          = isicol;
574bcd9e38bSBarry Smith     b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE;
57526fbe8dcSKarl Rupp 
576785e854fSJed Brown     ierr         = PetscMalloc1((n+1)*bs,&b->solve_work);CHKERRQ(ierr);
5776bce7ff8SHong Zhang     PetscFunctionReturn(0);
5786bce7ff8SHong Zhang   }
5796bce7ff8SHong Zhang 
5806bce7ff8SHong Zhang   /* general case perform the symbolic factorization */
5814e2b4712SSatish Balay   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
5824e2b4712SSatish Balay   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
5834e2b4712SSatish Balay 
5844e2b4712SSatish Balay   /* get new row pointers */
585854ce69bSBarry Smith   ierr     = PetscMalloc1(n+1,&ainew);CHKERRQ(ierr);
5864e2b4712SSatish Balay   ainew[0] = 0;
5874e2b4712SSatish Balay   /* don't know how many column pointers are needed so estimate */
588690b6cddSBarry Smith   jmax = (PetscInt)(f*ai[n] + 1);
589854ce69bSBarry Smith   ierr = PetscMalloc1(jmax,&ajnew);CHKERRQ(ierr);
5904e2b4712SSatish Balay   /* ajfill is level of fill for each fill entry */
591854ce69bSBarry Smith   ierr = PetscMalloc1(jmax,&ajfill);CHKERRQ(ierr);
5924e2b4712SSatish Balay   /* fill is a linked list of nonzeros in active row */
593854ce69bSBarry Smith   ierr = PetscMalloc1(n+1,&fill);CHKERRQ(ierr);
5944e2b4712SSatish Balay   /* im is level for each filled value */
595854ce69bSBarry Smith   ierr = PetscMalloc1(n+1,&im);CHKERRQ(ierr);
5964e2b4712SSatish Balay   /* dloc is location of diagonal in factor */
597854ce69bSBarry Smith   ierr    = PetscMalloc1(n+1,&dloc);CHKERRQ(ierr);
5984e2b4712SSatish Balay   dloc[0] = 0;
5994e2b4712SSatish Balay   for (prow=0; prow<n; prow++) {
600435faa5fSBarry Smith 
601435faa5fSBarry Smith     /* copy prow into linked list */
6024e2b4712SSatish Balay     nzf = nz = ai[r[prow]+1] - ai[r[prow]];
603e32f2f54SBarry Smith     if (!nz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[prow],prow);
6044e2b4712SSatish Balay     xi         = aj + ai[r[prow]];
6054e2b4712SSatish Balay     fill[n]    = n;
606435faa5fSBarry Smith     fill[prow] = -1;   /* marker for diagonal entry */
6074e2b4712SSatish Balay     while (nz--) {
6084e2b4712SSatish Balay       fm  = n;
6094e2b4712SSatish Balay       idx = ic[*xi++];
6104e2b4712SSatish Balay       do {
6114e2b4712SSatish Balay         m  = fm;
6124e2b4712SSatish Balay         fm = fill[m];
6134e2b4712SSatish Balay       } while (fm < idx);
6144e2b4712SSatish Balay       fill[m]   = idx;
6154e2b4712SSatish Balay       fill[idx] = fm;
6164e2b4712SSatish Balay       im[idx]   = 0;
6174e2b4712SSatish Balay     }
618435faa5fSBarry Smith 
619435faa5fSBarry Smith     /* make sure diagonal entry is included */
620435faa5fSBarry Smith     if (diagonal_fill && fill[prow] == -1) {
621435faa5fSBarry Smith       fm = n;
622435faa5fSBarry Smith       while (fill[fm] < prow) fm = fill[fm];
623435faa5fSBarry Smith       fill[prow] = fill[fm];    /* insert diagonal into linked list */
624435faa5fSBarry Smith       fill[fm]   = prow;
625435faa5fSBarry Smith       im[prow]   = 0;
626435faa5fSBarry Smith       nzf++;
627335d9088SBarry Smith       dcount++;
628435faa5fSBarry Smith     }
629435faa5fSBarry Smith 
6304e2b4712SSatish Balay     nzi = 0;
6314e2b4712SSatish Balay     row = fill[n];
6324e2b4712SSatish Balay     while (row < prow) {
6334e2b4712SSatish Balay       incrlev = im[row] + 1;
6344e2b4712SSatish Balay       nz      = dloc[row];
635435faa5fSBarry Smith       xi      = ajnew  + ainew[row] + nz + 1;
6364e2b4712SSatish Balay       flev    = ajfill + ainew[row] + nz + 1;
6374e2b4712SSatish Balay       nnz     = ainew[row+1] - ainew[row] - nz - 1;
6384e2b4712SSatish Balay       fm      = row;
6394e2b4712SSatish Balay       while (nnz-- > 0) {
6404e2b4712SSatish Balay         idx = *xi++;
6414e2b4712SSatish Balay         if (*flev + incrlev > levels) {
6424e2b4712SSatish Balay           flev++;
6434e2b4712SSatish Balay           continue;
6444e2b4712SSatish Balay         }
6454e2b4712SSatish Balay         do {
6464e2b4712SSatish Balay           m  = fm;
6474e2b4712SSatish Balay           fm = fill[m];
6484e2b4712SSatish Balay         } while (fm < idx);
6494e2b4712SSatish Balay         if (fm != idx) {
6504e2b4712SSatish Balay           im[idx]   = *flev + incrlev;
6514e2b4712SSatish Balay           fill[m]   = idx;
6524e2b4712SSatish Balay           fill[idx] = fm;
6534e2b4712SSatish Balay           fm        = idx;
6544e2b4712SSatish Balay           nzf++;
65526fbe8dcSKarl Rupp         } else if (im[idx] > *flev + incrlev) im[idx] = *flev+incrlev;
6564e2b4712SSatish Balay         flev++;
6574e2b4712SSatish Balay       }
6584e2b4712SSatish Balay       row = fill[row];
6594e2b4712SSatish Balay       nzi++;
6604e2b4712SSatish Balay     }
6614e2b4712SSatish Balay     /* copy new filled row into permanent storage */
6624e2b4712SSatish Balay     ainew[prow+1] = ainew[prow] + nzf;
6634e2b4712SSatish Balay     if (ainew[prow+1] > jmax) {
664ecf371e4SBarry Smith 
665ecf371e4SBarry Smith       /* estimate how much additional space we will need */
666ecf371e4SBarry Smith       /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */
667ecf371e4SBarry Smith       /* just double the memory each time */
668690b6cddSBarry Smith       PetscInt maxadd = jmax;
669ecf371e4SBarry Smith       /* maxadd = (int)(((f*ai[n]+1)*(n-prow+5))/n); */
6704e2b4712SSatish Balay       if (maxadd < nzf) maxadd = (n-prow)*(nzf+1);
6714e2b4712SSatish Balay       jmax += maxadd;
672ecf371e4SBarry Smith 
673ecf371e4SBarry Smith       /* allocate a longer ajnew and ajfill */
674785e854fSJed Brown       ierr   = PetscMalloc1(jmax,&xitmp);CHKERRQ(ierr);
6755d0c19d7SBarry Smith       ierr   = PetscMemcpy(xitmp,ajnew,ainew[prow]*sizeof(PetscInt));CHKERRQ(ierr);
676606d414cSSatish Balay       ierr   = PetscFree(ajnew);CHKERRQ(ierr);
6775d0c19d7SBarry Smith       ajnew  = xitmp;
678785e854fSJed Brown       ierr   = PetscMalloc1(jmax,&xitmp);CHKERRQ(ierr);
6795d0c19d7SBarry Smith       ierr   = PetscMemcpy(xitmp,ajfill,ainew[prow]*sizeof(PetscInt));CHKERRQ(ierr);
680606d414cSSatish Balay       ierr   = PetscFree(ajfill);CHKERRQ(ierr);
6815d0c19d7SBarry Smith       ajfill = xitmp;
682eb150c5cSKris Buschelman       reallocate++;   /* count how many reallocations are needed */
6834e2b4712SSatish Balay     }
6845d0c19d7SBarry Smith     xitmp      = ajnew + ainew[prow];
6854e2b4712SSatish Balay     flev       = ajfill + ainew[prow];
6864e2b4712SSatish Balay     dloc[prow] = nzi;
6874e2b4712SSatish Balay     fm         = fill[n];
6884e2b4712SSatish Balay     while (nzf--) {
6895d0c19d7SBarry Smith       *xitmp++ = fm;
6904e2b4712SSatish Balay       *flev++  = im[fm];
6914e2b4712SSatish Balay       fm       = fill[fm];
6924e2b4712SSatish Balay     }
693435faa5fSBarry Smith     /* make sure row has diagonal entry */
694f23aa3ddSBarry Smith     if (ajnew[ainew[prow]+dloc[prow]] != prow) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Row %D has missing diagonal in factored matrix\n\
6952401956bSBarry Smith                                                         try running with -pc_factor_nonzeros_along_diagonal or -pc_factor_diagonal_fill",prow);
696435faa5fSBarry Smith   }
697606d414cSSatish Balay   ierr = PetscFree(ajfill);CHKERRQ(ierr);
6984e2b4712SSatish Balay   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
6994e2b4712SSatish Balay   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
700606d414cSSatish Balay   ierr = PetscFree(fill);CHKERRQ(ierr);
701606d414cSSatish Balay   ierr = PetscFree(im);CHKERRQ(ierr);
7024e2b4712SSatish Balay 
7036cf91177SBarry Smith #if defined(PETSC_USE_INFO)
7044e2b4712SSatish Balay   {
705329f5518SBarry Smith     PetscReal af = ((PetscReal)ainew[n])/((PetscReal)ai[n]);
70657622a8eSBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocate,(double)f,(double)af);CHKERRQ(ierr);
70757622a8eSBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);CHKERRQ(ierr);
70857622a8eSBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%g);\n",(double)af);CHKERRQ(ierr);
709ae15b995SBarry Smith     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
710335d9088SBarry Smith     if (diagonal_fill) {
711ae15b995SBarry Smith       ierr = PetscInfo1(A,"Detected and replaced %D missing diagonals\n",dcount);CHKERRQ(ierr);
712335d9088SBarry Smith     }
7134e2b4712SSatish Balay   }
71463ba0a88SBarry Smith #endif
7154e2b4712SSatish Balay 
7164e2b4712SSatish Balay   /* put together the new matrix */
7170298fd71SBarry Smith   ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(fact,bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr);
7183bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)fact,(PetscObject)isicol);CHKERRQ(ierr);
719ae3d28f0SHong Zhang   b    = (Mat_SeqBAIJ*)fact->data;
72026fbe8dcSKarl Rupp 
721e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
722e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
7237c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
72426fbe8dcSKarl Rupp 
725785e854fSJed Brown   ierr = PetscMalloc1(bs2*ainew[n],&b->a);CHKERRQ(ierr);
72626fbe8dcSKarl Rupp 
7274e2b4712SSatish Balay   b->j          = ajnew;
7284e2b4712SSatish Balay   b->i          = ainew;
7294e2b4712SSatish Balay   for (i=0; i<n; i++) dloc[i] += ainew[i];
7304e2b4712SSatish Balay   b->diag          = dloc;
7317f53bb6cSHong Zhang   b->free_diag     = PETSC_TRUE;
7324e2b4712SSatish Balay   b->ilen          = 0;
7334e2b4712SSatish Balay   b->imax          = 0;
7344e2b4712SSatish Balay   b->row           = isrow;
7354e2b4712SSatish Balay   b->col           = iscol;
736bcd9e38bSBarry Smith   b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE;
73726fbe8dcSKarl Rupp 
738c38d4ed2SBarry Smith   ierr    = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
739c38d4ed2SBarry Smith   ierr    = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
740e51c0b9cSSatish Balay   b->icol = isicol;
741854ce69bSBarry Smith   ierr    = PetscMalloc1(bs*n+bs,&b->solve_work);CHKERRQ(ierr);
7424e2b4712SSatish Balay   /* In b structure:  Free imax, ilen, old a, old j.
7434e2b4712SSatish Balay      Allocate dloc, solve_work, new a, new j */
7443bb1ff40SBarry Smith   ierr     = PetscLogObjectMemory((PetscObject)fact,(ainew[n]-n)*(sizeof(PetscInt))+bs2*ainew[n]*sizeof(PetscScalar));CHKERRQ(ierr);
7454e2b4712SSatish Balay   b->maxnz = b->nz = ainew[n];
7464e2b4712SSatish Balay 
747ae3d28f0SHong Zhang   fact->info.factor_mallocs    = reallocate;
748ae3d28f0SHong Zhang   fact->info.fill_ratio_given  = f;
749ae3d28f0SHong Zhang   fact->info.fill_ratio_needed = ((PetscReal)ainew[n])/((PetscReal)ai[prow]);
7506bce7ff8SHong Zhang 
7518b1456e3SHong Zhang   ierr = MatSeqBAIJSetNumericFactorization_inplace(fact,both_identity);CHKERRQ(ierr);
7528661488fSKris Buschelman   PetscFunctionReturn(0);
7538661488fSKris Buschelman }
7548661488fSKris Buschelman 
755732ee342SKris Buschelman #undef __FUNCT__
7567e7071cdSKris Buschelman #define __FUNCT__ "MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE"
757dfbe8321SBarry Smith PetscErrorCode MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE(Mat A)
7587e7071cdSKris Buschelman {
75912272027SHong Zhang   /* Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; */
76012272027SHong Zhang   /* int i,*AJ=a->j,nz=a->nz; */
7615fd66863SKarl Rupp 
7625a9542e3SKris Buschelman   PetscFunctionBegin;
7637cf1b8d3SKris Buschelman   /* Undo Column scaling */
7647cf1b8d3SKris Buschelman   /*    while (nz--) { */
7657cf1b8d3SKris Buschelman   /*      AJ[i] = AJ[i]/4; */
7667cf1b8d3SKris Buschelman   /*    } */
767c115a38dSKris Buschelman   /* This should really invoke a push/pop logic, but we don't have that yet. */
7680298fd71SBarry Smith   A->ops->setunfactored = NULL;
7697cf1b8d3SKris Buschelman   PetscFunctionReturn(0);
7707cf1b8d3SKris Buschelman }
7717cf1b8d3SKris Buschelman 
7727cf1b8d3SKris Buschelman #undef __FUNCT__
7737cf1b8d3SKris Buschelman #define __FUNCT__ "MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE_usj"
774dfbe8321SBarry Smith PetscErrorCode MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE_usj(Mat A)
7757cf1b8d3SKris Buschelman {
7767cf1b8d3SKris Buschelman   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
777b24ad042SBarry Smith   PetscInt       *AJ=a->j,nz=a->nz;
7782aa5897fSKris Buschelman   unsigned short *aj=(unsigned short*)AJ;
7795fd66863SKarl Rupp 
7805a9542e3SKris Buschelman   PetscFunctionBegin;
7810b9da03eSKris Buschelman   /* Is this really necessary? */
78220235379SKris Buschelman   while (nz--) {
7830b9da03eSKris Buschelman     AJ[nz] = (int)((unsigned int)aj[nz]); /* First extend, then convert to signed. */
7847e7071cdSKris Buschelman   }
7850298fd71SBarry Smith   A->ops->setunfactored = NULL;
7867e7071cdSKris Buschelman   PetscFunctionReturn(0);
7877e7071cdSKris Buschelman }
7887e7071cdSKris Buschelman 
789732ee342SKris Buschelman 
790