xref: /petsc/src/mat/impls/baij/seq/baijfact2.c (revision ba7f0461b1b47cff793348220e9cffc43db9f100)
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 
14766f9fbaSBarry Smith /*
15766f9fbaSBarry Smith    This is not much faster than MatLUFactorNumeric_SeqBAIJ_N() but the solve is faster at least sometimes
16766f9fbaSBarry Smith */
1729a97285SShri Abhyankar PetscErrorCode MatLUFactorNumeric_SeqBAIJ_15_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info)
182b0b2ea7SShri Abhyankar {
192b0b2ea7SShri Abhyankar   Mat             C =B;
202b0b2ea7SShri Abhyankar   Mat_SeqBAIJ     *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data;
212b0b2ea7SShri Abhyankar   PetscErrorCode  ierr;
22766f9fbaSBarry Smith   PetscInt        i,j,k,ipvt[15];
23766f9fbaSBarry Smith   const PetscInt  n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*ajtmp,*bjtmp,*bdiag=b->diag,*pj;
24766f9fbaSBarry Smith   PetscInt        nz,nzL,row;
25766f9fbaSBarry Smith   MatScalar       *rtmp,*pc,*mwork,*pv,*vv,work[225];
26766f9fbaSBarry Smith   const MatScalar *v,*aa=a->a;
272b0b2ea7SShri Abhyankar   PetscInt        bs2 = a->bs2,bs=A->rmap->bs,flg;
280fa040f9SShri Abhyankar   PetscInt        sol_ver;
29a455e926SHong Zhang   PetscBool       allowzeropivot,zeropivotdetected;
302b0b2ea7SShri Abhyankar 
312b0b2ea7SShri Abhyankar   PetscFunctionBegin;
320164db54SHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
33c5929fdfSBarry Smith   ierr = PetscOptionsGetInt(NULL,((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);
105a455e926SHong Zhang     ierr = PetscKernel_A_gets_inverse_A_15(pv,ipvt,work,info->shiftamount,allowzeropivot,&zeropivotdetected);CHKERRQ(ierr);
1067b6c816cSBarry Smith     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
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 
1274dd39f65SShri Abhyankar PetscErrorCode MatLUFactorNumeric_SeqBAIJ_N(Mat B,Mat A,const MatFactorInfo *info)
1286bce7ff8SHong Zhang {
1296bce7ff8SHong Zhang   Mat            C     =B;
1306bce7ff8SHong Zhang   Mat_SeqBAIJ    *a    =(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data;
1316bce7ff8SHong Zhang   IS             isrow = b->row,isicol = b->icol;
1326bce7ff8SHong Zhang   PetscErrorCode ierr;
1335a586d82SBarry Smith   const PetscInt *r,*ic;
1346bce7ff8SHong Zhang   PetscInt       i,j,k,n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
1356bce7ff8SHong Zhang   PetscInt       *ajtmp,*bjtmp,nz,nzL,row,*bdiag=b->diag,*pj;
136b588c5a2SHong Zhang   MatScalar      *rtmp,*pc,*mwork,*v,*pv,*aa=a->a;
137914a18a2SHong Zhang   PetscInt       bs=A->rmap->bs,bs2 = a->bs2,*v_pivots,flg;
138914a18a2SHong Zhang   MatScalar      *v_work;
139ace3abfcSBarry Smith   PetscBool      col_identity,row_identity,both_identity;
1405f8bbccaSHong Zhang   PetscBool      allowzeropivot,zeropivotdetected;
1416bce7ff8SHong Zhang 
1426bce7ff8SHong Zhang   PetscFunctionBegin;
1436bce7ff8SHong Zhang   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
1446bce7ff8SHong Zhang   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
1455f8bbccaSHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
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];
2141a83e813SShri Abhyankar     ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr);
2155f8bbccaSHong Zhang 
2165f8bbccaSHong Zhang     ierr = PetscKernel_A_gets_inverse_A(bs,pv,v_pivots,v_work,allowzeropivot,&zeropivotdetected);CHKERRQ(ierr);
2177b6c816cSBarry Smith     if (zeropivotdetected) B->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
2181a83e813SShri Abhyankar 
2191a83e813SShri Abhyankar     /* U part */
2201a83e813SShri Abhyankar     pv = b->a + bs2*(bdiag[i+1]+1);
2211a83e813SShri Abhyankar     pj = b->j + bdiag[i+1]+1;
2221a83e813SShri Abhyankar     nz = bdiag[i] - bdiag[i+1] - 1;
2231a83e813SShri Abhyankar     for (j=0; j<nz; j++) {
2241a83e813SShri Abhyankar       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
2251a83e813SShri Abhyankar     }
2261a83e813SShri Abhyankar   }
2271a83e813SShri Abhyankar 
2281a83e813SShri Abhyankar   ierr = PetscFree(rtmp);CHKERRQ(ierr);
229fca92195SBarry Smith   ierr = PetscFree3(v_work,mwork,v_pivots);CHKERRQ(ierr);
2301a83e813SShri Abhyankar   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
2311a83e813SShri Abhyankar   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
2321a83e813SShri Abhyankar 
233ae3d28f0SHong Zhang   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
234ae3d28f0SHong Zhang   ierr = ISIdentity(isicol,&col_identity);CHKERRQ(ierr);
23526fbe8dcSKarl Rupp 
236ace3abfcSBarry Smith   both_identity = (PetscBool) (row_identity && col_identity);
237ae3d28f0SHong Zhang   if (both_identity) {
238*ba7f0461SHong Zhang     switch (bs) {
239*ba7f0461SHong Zhang     case 11:
240*ba7f0461SHong Zhang       C->ops->solve = MatSolve_SeqBAIJ_11_NaturalOrdering;
241*ba7f0461SHong Zhang       break;
242*ba7f0461SHong Zhang     case 12:
243*ba7f0461SHong Zhang       C->ops->solve = MatSolve_SeqBAIJ_12_NaturalOrdering;
244*ba7f0461SHong Zhang       break;
245*ba7f0461SHong Zhang     case 13:
246*ba7f0461SHong Zhang       C->ops->solve = MatSolve_SeqBAIJ_13_NaturalOrdering;
247*ba7f0461SHong Zhang       break;
248*ba7f0461SHong Zhang     case 14:
249*ba7f0461SHong Zhang       C->ops->solve = MatSolve_SeqBAIJ_14_NaturalOrdering;
250*ba7f0461SHong Zhang       break;
251*ba7f0461SHong Zhang     default:
2524dd39f65SShri Abhyankar       C->ops->solve = MatSolve_SeqBAIJ_N_NaturalOrdering;
253*ba7f0461SHong Zhang       break;
254*ba7f0461SHong Zhang     }
255ae3d28f0SHong Zhang   } else {
2564dd39f65SShri Abhyankar     C->ops->solve = MatSolve_SeqBAIJ_N;
257ae3d28f0SHong Zhang   }
2584dd39f65SShri Abhyankar   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_N;
259ae3d28f0SHong Zhang 
2601a83e813SShri Abhyankar   C->assembled = PETSC_TRUE;
26126fbe8dcSKarl Rupp 
262766f9fbaSBarry Smith   ierr = PetscLogFlops(1.333333333333*bs*bs2*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
2631a83e813SShri Abhyankar   PetscFunctionReturn(0);
2641a83e813SShri Abhyankar }
2651a83e813SShri Abhyankar 
2666bce7ff8SHong Zhang /*
2676bce7ff8SHong Zhang    ilu(0) with natural ordering under new data structure.
2684dd39f65SShri Abhyankar    See MatILUFactorSymbolic_SeqAIJ_ilu0() for detailed description
2694dd39f65SShri Abhyankar    because this code is almost identical to MatILUFactorSymbolic_SeqAIJ_ilu0_inplace().
2706bce7ff8SHong Zhang */
271c0c7eb62SShri Abhyankar 
2724dd39f65SShri Abhyankar PetscErrorCode MatILUFactorSymbolic_SeqBAIJ_ilu0(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
2736bce7ff8SHong Zhang {
2746bce7ff8SHong Zhang 
2756bce7ff8SHong Zhang   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b;
2766bce7ff8SHong Zhang   PetscErrorCode ierr;
27716a2bf60SHong Zhang   PetscInt       n=a->mbs,*ai=a->i,*aj,*adiag=a->diag,bs2 = a->bs2;
27835aa4fcfSShri Abhyankar   PetscInt       i,j,nz,*bi,*bj,*bdiag,bi_temp;
27935aa4fcfSShri Abhyankar 
28035aa4fcfSShri Abhyankar   PetscFunctionBegin;
28135aa4fcfSShri Abhyankar   ierr = MatDuplicateNoCreate_SeqBAIJ(fact,A,MAT_DO_NOT_COPY_VALUES,PETSC_FALSE);CHKERRQ(ierr);
28235aa4fcfSShri Abhyankar   b    = (Mat_SeqBAIJ*)(fact)->data;
28335aa4fcfSShri Abhyankar 
28435aa4fcfSShri Abhyankar   /* allocate matrix arrays for new data structure */
285dcca6d9dSJed Brown   ierr = PetscMalloc3(bs2*ai[n]+1,&b->a,ai[n]+1,&b->j,n+1,&b->i);CHKERRQ(ierr);
2863bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)fact,ai[n]*(bs2*sizeof(PetscScalar)+sizeof(PetscInt))+(n+1)*sizeof(PetscInt));CHKERRQ(ierr);
28726fbe8dcSKarl Rupp 
28835aa4fcfSShri Abhyankar   b->singlemalloc    = PETSC_TRUE;
289379be0ddSLisandro Dalcin   b->free_a          = PETSC_TRUE;
290379be0ddSLisandro Dalcin   b->free_ij         = PETSC_TRUE;
2911e40a84eSLisandro Dalcin   fact->preallocated = PETSC_TRUE;
2921e40a84eSLisandro Dalcin   fact->assembled    = PETSC_TRUE;
29335aa4fcfSShri Abhyankar   if (!b->diag) {
294854ce69bSBarry Smith     ierr = PetscMalloc1(n+1,&b->diag);CHKERRQ(ierr);
2953bb1ff40SBarry Smith     ierr = PetscLogObjectMemory((PetscObject)fact,(n+1)*sizeof(PetscInt));CHKERRQ(ierr);
29635aa4fcfSShri Abhyankar   }
29735aa4fcfSShri Abhyankar   bdiag = b->diag;
29835aa4fcfSShri Abhyankar 
29935aa4fcfSShri Abhyankar   if (n > 0) {
30035aa4fcfSShri Abhyankar     ierr = PetscMemzero(b->a,bs2*ai[n]*sizeof(MatScalar));CHKERRQ(ierr);
30135aa4fcfSShri Abhyankar   }
30235aa4fcfSShri Abhyankar 
30335aa4fcfSShri Abhyankar   /* set bi and bj with new data structure */
30435aa4fcfSShri Abhyankar   bi = b->i;
30535aa4fcfSShri Abhyankar   bj = b->j;
30635aa4fcfSShri Abhyankar 
30735aa4fcfSShri Abhyankar   /* L part */
30835aa4fcfSShri Abhyankar   bi[0] = 0;
30935aa4fcfSShri Abhyankar   for (i=0; i<n; i++) {
31035aa4fcfSShri Abhyankar     nz      = adiag[i] - ai[i];
31135aa4fcfSShri Abhyankar     bi[i+1] = bi[i] + nz;
31235aa4fcfSShri Abhyankar     aj      = a->j + ai[i];
31335aa4fcfSShri Abhyankar     for (j=0; j<nz; j++) {
31435aa4fcfSShri Abhyankar       *bj = aj[j]; bj++;
31535aa4fcfSShri Abhyankar     }
31635aa4fcfSShri Abhyankar   }
31735aa4fcfSShri Abhyankar 
31835aa4fcfSShri Abhyankar   /* U part */
31935aa4fcfSShri Abhyankar   bi_temp  = bi[n];
32035aa4fcfSShri Abhyankar   bdiag[n] = bi[n]-1;
32135aa4fcfSShri Abhyankar   for (i=n-1; i>=0; i--) {
32235aa4fcfSShri Abhyankar     nz      = ai[i+1] - adiag[i] - 1;
32335aa4fcfSShri Abhyankar     bi_temp = bi_temp + nz + 1;
32435aa4fcfSShri Abhyankar     aj      = a->j + adiag[i] + 1;
32535aa4fcfSShri Abhyankar     for (j=0; j<nz; j++) {
32635aa4fcfSShri Abhyankar       *bj = aj[j]; bj++;
32735aa4fcfSShri Abhyankar     }
32835aa4fcfSShri Abhyankar     /* diag[i] */
32935aa4fcfSShri Abhyankar     *bj      = i; bj++;
33035aa4fcfSShri Abhyankar     bdiag[i] = bi_temp - 1;
33135aa4fcfSShri Abhyankar   }
33235aa4fcfSShri Abhyankar   PetscFunctionReturn(0);
33335aa4fcfSShri Abhyankar }
33435aa4fcfSShri Abhyankar 
3354dd39f65SShri Abhyankar PetscErrorCode MatILUFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
33616a2bf60SHong Zhang {
33716a2bf60SHong Zhang   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data,*b;
33816a2bf60SHong Zhang   IS                 isicol;
33916a2bf60SHong Zhang   PetscErrorCode     ierr;
34016a2bf60SHong Zhang   const PetscInt     *r,*ic;
3417fa3a6a0SHong Zhang   PetscInt           n=a->mbs,*ai=a->i,*aj=a->j,d;
34216a2bf60SHong Zhang   PetscInt           *bi,*cols,nnz,*cols_lvl;
34316a2bf60SHong Zhang   PetscInt           *bdiag,prow,fm,nzbd,reallocs=0,dcount=0;
34416a2bf60SHong Zhang   PetscInt           i,levels,diagonal_fill;
345ace3abfcSBarry Smith   PetscBool          col_identity,row_identity,both_identity;
34616a2bf60SHong Zhang   PetscReal          f;
3470298fd71SBarry Smith   PetscInt           nlnk,*lnk,*lnk_lvl=NULL;
34816a2bf60SHong Zhang   PetscBT            lnkbt;
34916a2bf60SHong Zhang   PetscInt           nzi,*bj,**bj_ptr,**bjlvl_ptr;
3500298fd71SBarry Smith   PetscFreeSpaceList free_space    =NULL,current_space=NULL;
3510298fd71SBarry Smith   PetscFreeSpaceList free_space_lvl=NULL,current_space_lvl=NULL;
352ace3abfcSBarry Smith   PetscBool          missing;
3537fa3a6a0SHong Zhang   PetscInt           bs=A->rmap->bs,bs2=a->bs2;
35416a2bf60SHong Zhang 
35516a2bf60SHong Zhang   PetscFunctionBegin;
356e32f2f54SBarry 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);
3576ba06ab7SHong Zhang   if (bs>1) {  /* check shifttype */
3586c4ed002SBarry 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");
3596ba06ab7SHong Zhang   }
3606ba06ab7SHong Zhang 
36116a2bf60SHong Zhang   ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr);
362e32f2f54SBarry Smith   if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
36316a2bf60SHong Zhang 
36416a2bf60SHong Zhang   f             = info->fill;
36516a2bf60SHong Zhang   levels        = (PetscInt)info->levels;
36616a2bf60SHong Zhang   diagonal_fill = (PetscInt)info->diagonal_fill;
36726fbe8dcSKarl Rupp 
36816a2bf60SHong Zhang   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
36916a2bf60SHong Zhang 
37016a2bf60SHong Zhang   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
37116a2bf60SHong Zhang   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
37226fbe8dcSKarl Rupp 
373ace3abfcSBarry Smith   both_identity = (PetscBool) (row_identity && col_identity);
37416a2bf60SHong Zhang 
3757fa3a6a0SHong Zhang   if (!levels && both_identity) {
37616a2bf60SHong Zhang     /* special case: ilu(0) with natural ordering */
3774dd39f65SShri Abhyankar     ierr = MatILUFactorSymbolic_SeqBAIJ_ilu0(fact,A,isrow,iscol,info);CHKERRQ(ierr);
3784dd39f65SShri Abhyankar     ierr = MatSeqBAIJSetNumericFactorization(fact,both_identity);CHKERRQ(ierr);
37935aa4fcfSShri Abhyankar 
380d5f3da31SBarry Smith     fact->factortype               = MAT_FACTOR_ILU;
38135aa4fcfSShri Abhyankar     (fact)->info.factor_mallocs    = 0;
38235aa4fcfSShri Abhyankar     (fact)->info.fill_ratio_given  = info->fill;
38335aa4fcfSShri Abhyankar     (fact)->info.fill_ratio_needed = 1.0;
38426fbe8dcSKarl Rupp 
38535aa4fcfSShri Abhyankar     b                = (Mat_SeqBAIJ*)(fact)->data;
38635aa4fcfSShri Abhyankar     b->row           = isrow;
38735aa4fcfSShri Abhyankar     b->col           = iscol;
38835aa4fcfSShri Abhyankar     b->icol          = isicol;
38935aa4fcfSShri Abhyankar     ierr             = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
39035aa4fcfSShri Abhyankar     ierr             = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
39135aa4fcfSShri Abhyankar     b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE;
39226fbe8dcSKarl Rupp 
393785e854fSJed Brown     ierr = PetscMalloc1((n+1)*bs,&b->solve_work);CHKERRQ(ierr);
39435aa4fcfSShri Abhyankar     PetscFunctionReturn(0);
39535aa4fcfSShri Abhyankar   }
39635aa4fcfSShri Abhyankar 
39735aa4fcfSShri Abhyankar   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
39835aa4fcfSShri Abhyankar   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
39935aa4fcfSShri Abhyankar 
40035aa4fcfSShri Abhyankar   /* get new row pointers */
401854ce69bSBarry Smith   ierr  = PetscMalloc1(n+1,&bi);CHKERRQ(ierr);
40235aa4fcfSShri Abhyankar   bi[0] = 0;
40335aa4fcfSShri Abhyankar   /* bdiag is location of diagonal in factor */
404854ce69bSBarry Smith   ierr     = PetscMalloc1(n+1,&bdiag);CHKERRQ(ierr);
40535aa4fcfSShri Abhyankar   bdiag[0] = 0;
40635aa4fcfSShri Abhyankar 
407dcca6d9dSJed Brown   ierr = PetscMalloc2(n,&bj_ptr,n,&bjlvl_ptr);CHKERRQ(ierr);
40835aa4fcfSShri Abhyankar 
40935aa4fcfSShri Abhyankar   /* create a linked list for storing column indices of the active row */
41035aa4fcfSShri Abhyankar   nlnk = n + 1;
41135aa4fcfSShri Abhyankar   ierr = PetscIncompleteLLCreate(n,n,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
41235aa4fcfSShri Abhyankar 
41335aa4fcfSShri Abhyankar   /* initial FreeSpace size is f*(ai[n]+1) */
414f91af8c7SBarry Smith   ierr              = PetscFreeSpaceGet(PetscRealIntMultTruncate(f,ai[n]+1),&free_space);CHKERRQ(ierr);
41535aa4fcfSShri Abhyankar   current_space     = free_space;
416f91af8c7SBarry Smith   ierr              = PetscFreeSpaceGet(PetscRealIntMultTruncate(f,ai[n]+1),&free_space_lvl);CHKERRQ(ierr);
41735aa4fcfSShri Abhyankar   current_space_lvl = free_space_lvl;
41835aa4fcfSShri Abhyankar 
41935aa4fcfSShri Abhyankar   for (i=0; i<n; i++) {
42035aa4fcfSShri Abhyankar     nzi = 0;
42135aa4fcfSShri Abhyankar     /* copy current row into linked list */
42235aa4fcfSShri Abhyankar     nnz = ai[r[i]+1] - ai[r[i]];
423e32f2f54SBarry 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);
42435aa4fcfSShri Abhyankar     cols   = aj + ai[r[i]];
42535aa4fcfSShri Abhyankar     lnk[i] = -1; /* marker to indicate if diagonal exists */
42635aa4fcfSShri Abhyankar     ierr   = PetscIncompleteLLInit(nnz,cols,n,ic,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
42735aa4fcfSShri Abhyankar     nzi   += nlnk;
42835aa4fcfSShri Abhyankar 
42935aa4fcfSShri Abhyankar     /* make sure diagonal entry is included */
43035aa4fcfSShri Abhyankar     if (diagonal_fill && lnk[i] == -1) {
43135aa4fcfSShri Abhyankar       fm = n;
43235aa4fcfSShri Abhyankar       while (lnk[fm] < i) fm = lnk[fm];
43335aa4fcfSShri Abhyankar       lnk[i]     = lnk[fm]; /* insert diagonal into linked list */
43435aa4fcfSShri Abhyankar       lnk[fm]    = i;
43535aa4fcfSShri Abhyankar       lnk_lvl[i] = 0;
43635aa4fcfSShri Abhyankar       nzi++; dcount++;
43735aa4fcfSShri Abhyankar     }
43835aa4fcfSShri Abhyankar 
43935aa4fcfSShri Abhyankar     /* add pivot rows into the active row */
44035aa4fcfSShri Abhyankar     nzbd = 0;
44135aa4fcfSShri Abhyankar     prow = lnk[n];
44235aa4fcfSShri Abhyankar     while (prow < i) {
44335aa4fcfSShri Abhyankar       nnz      = bdiag[prow];
44435aa4fcfSShri Abhyankar       cols     = bj_ptr[prow] + nnz + 1;
44535aa4fcfSShri Abhyankar       cols_lvl = bjlvl_ptr[prow] + nnz + 1;
44635aa4fcfSShri Abhyankar       nnz      = bi[prow+1] - bi[prow] - nnz - 1;
44726fbe8dcSKarl Rupp 
44835aa4fcfSShri Abhyankar       ierr = PetscILULLAddSorted(nnz,cols,levels,cols_lvl,prow,nlnk,lnk,lnk_lvl,lnkbt,prow);CHKERRQ(ierr);
44935aa4fcfSShri Abhyankar       nzi += nlnk;
45035aa4fcfSShri Abhyankar       prow = lnk[prow];
45135aa4fcfSShri Abhyankar       nzbd++;
45235aa4fcfSShri Abhyankar     }
45335aa4fcfSShri Abhyankar     bdiag[i] = nzbd;
45435aa4fcfSShri Abhyankar     bi[i+1]  = bi[i] + nzi;
45535aa4fcfSShri Abhyankar 
45635aa4fcfSShri Abhyankar     /* if free space is not available, make more free space */
45735aa4fcfSShri Abhyankar     if (current_space->local_remaining<nzi) {
458f91af8c7SBarry Smith       nnz  = PetscIntMultTruncate(2,PetscIntMultTruncate(nzi,(n - i))); /* estimated and max additional space needed */
45935aa4fcfSShri Abhyankar       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
46035aa4fcfSShri Abhyankar       ierr = PetscFreeSpaceGet(nnz,&current_space_lvl);CHKERRQ(ierr);
46135aa4fcfSShri Abhyankar       reallocs++;
46235aa4fcfSShri Abhyankar     }
46335aa4fcfSShri Abhyankar 
46435aa4fcfSShri Abhyankar     /* copy data into free_space and free_space_lvl, then initialize lnk */
46535aa4fcfSShri Abhyankar     ierr = PetscIncompleteLLClean(n,n,nzi,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
46626fbe8dcSKarl Rupp 
46735aa4fcfSShri Abhyankar     bj_ptr[i]    = current_space->array;
46835aa4fcfSShri Abhyankar     bjlvl_ptr[i] = current_space_lvl->array;
46935aa4fcfSShri Abhyankar 
47035aa4fcfSShri Abhyankar     /* make sure the active row i has diagonal entry */
47165e19b50SBarry 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);
47235aa4fcfSShri Abhyankar 
47335aa4fcfSShri Abhyankar     current_space->array           += nzi;
47435aa4fcfSShri Abhyankar     current_space->local_used      += nzi;
47535aa4fcfSShri Abhyankar     current_space->local_remaining -= nzi;
47626fbe8dcSKarl Rupp 
47735aa4fcfSShri Abhyankar     current_space_lvl->array           += nzi;
47835aa4fcfSShri Abhyankar     current_space_lvl->local_used      += nzi;
47935aa4fcfSShri Abhyankar     current_space_lvl->local_remaining -= nzi;
48035aa4fcfSShri Abhyankar   }
48135aa4fcfSShri Abhyankar 
48235aa4fcfSShri Abhyankar   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
48335aa4fcfSShri Abhyankar   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
48435aa4fcfSShri Abhyankar 
48535aa4fcfSShri Abhyankar   /* copy free_space into bj and free free_space; set bi, bj, bdiag in new datastructure; */
486854ce69bSBarry Smith   ierr = PetscMalloc1(bi[n]+1,&bj);CHKERRQ(ierr);
4872ce24eb6SHong Zhang   ierr = PetscFreeSpaceContiguous_LU(&free_space,bj,n,bi,bdiag);CHKERRQ(ierr);
48835aa4fcfSShri Abhyankar 
48935aa4fcfSShri Abhyankar   ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
49035aa4fcfSShri Abhyankar   ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
491fca92195SBarry Smith   ierr = PetscFree2(bj_ptr,bjlvl_ptr);CHKERRQ(ierr);
49235aa4fcfSShri Abhyankar 
49335aa4fcfSShri Abhyankar #if defined(PETSC_USE_INFO)
49435aa4fcfSShri Abhyankar   {
495aef85c9fSShri Abhyankar     PetscReal af = ((PetscReal)(bdiag[0]+1))/((PetscReal)ai[n]);
49657622a8eSBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)f,(double)af);CHKERRQ(ierr);
49757622a8eSBarry Smith     ierr = PetscInfo1(A,"Run with -[sub_]pc_factor_fill %g or use \n",(double)af);CHKERRQ(ierr);
49857622a8eSBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill([sub]pc,%g);\n",(double)af);CHKERRQ(ierr);
49935aa4fcfSShri Abhyankar     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
50035aa4fcfSShri Abhyankar     if (diagonal_fill) {
501955c1f14SBarry Smith       ierr = PetscInfo1(A,"Detected and replaced %D missing diagonals\n",dcount);CHKERRQ(ierr);
50235aa4fcfSShri Abhyankar     }
50335aa4fcfSShri Abhyankar   }
50435aa4fcfSShri Abhyankar #endif
50535aa4fcfSShri Abhyankar 
50635aa4fcfSShri Abhyankar   /* put together the new matrix */
507367daffbSBarry Smith   ierr = MatSeqBAIJSetPreallocation(fact,bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr);
5083bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)fact,(PetscObject)isicol);CHKERRQ(ierr);
50926fbe8dcSKarl Rupp 
51035aa4fcfSShri Abhyankar   b               = (Mat_SeqBAIJ*)(fact)->data;
51135aa4fcfSShri Abhyankar   b->free_a       = PETSC_TRUE;
51235aa4fcfSShri Abhyankar   b->free_ij      = PETSC_TRUE;
51335aa4fcfSShri Abhyankar   b->singlemalloc = PETSC_FALSE;
51426fbe8dcSKarl Rupp 
515854ce69bSBarry Smith   ierr = PetscMalloc1(bs2*(bdiag[0]+1),&b->a);CHKERRQ(ierr);
51626fbe8dcSKarl Rupp 
51735aa4fcfSShri Abhyankar   b->j          = bj;
51835aa4fcfSShri Abhyankar   b->i          = bi;
51935aa4fcfSShri Abhyankar   b->diag       = bdiag;
52035aa4fcfSShri Abhyankar   b->free_diag  = PETSC_TRUE;
52135aa4fcfSShri Abhyankar   b->ilen       = 0;
52235aa4fcfSShri Abhyankar   b->imax       = 0;
52335aa4fcfSShri Abhyankar   b->row        = isrow;
52435aa4fcfSShri Abhyankar   b->col        = iscol;
52535aa4fcfSShri Abhyankar   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
52635aa4fcfSShri Abhyankar   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
52735aa4fcfSShri Abhyankar   b->icol       = isicol;
52826fbe8dcSKarl Rupp 
529854ce69bSBarry Smith   ierr = PetscMalloc1(bs*n+bs,&b->solve_work);CHKERRQ(ierr);
53035aa4fcfSShri Abhyankar   /* In b structure:  Free imax, ilen, old a, old j.
53135aa4fcfSShri Abhyankar      Allocate bdiag, solve_work, new a, new j */
5323bb1ff40SBarry Smith   ierr     = PetscLogObjectMemory((PetscObject)fact,(bdiag[0]+1) * (sizeof(PetscInt)+bs2*sizeof(PetscScalar)));CHKERRQ(ierr);
53335aa4fcfSShri Abhyankar   b->maxnz = b->nz = bdiag[0]+1;
53426fbe8dcSKarl Rupp 
535ae3d28f0SHong Zhang   fact->info.factor_mallocs    = reallocs;
536ae3d28f0SHong Zhang   fact->info.fill_ratio_given  = f;
537ae3d28f0SHong Zhang   fact->info.fill_ratio_needed = ((PetscReal)(bdiag[0]+1))/((PetscReal)ai[n]);
53826fbe8dcSKarl Rupp 
5394dd39f65SShri Abhyankar   ierr = MatSeqBAIJSetNumericFactorization(fact,both_identity);CHKERRQ(ierr);
54035aa4fcfSShri Abhyankar   PetscFunctionReturn(0);
54135aa4fcfSShri Abhyankar }
54235aa4fcfSShri Abhyankar 
5434e2b4712SSatish Balay /*
5444e2b4712SSatish Balay      This code is virtually identical to MatILUFactorSymbolic_SeqAIJ
5454e2b4712SSatish Balay    except that the data structure of Mat_SeqAIJ is slightly different.
5464e2b4712SSatish Balay    Not a good example of code reuse.
5474e2b4712SSatish Balay */
54806e38f1dSHong Zhang PetscErrorCode MatILUFactorSymbolic_SeqBAIJ_inplace(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
5494e2b4712SSatish Balay {
5504e2b4712SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b;
5514e2b4712SSatish Balay   IS             isicol;
5526849ba73SBarry Smith   PetscErrorCode ierr;
5535d0c19d7SBarry Smith   const PetscInt *r,*ic,*ai = a->i,*aj = a->j,*xi;
5545d0c19d7SBarry Smith   PetscInt       prow,n = a->mbs,*ainew,*ajnew,jmax,*fill,nz,*im,*ajfill,*flev,*xitmp;
555a96a251dSBarry Smith   PetscInt       *dloc,idx,row,m,fm,nzf,nzi,reallocate = 0,dcount = 0;
556d0f46423SBarry Smith   PetscInt       incrlev,nnz,i,bs = A->rmap->bs,bs2 = a->bs2,levels,diagonal_fill,dd;
557ace3abfcSBarry Smith   PetscBool      col_identity,row_identity,both_identity,flg;
558329f5518SBarry Smith   PetscReal      f;
5594e2b4712SSatish Balay 
5604e2b4712SSatish Balay   PetscFunctionBegin;
5616bce7ff8SHong Zhang   ierr = MatMissingDiagonal_SeqBAIJ(A,&flg,&dd);CHKERRQ(ierr);
562e32f2f54SBarry Smith   if (flg) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix A is missing diagonal entry in row %D",dd);
5636bce7ff8SHong Zhang 
564435faa5fSBarry Smith   f             = info->fill;
565690b6cddSBarry Smith   levels        = (PetscInt)info->levels;
566690b6cddSBarry Smith   diagonal_fill = (PetscInt)info->diagonal_fill;
56726fbe8dcSKarl Rupp 
5684c49b128SBarry Smith   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
56916a2bf60SHong Zhang 
570667159a5SBarry Smith   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
571667159a5SBarry Smith   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
572ace3abfcSBarry Smith   both_identity = (PetscBool) (row_identity && col_identity);
573309c388cSBarry Smith 
57441df41f0SMatthew Knepley   if (!levels && both_identity) {  /* special case copy the nonzero structure */
57516a2bf60SHong Zhang     ierr = MatDuplicateNoCreate_SeqBAIJ(fact,A,MAT_DO_NOT_COPY_VALUES,PETSC_TRUE);CHKERRQ(ierr);
5768b1456e3SHong Zhang     ierr = MatSeqBAIJSetNumericFactorization_inplace(fact,both_identity);CHKERRQ(ierr);
5776bce7ff8SHong Zhang 
578d5f3da31SBarry Smith     fact->factortype = MAT_FACTOR_ILU;
579ae3d28f0SHong Zhang     b                = (Mat_SeqBAIJ*)fact->data;
580bb3d539aSBarry Smith     b->row           = isrow;
581bb3d539aSBarry Smith     b->col           = iscol;
582bb3d539aSBarry Smith     ierr             = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
583bb3d539aSBarry Smith     ierr             = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
584bb3d539aSBarry Smith     b->icol          = isicol;
585bcd9e38bSBarry Smith     b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE;
58626fbe8dcSKarl Rupp 
587785e854fSJed Brown     ierr         = PetscMalloc1((n+1)*bs,&b->solve_work);CHKERRQ(ierr);
5886bce7ff8SHong Zhang     PetscFunctionReturn(0);
5896bce7ff8SHong Zhang   }
5906bce7ff8SHong Zhang 
5916bce7ff8SHong Zhang   /* general case perform the symbolic factorization */
5924e2b4712SSatish Balay   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
5934e2b4712SSatish Balay   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
5944e2b4712SSatish Balay 
5954e2b4712SSatish Balay   /* get new row pointers */
596854ce69bSBarry Smith   ierr     = PetscMalloc1(n+1,&ainew);CHKERRQ(ierr);
5974e2b4712SSatish Balay   ainew[0] = 0;
5984e2b4712SSatish Balay   /* don't know how many column pointers are needed so estimate */
599690b6cddSBarry Smith   jmax = (PetscInt)(f*ai[n] + 1);
600854ce69bSBarry Smith   ierr = PetscMalloc1(jmax,&ajnew);CHKERRQ(ierr);
6014e2b4712SSatish Balay   /* ajfill is level of fill for each fill entry */
602854ce69bSBarry Smith   ierr = PetscMalloc1(jmax,&ajfill);CHKERRQ(ierr);
6034e2b4712SSatish Balay   /* fill is a linked list of nonzeros in active row */
604854ce69bSBarry Smith   ierr = PetscMalloc1(n+1,&fill);CHKERRQ(ierr);
6054e2b4712SSatish Balay   /* im is level for each filled value */
606854ce69bSBarry Smith   ierr = PetscMalloc1(n+1,&im);CHKERRQ(ierr);
6074e2b4712SSatish Balay   /* dloc is location of diagonal in factor */
608854ce69bSBarry Smith   ierr    = PetscMalloc1(n+1,&dloc);CHKERRQ(ierr);
6094e2b4712SSatish Balay   dloc[0] = 0;
6104e2b4712SSatish Balay   for (prow=0; prow<n; prow++) {
611435faa5fSBarry Smith 
612435faa5fSBarry Smith     /* copy prow into linked list */
6134e2b4712SSatish Balay     nzf = nz = ai[r[prow]+1] - ai[r[prow]];
614e32f2f54SBarry 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);
6154e2b4712SSatish Balay     xi         = aj + ai[r[prow]];
6164e2b4712SSatish Balay     fill[n]    = n;
617435faa5fSBarry Smith     fill[prow] = -1;   /* marker for diagonal entry */
6184e2b4712SSatish Balay     while (nz--) {
6194e2b4712SSatish Balay       fm  = n;
6204e2b4712SSatish Balay       idx = ic[*xi++];
6214e2b4712SSatish Balay       do {
6224e2b4712SSatish Balay         m  = fm;
6234e2b4712SSatish Balay         fm = fill[m];
6244e2b4712SSatish Balay       } while (fm < idx);
6254e2b4712SSatish Balay       fill[m]   = idx;
6264e2b4712SSatish Balay       fill[idx] = fm;
6274e2b4712SSatish Balay       im[idx]   = 0;
6284e2b4712SSatish Balay     }
629435faa5fSBarry Smith 
630435faa5fSBarry Smith     /* make sure diagonal entry is included */
631435faa5fSBarry Smith     if (diagonal_fill && fill[prow] == -1) {
632435faa5fSBarry Smith       fm = n;
633435faa5fSBarry Smith       while (fill[fm] < prow) fm = fill[fm];
634435faa5fSBarry Smith       fill[prow] = fill[fm];    /* insert diagonal into linked list */
635435faa5fSBarry Smith       fill[fm]   = prow;
636435faa5fSBarry Smith       im[prow]   = 0;
637435faa5fSBarry Smith       nzf++;
638335d9088SBarry Smith       dcount++;
639435faa5fSBarry Smith     }
640435faa5fSBarry Smith 
6414e2b4712SSatish Balay     nzi = 0;
6424e2b4712SSatish Balay     row = fill[n];
6434e2b4712SSatish Balay     while (row < prow) {
6444e2b4712SSatish Balay       incrlev = im[row] + 1;
6454e2b4712SSatish Balay       nz      = dloc[row];
646435faa5fSBarry Smith       xi      = ajnew  + ainew[row] + nz + 1;
6474e2b4712SSatish Balay       flev    = ajfill + ainew[row] + nz + 1;
6484e2b4712SSatish Balay       nnz     = ainew[row+1] - ainew[row] - nz - 1;
6494e2b4712SSatish Balay       fm      = row;
6504e2b4712SSatish Balay       while (nnz-- > 0) {
6514e2b4712SSatish Balay         idx = *xi++;
6524e2b4712SSatish Balay         if (*flev + incrlev > levels) {
6534e2b4712SSatish Balay           flev++;
6544e2b4712SSatish Balay           continue;
6554e2b4712SSatish Balay         }
6564e2b4712SSatish Balay         do {
6574e2b4712SSatish Balay           m  = fm;
6584e2b4712SSatish Balay           fm = fill[m];
6594e2b4712SSatish Balay         } while (fm < idx);
6604e2b4712SSatish Balay         if (fm != idx) {
6614e2b4712SSatish Balay           im[idx]   = *flev + incrlev;
6624e2b4712SSatish Balay           fill[m]   = idx;
6634e2b4712SSatish Balay           fill[idx] = fm;
6644e2b4712SSatish Balay           fm        = idx;
6654e2b4712SSatish Balay           nzf++;
66626fbe8dcSKarl Rupp         } else if (im[idx] > *flev + incrlev) im[idx] = *flev+incrlev;
6674e2b4712SSatish Balay         flev++;
6684e2b4712SSatish Balay       }
6694e2b4712SSatish Balay       row = fill[row];
6704e2b4712SSatish Balay       nzi++;
6714e2b4712SSatish Balay     }
6724e2b4712SSatish Balay     /* copy new filled row into permanent storage */
6734e2b4712SSatish Balay     ainew[prow+1] = ainew[prow] + nzf;
6744e2b4712SSatish Balay     if (ainew[prow+1] > jmax) {
675ecf371e4SBarry Smith 
676ecf371e4SBarry Smith       /* estimate how much additional space we will need */
677ecf371e4SBarry Smith       /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */
678ecf371e4SBarry Smith       /* just double the memory each time */
679690b6cddSBarry Smith       PetscInt maxadd = jmax;
680ecf371e4SBarry Smith       /* maxadd = (int)(((f*ai[n]+1)*(n-prow+5))/n); */
6814e2b4712SSatish Balay       if (maxadd < nzf) maxadd = (n-prow)*(nzf+1);
6824e2b4712SSatish Balay       jmax += maxadd;
683ecf371e4SBarry Smith 
684ecf371e4SBarry Smith       /* allocate a longer ajnew and ajfill */
685785e854fSJed Brown       ierr   = PetscMalloc1(jmax,&xitmp);CHKERRQ(ierr);
6865d0c19d7SBarry Smith       ierr   = PetscMemcpy(xitmp,ajnew,ainew[prow]*sizeof(PetscInt));CHKERRQ(ierr);
687606d414cSSatish Balay       ierr   = PetscFree(ajnew);CHKERRQ(ierr);
6885d0c19d7SBarry Smith       ajnew  = xitmp;
689785e854fSJed Brown       ierr   = PetscMalloc1(jmax,&xitmp);CHKERRQ(ierr);
6905d0c19d7SBarry Smith       ierr   = PetscMemcpy(xitmp,ajfill,ainew[prow]*sizeof(PetscInt));CHKERRQ(ierr);
691606d414cSSatish Balay       ierr   = PetscFree(ajfill);CHKERRQ(ierr);
6925d0c19d7SBarry Smith       ajfill = xitmp;
693eb150c5cSKris Buschelman       reallocate++;   /* count how many reallocations are needed */
6944e2b4712SSatish Balay     }
6955d0c19d7SBarry Smith     xitmp      = ajnew + ainew[prow];
6964e2b4712SSatish Balay     flev       = ajfill + ainew[prow];
6974e2b4712SSatish Balay     dloc[prow] = nzi;
6984e2b4712SSatish Balay     fm         = fill[n];
6994e2b4712SSatish Balay     while (nzf--) {
7005d0c19d7SBarry Smith       *xitmp++ = fm;
7014e2b4712SSatish Balay       *flev++  = im[fm];
7024e2b4712SSatish Balay       fm       = fill[fm];
7034e2b4712SSatish Balay     }
704435faa5fSBarry Smith     /* make sure row has diagonal entry */
705f23aa3ddSBarry 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\
7062401956bSBarry Smith                                                         try running with -pc_factor_nonzeros_along_diagonal or -pc_factor_diagonal_fill",prow);
707435faa5fSBarry Smith   }
708606d414cSSatish Balay   ierr = PetscFree(ajfill);CHKERRQ(ierr);
7094e2b4712SSatish Balay   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
7104e2b4712SSatish Balay   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
711606d414cSSatish Balay   ierr = PetscFree(fill);CHKERRQ(ierr);
712606d414cSSatish Balay   ierr = PetscFree(im);CHKERRQ(ierr);
7134e2b4712SSatish Balay 
7146cf91177SBarry Smith #if defined(PETSC_USE_INFO)
7154e2b4712SSatish Balay   {
716329f5518SBarry Smith     PetscReal af = ((PetscReal)ainew[n])/((PetscReal)ai[n]);
71757622a8eSBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocate,(double)f,(double)af);CHKERRQ(ierr);
71857622a8eSBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);CHKERRQ(ierr);
71957622a8eSBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%g);\n",(double)af);CHKERRQ(ierr);
720ae15b995SBarry Smith     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
721335d9088SBarry Smith     if (diagonal_fill) {
722ae15b995SBarry Smith       ierr = PetscInfo1(A,"Detected and replaced %D missing diagonals\n",dcount);CHKERRQ(ierr);
723335d9088SBarry Smith     }
7244e2b4712SSatish Balay   }
72563ba0a88SBarry Smith #endif
7264e2b4712SSatish Balay 
7274e2b4712SSatish Balay   /* put together the new matrix */
728367daffbSBarry Smith   ierr = MatSeqBAIJSetPreallocation(fact,bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr);
7293bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)fact,(PetscObject)isicol);CHKERRQ(ierr);
730ae3d28f0SHong Zhang   b    = (Mat_SeqBAIJ*)fact->data;
73126fbe8dcSKarl Rupp 
732e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
733e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
7347c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
73526fbe8dcSKarl Rupp 
736785e854fSJed Brown   ierr = PetscMalloc1(bs2*ainew[n],&b->a);CHKERRQ(ierr);
73726fbe8dcSKarl Rupp 
7384e2b4712SSatish Balay   b->j          = ajnew;
7394e2b4712SSatish Balay   b->i          = ainew;
7404e2b4712SSatish Balay   for (i=0; i<n; i++) dloc[i] += ainew[i];
7414e2b4712SSatish Balay   b->diag          = dloc;
7427f53bb6cSHong Zhang   b->free_diag     = PETSC_TRUE;
7434e2b4712SSatish Balay   b->ilen          = 0;
7444e2b4712SSatish Balay   b->imax          = 0;
7454e2b4712SSatish Balay   b->row           = isrow;
7464e2b4712SSatish Balay   b->col           = iscol;
747bcd9e38bSBarry Smith   b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE;
74826fbe8dcSKarl Rupp 
749c38d4ed2SBarry Smith   ierr    = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
750c38d4ed2SBarry Smith   ierr    = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
751e51c0b9cSSatish Balay   b->icol = isicol;
752854ce69bSBarry Smith   ierr    = PetscMalloc1(bs*n+bs,&b->solve_work);CHKERRQ(ierr);
7534e2b4712SSatish Balay   /* In b structure:  Free imax, ilen, old a, old j.
7544e2b4712SSatish Balay      Allocate dloc, solve_work, new a, new j */
7553bb1ff40SBarry Smith   ierr     = PetscLogObjectMemory((PetscObject)fact,(ainew[n]-n)*(sizeof(PetscInt))+bs2*ainew[n]*sizeof(PetscScalar));CHKERRQ(ierr);
7564e2b4712SSatish Balay   b->maxnz = b->nz = ainew[n];
7574e2b4712SSatish Balay 
758ae3d28f0SHong Zhang   fact->info.factor_mallocs    = reallocate;
759ae3d28f0SHong Zhang   fact->info.fill_ratio_given  = f;
760ae3d28f0SHong Zhang   fact->info.fill_ratio_needed = ((PetscReal)ainew[n])/((PetscReal)ai[prow]);
7616bce7ff8SHong Zhang 
7628b1456e3SHong Zhang   ierr = MatSeqBAIJSetNumericFactorization_inplace(fact,both_identity);CHKERRQ(ierr);
7638661488fSKris Buschelman   PetscFunctionReturn(0);
7648661488fSKris Buschelman }
7658661488fSKris Buschelman 
766dfbe8321SBarry Smith PetscErrorCode MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE(Mat A)
7677e7071cdSKris Buschelman {
76812272027SHong Zhang   /* Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; */
76912272027SHong Zhang   /* int i,*AJ=a->j,nz=a->nz; */
7705fd66863SKarl Rupp 
7715a9542e3SKris Buschelman   PetscFunctionBegin;
7727cf1b8d3SKris Buschelman   /* Undo Column scaling */
7737cf1b8d3SKris Buschelman   /*    while (nz--) { */
7747cf1b8d3SKris Buschelman   /*      AJ[i] = AJ[i]/4; */
7757cf1b8d3SKris Buschelman   /*    } */
776c115a38dSKris Buschelman   /* This should really invoke a push/pop logic, but we don't have that yet. */
7770298fd71SBarry Smith   A->ops->setunfactored = NULL;
7787cf1b8d3SKris Buschelman   PetscFunctionReturn(0);
7797cf1b8d3SKris Buschelman }
7807cf1b8d3SKris Buschelman 
781dfbe8321SBarry Smith PetscErrorCode MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE_usj(Mat A)
7827cf1b8d3SKris Buschelman {
7837cf1b8d3SKris Buschelman   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
784b24ad042SBarry Smith   PetscInt       *AJ=a->j,nz=a->nz;
7852aa5897fSKris Buschelman   unsigned short *aj=(unsigned short*)AJ;
7865fd66863SKarl Rupp 
7875a9542e3SKris Buschelman   PetscFunctionBegin;
7880b9da03eSKris Buschelman   /* Is this really necessary? */
78920235379SKris Buschelman   while (nz--) {
7900b9da03eSKris Buschelman     AJ[nz] = (int)((unsigned int)aj[nz]); /* First extend, then convert to signed. */
7917e7071cdSKris Buschelman   }
7920298fd71SBarry Smith   A->ops->setunfactored = NULL;
7937e7071cdSKris Buschelman   PetscFunctionReturn(0);
7947e7071cdSKris Buschelman }
7957e7071cdSKris Buschelman 
796732ee342SKris Buschelman 
797