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) { 238ba7f0461SHong Zhang switch (bs) { 239*96e086a2SDaniel Kokron case 9: 240*96e086a2SDaniel Kokron C->ops->solve = MatSolve_SeqBAIJ_9_NaturalOrdering; 241*96e086a2SDaniel Kokron break; 242ba7f0461SHong Zhang case 11: 243ba7f0461SHong Zhang C->ops->solve = MatSolve_SeqBAIJ_11_NaturalOrdering; 244ba7f0461SHong Zhang break; 245ba7f0461SHong Zhang case 12: 246ba7f0461SHong Zhang C->ops->solve = MatSolve_SeqBAIJ_12_NaturalOrdering; 247ba7f0461SHong Zhang break; 248ba7f0461SHong Zhang case 13: 249ba7f0461SHong Zhang C->ops->solve = MatSolve_SeqBAIJ_13_NaturalOrdering; 250ba7f0461SHong Zhang break; 251ba7f0461SHong Zhang case 14: 252ba7f0461SHong Zhang C->ops->solve = MatSolve_SeqBAIJ_14_NaturalOrdering; 253ba7f0461SHong Zhang break; 254ba7f0461SHong Zhang default: 2554dd39f65SShri Abhyankar C->ops->solve = MatSolve_SeqBAIJ_N_NaturalOrdering; 256ba7f0461SHong Zhang break; 257ba7f0461SHong Zhang } 258ae3d28f0SHong Zhang } else { 2594dd39f65SShri Abhyankar C->ops->solve = MatSolve_SeqBAIJ_N; 260ae3d28f0SHong Zhang } 2614dd39f65SShri Abhyankar C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_N; 262ae3d28f0SHong Zhang 2631a83e813SShri Abhyankar C->assembled = PETSC_TRUE; 26426fbe8dcSKarl Rupp 265766f9fbaSBarry Smith ierr = PetscLogFlops(1.333333333333*bs*bs2*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 2661a83e813SShri Abhyankar PetscFunctionReturn(0); 2671a83e813SShri Abhyankar } 2681a83e813SShri Abhyankar 2696bce7ff8SHong Zhang /* 2706bce7ff8SHong Zhang ilu(0) with natural ordering under new data structure. 2714dd39f65SShri Abhyankar See MatILUFactorSymbolic_SeqAIJ_ilu0() for detailed description 2724dd39f65SShri Abhyankar because this code is almost identical to MatILUFactorSymbolic_SeqAIJ_ilu0_inplace(). 2736bce7ff8SHong Zhang */ 274c0c7eb62SShri Abhyankar 2754dd39f65SShri Abhyankar PetscErrorCode MatILUFactorSymbolic_SeqBAIJ_ilu0(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info) 2766bce7ff8SHong Zhang { 2776bce7ff8SHong Zhang 2786bce7ff8SHong Zhang Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b; 2796bce7ff8SHong Zhang PetscErrorCode ierr; 28016a2bf60SHong Zhang PetscInt n=a->mbs,*ai=a->i,*aj,*adiag=a->diag,bs2 = a->bs2; 28135aa4fcfSShri Abhyankar PetscInt i,j,nz,*bi,*bj,*bdiag,bi_temp; 28235aa4fcfSShri Abhyankar 28335aa4fcfSShri Abhyankar PetscFunctionBegin; 28435aa4fcfSShri Abhyankar ierr = MatDuplicateNoCreate_SeqBAIJ(fact,A,MAT_DO_NOT_COPY_VALUES,PETSC_FALSE);CHKERRQ(ierr); 28535aa4fcfSShri Abhyankar b = (Mat_SeqBAIJ*)(fact)->data; 28635aa4fcfSShri Abhyankar 28735aa4fcfSShri Abhyankar /* allocate matrix arrays for new data structure */ 288dcca6d9dSJed Brown ierr = PetscMalloc3(bs2*ai[n]+1,&b->a,ai[n]+1,&b->j,n+1,&b->i);CHKERRQ(ierr); 2893bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)fact,ai[n]*(bs2*sizeof(PetscScalar)+sizeof(PetscInt))+(n+1)*sizeof(PetscInt));CHKERRQ(ierr); 29026fbe8dcSKarl Rupp 29135aa4fcfSShri Abhyankar b->singlemalloc = PETSC_TRUE; 292379be0ddSLisandro Dalcin b->free_a = PETSC_TRUE; 293379be0ddSLisandro Dalcin b->free_ij = PETSC_TRUE; 2941e40a84eSLisandro Dalcin fact->preallocated = PETSC_TRUE; 2951e40a84eSLisandro Dalcin fact->assembled = PETSC_TRUE; 29635aa4fcfSShri Abhyankar if (!b->diag) { 297854ce69bSBarry Smith ierr = PetscMalloc1(n+1,&b->diag);CHKERRQ(ierr); 2983bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)fact,(n+1)*sizeof(PetscInt));CHKERRQ(ierr); 29935aa4fcfSShri Abhyankar } 30035aa4fcfSShri Abhyankar bdiag = b->diag; 30135aa4fcfSShri Abhyankar 30235aa4fcfSShri Abhyankar if (n > 0) { 30335aa4fcfSShri Abhyankar ierr = PetscMemzero(b->a,bs2*ai[n]*sizeof(MatScalar));CHKERRQ(ierr); 30435aa4fcfSShri Abhyankar } 30535aa4fcfSShri Abhyankar 30635aa4fcfSShri Abhyankar /* set bi and bj with new data structure */ 30735aa4fcfSShri Abhyankar bi = b->i; 30835aa4fcfSShri Abhyankar bj = b->j; 30935aa4fcfSShri Abhyankar 31035aa4fcfSShri Abhyankar /* L part */ 31135aa4fcfSShri Abhyankar bi[0] = 0; 31235aa4fcfSShri Abhyankar for (i=0; i<n; i++) { 31335aa4fcfSShri Abhyankar nz = adiag[i] - ai[i]; 31435aa4fcfSShri Abhyankar bi[i+1] = bi[i] + nz; 31535aa4fcfSShri Abhyankar aj = a->j + ai[i]; 31635aa4fcfSShri Abhyankar for (j=0; j<nz; j++) { 31735aa4fcfSShri Abhyankar *bj = aj[j]; bj++; 31835aa4fcfSShri Abhyankar } 31935aa4fcfSShri Abhyankar } 32035aa4fcfSShri Abhyankar 32135aa4fcfSShri Abhyankar /* U part */ 32235aa4fcfSShri Abhyankar bi_temp = bi[n]; 32335aa4fcfSShri Abhyankar bdiag[n] = bi[n]-1; 32435aa4fcfSShri Abhyankar for (i=n-1; i>=0; i--) { 32535aa4fcfSShri Abhyankar nz = ai[i+1] - adiag[i] - 1; 32635aa4fcfSShri Abhyankar bi_temp = bi_temp + nz + 1; 32735aa4fcfSShri Abhyankar aj = a->j + adiag[i] + 1; 32835aa4fcfSShri Abhyankar for (j=0; j<nz; j++) { 32935aa4fcfSShri Abhyankar *bj = aj[j]; bj++; 33035aa4fcfSShri Abhyankar } 33135aa4fcfSShri Abhyankar /* diag[i] */ 33235aa4fcfSShri Abhyankar *bj = i; bj++; 33335aa4fcfSShri Abhyankar bdiag[i] = bi_temp - 1; 33435aa4fcfSShri Abhyankar } 33535aa4fcfSShri Abhyankar PetscFunctionReturn(0); 33635aa4fcfSShri Abhyankar } 33735aa4fcfSShri Abhyankar 3384dd39f65SShri Abhyankar PetscErrorCode MatILUFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info) 33916a2bf60SHong Zhang { 34016a2bf60SHong Zhang Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b; 34116a2bf60SHong Zhang IS isicol; 34216a2bf60SHong Zhang PetscErrorCode ierr; 34316a2bf60SHong Zhang const PetscInt *r,*ic; 3447fa3a6a0SHong Zhang PetscInt n=a->mbs,*ai=a->i,*aj=a->j,d; 34516a2bf60SHong Zhang PetscInt *bi,*cols,nnz,*cols_lvl; 34616a2bf60SHong Zhang PetscInt *bdiag,prow,fm,nzbd,reallocs=0,dcount=0; 34716a2bf60SHong Zhang PetscInt i,levels,diagonal_fill; 348ace3abfcSBarry Smith PetscBool col_identity,row_identity,both_identity; 34916a2bf60SHong Zhang PetscReal f; 3500298fd71SBarry Smith PetscInt nlnk,*lnk,*lnk_lvl=NULL; 35116a2bf60SHong Zhang PetscBT lnkbt; 35216a2bf60SHong Zhang PetscInt nzi,*bj,**bj_ptr,**bjlvl_ptr; 3530298fd71SBarry Smith PetscFreeSpaceList free_space =NULL,current_space=NULL; 3540298fd71SBarry Smith PetscFreeSpaceList free_space_lvl=NULL,current_space_lvl=NULL; 355ace3abfcSBarry Smith PetscBool missing; 3567fa3a6a0SHong Zhang PetscInt bs=A->rmap->bs,bs2=a->bs2; 35716a2bf60SHong Zhang 35816a2bf60SHong Zhang PetscFunctionBegin; 359e32f2f54SBarry 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); 3606ba06ab7SHong Zhang if (bs>1) { /* check shifttype */ 3616c4ed002SBarry 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"); 3626ba06ab7SHong Zhang } 3636ba06ab7SHong Zhang 36416a2bf60SHong Zhang ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr); 365e32f2f54SBarry Smith if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d); 36616a2bf60SHong Zhang 36716a2bf60SHong Zhang f = info->fill; 36816a2bf60SHong Zhang levels = (PetscInt)info->levels; 36916a2bf60SHong Zhang diagonal_fill = (PetscInt)info->diagonal_fill; 37026fbe8dcSKarl Rupp 37116a2bf60SHong Zhang ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 37216a2bf60SHong Zhang 37316a2bf60SHong Zhang ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 37416a2bf60SHong Zhang ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr); 37526fbe8dcSKarl Rupp 376ace3abfcSBarry Smith both_identity = (PetscBool) (row_identity && col_identity); 37716a2bf60SHong Zhang 3787fa3a6a0SHong Zhang if (!levels && both_identity) { 37916a2bf60SHong Zhang /* special case: ilu(0) with natural ordering */ 3804dd39f65SShri Abhyankar ierr = MatILUFactorSymbolic_SeqBAIJ_ilu0(fact,A,isrow,iscol,info);CHKERRQ(ierr); 3814dd39f65SShri Abhyankar ierr = MatSeqBAIJSetNumericFactorization(fact,both_identity);CHKERRQ(ierr); 38235aa4fcfSShri Abhyankar 383d5f3da31SBarry Smith fact->factortype = MAT_FACTOR_ILU; 38435aa4fcfSShri Abhyankar (fact)->info.factor_mallocs = 0; 38535aa4fcfSShri Abhyankar (fact)->info.fill_ratio_given = info->fill; 38635aa4fcfSShri Abhyankar (fact)->info.fill_ratio_needed = 1.0; 38726fbe8dcSKarl Rupp 38835aa4fcfSShri Abhyankar b = (Mat_SeqBAIJ*)(fact)->data; 38935aa4fcfSShri Abhyankar b->row = isrow; 39035aa4fcfSShri Abhyankar b->col = iscol; 39135aa4fcfSShri Abhyankar b->icol = isicol; 39235aa4fcfSShri Abhyankar ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 39335aa4fcfSShri Abhyankar ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 39435aa4fcfSShri Abhyankar b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE; 39526fbe8dcSKarl Rupp 396785e854fSJed Brown ierr = PetscMalloc1((n+1)*bs,&b->solve_work);CHKERRQ(ierr); 39735aa4fcfSShri Abhyankar PetscFunctionReturn(0); 39835aa4fcfSShri Abhyankar } 39935aa4fcfSShri Abhyankar 40035aa4fcfSShri Abhyankar ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 40135aa4fcfSShri Abhyankar ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 40235aa4fcfSShri Abhyankar 40335aa4fcfSShri Abhyankar /* get new row pointers */ 404854ce69bSBarry Smith ierr = PetscMalloc1(n+1,&bi);CHKERRQ(ierr); 40535aa4fcfSShri Abhyankar bi[0] = 0; 40635aa4fcfSShri Abhyankar /* bdiag is location of diagonal in factor */ 407854ce69bSBarry Smith ierr = PetscMalloc1(n+1,&bdiag);CHKERRQ(ierr); 40835aa4fcfSShri Abhyankar bdiag[0] = 0; 40935aa4fcfSShri Abhyankar 410dcca6d9dSJed Brown ierr = PetscMalloc2(n,&bj_ptr,n,&bjlvl_ptr);CHKERRQ(ierr); 41135aa4fcfSShri Abhyankar 41235aa4fcfSShri Abhyankar /* create a linked list for storing column indices of the active row */ 41335aa4fcfSShri Abhyankar nlnk = n + 1; 41435aa4fcfSShri Abhyankar ierr = PetscIncompleteLLCreate(n,n,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 41535aa4fcfSShri Abhyankar 41635aa4fcfSShri Abhyankar /* initial FreeSpace size is f*(ai[n]+1) */ 417f91af8c7SBarry Smith ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(f,ai[n]+1),&free_space);CHKERRQ(ierr); 41835aa4fcfSShri Abhyankar current_space = free_space; 419f91af8c7SBarry Smith ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(f,ai[n]+1),&free_space_lvl);CHKERRQ(ierr); 42035aa4fcfSShri Abhyankar current_space_lvl = free_space_lvl; 42135aa4fcfSShri Abhyankar 42235aa4fcfSShri Abhyankar for (i=0; i<n; i++) { 42335aa4fcfSShri Abhyankar nzi = 0; 42435aa4fcfSShri Abhyankar /* copy current row into linked list */ 42535aa4fcfSShri Abhyankar nnz = ai[r[i]+1] - ai[r[i]]; 426e32f2f54SBarry 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); 42735aa4fcfSShri Abhyankar cols = aj + ai[r[i]]; 42835aa4fcfSShri Abhyankar lnk[i] = -1; /* marker to indicate if diagonal exists */ 42935aa4fcfSShri Abhyankar ierr = PetscIncompleteLLInit(nnz,cols,n,ic,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 43035aa4fcfSShri Abhyankar nzi += nlnk; 43135aa4fcfSShri Abhyankar 43235aa4fcfSShri Abhyankar /* make sure diagonal entry is included */ 43335aa4fcfSShri Abhyankar if (diagonal_fill && lnk[i] == -1) { 43435aa4fcfSShri Abhyankar fm = n; 43535aa4fcfSShri Abhyankar while (lnk[fm] < i) fm = lnk[fm]; 43635aa4fcfSShri Abhyankar lnk[i] = lnk[fm]; /* insert diagonal into linked list */ 43735aa4fcfSShri Abhyankar lnk[fm] = i; 43835aa4fcfSShri Abhyankar lnk_lvl[i] = 0; 43935aa4fcfSShri Abhyankar nzi++; dcount++; 44035aa4fcfSShri Abhyankar } 44135aa4fcfSShri Abhyankar 44235aa4fcfSShri Abhyankar /* add pivot rows into the active row */ 44335aa4fcfSShri Abhyankar nzbd = 0; 44435aa4fcfSShri Abhyankar prow = lnk[n]; 44535aa4fcfSShri Abhyankar while (prow < i) { 44635aa4fcfSShri Abhyankar nnz = bdiag[prow]; 44735aa4fcfSShri Abhyankar cols = bj_ptr[prow] + nnz + 1; 44835aa4fcfSShri Abhyankar cols_lvl = bjlvl_ptr[prow] + nnz + 1; 44935aa4fcfSShri Abhyankar nnz = bi[prow+1] - bi[prow] - nnz - 1; 45026fbe8dcSKarl Rupp 45135aa4fcfSShri Abhyankar ierr = PetscILULLAddSorted(nnz,cols,levels,cols_lvl,prow,nlnk,lnk,lnk_lvl,lnkbt,prow);CHKERRQ(ierr); 45235aa4fcfSShri Abhyankar nzi += nlnk; 45335aa4fcfSShri Abhyankar prow = lnk[prow]; 45435aa4fcfSShri Abhyankar nzbd++; 45535aa4fcfSShri Abhyankar } 45635aa4fcfSShri Abhyankar bdiag[i] = nzbd; 45735aa4fcfSShri Abhyankar bi[i+1] = bi[i] + nzi; 45835aa4fcfSShri Abhyankar 45935aa4fcfSShri Abhyankar /* if free space is not available, make more free space */ 46035aa4fcfSShri Abhyankar if (current_space->local_remaining<nzi) { 461f91af8c7SBarry Smith nnz = PetscIntMultTruncate(2,PetscIntMultTruncate(nzi,(n - i))); /* estimated and max additional space needed */ 46235aa4fcfSShri Abhyankar ierr = PetscFreeSpaceGet(nnz,¤t_space);CHKERRQ(ierr); 46335aa4fcfSShri Abhyankar ierr = PetscFreeSpaceGet(nnz,¤t_space_lvl);CHKERRQ(ierr); 46435aa4fcfSShri Abhyankar reallocs++; 46535aa4fcfSShri Abhyankar } 46635aa4fcfSShri Abhyankar 46735aa4fcfSShri Abhyankar /* copy data into free_space and free_space_lvl, then initialize lnk */ 46835aa4fcfSShri Abhyankar ierr = PetscIncompleteLLClean(n,n,nzi,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr); 46926fbe8dcSKarl Rupp 47035aa4fcfSShri Abhyankar bj_ptr[i] = current_space->array; 47135aa4fcfSShri Abhyankar bjlvl_ptr[i] = current_space_lvl->array; 47235aa4fcfSShri Abhyankar 47335aa4fcfSShri Abhyankar /* make sure the active row i has diagonal entry */ 47465e19b50SBarry 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); 47535aa4fcfSShri Abhyankar 47635aa4fcfSShri Abhyankar current_space->array += nzi; 47735aa4fcfSShri Abhyankar current_space->local_used += nzi; 47835aa4fcfSShri Abhyankar current_space->local_remaining -= nzi; 47926fbe8dcSKarl Rupp 48035aa4fcfSShri Abhyankar current_space_lvl->array += nzi; 48135aa4fcfSShri Abhyankar current_space_lvl->local_used += nzi; 48235aa4fcfSShri Abhyankar current_space_lvl->local_remaining -= nzi; 48335aa4fcfSShri Abhyankar } 48435aa4fcfSShri Abhyankar 48535aa4fcfSShri Abhyankar ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 48635aa4fcfSShri Abhyankar ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 48735aa4fcfSShri Abhyankar 48835aa4fcfSShri Abhyankar /* copy free_space into bj and free free_space; set bi, bj, bdiag in new datastructure; */ 489854ce69bSBarry Smith ierr = PetscMalloc1(bi[n]+1,&bj);CHKERRQ(ierr); 4902ce24eb6SHong Zhang ierr = PetscFreeSpaceContiguous_LU(&free_space,bj,n,bi,bdiag);CHKERRQ(ierr); 49135aa4fcfSShri Abhyankar 49235aa4fcfSShri Abhyankar ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 49335aa4fcfSShri Abhyankar ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr); 494fca92195SBarry Smith ierr = PetscFree2(bj_ptr,bjlvl_ptr);CHKERRQ(ierr); 49535aa4fcfSShri Abhyankar 49635aa4fcfSShri Abhyankar #if defined(PETSC_USE_INFO) 49735aa4fcfSShri Abhyankar { 498aef85c9fSShri Abhyankar PetscReal af = ((PetscReal)(bdiag[0]+1))/((PetscReal)ai[n]); 49957622a8eSBarry Smith ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)f,(double)af);CHKERRQ(ierr); 50057622a8eSBarry Smith ierr = PetscInfo1(A,"Run with -[sub_]pc_factor_fill %g or use \n",(double)af);CHKERRQ(ierr); 50157622a8eSBarry Smith ierr = PetscInfo1(A,"PCFactorSetFill([sub]pc,%g);\n",(double)af);CHKERRQ(ierr); 50235aa4fcfSShri Abhyankar ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr); 50335aa4fcfSShri Abhyankar if (diagonal_fill) { 504955c1f14SBarry Smith ierr = PetscInfo1(A,"Detected and replaced %D missing diagonals\n",dcount);CHKERRQ(ierr); 50535aa4fcfSShri Abhyankar } 50635aa4fcfSShri Abhyankar } 50735aa4fcfSShri Abhyankar #endif 50835aa4fcfSShri Abhyankar 50935aa4fcfSShri Abhyankar /* put together the new matrix */ 510367daffbSBarry Smith ierr = MatSeqBAIJSetPreallocation(fact,bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr); 5113bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)fact,(PetscObject)isicol);CHKERRQ(ierr); 51226fbe8dcSKarl Rupp 51335aa4fcfSShri Abhyankar b = (Mat_SeqBAIJ*)(fact)->data; 51435aa4fcfSShri Abhyankar b->free_a = PETSC_TRUE; 51535aa4fcfSShri Abhyankar b->free_ij = PETSC_TRUE; 51635aa4fcfSShri Abhyankar b->singlemalloc = PETSC_FALSE; 51726fbe8dcSKarl Rupp 518854ce69bSBarry Smith ierr = PetscMalloc1(bs2*(bdiag[0]+1),&b->a);CHKERRQ(ierr); 51926fbe8dcSKarl Rupp 52035aa4fcfSShri Abhyankar b->j = bj; 52135aa4fcfSShri Abhyankar b->i = bi; 52235aa4fcfSShri Abhyankar b->diag = bdiag; 52335aa4fcfSShri Abhyankar b->free_diag = PETSC_TRUE; 52435aa4fcfSShri Abhyankar b->ilen = 0; 52535aa4fcfSShri Abhyankar b->imax = 0; 52635aa4fcfSShri Abhyankar b->row = isrow; 52735aa4fcfSShri Abhyankar b->col = iscol; 52835aa4fcfSShri Abhyankar ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 52935aa4fcfSShri Abhyankar ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 53035aa4fcfSShri Abhyankar b->icol = isicol; 53126fbe8dcSKarl Rupp 532854ce69bSBarry Smith ierr = PetscMalloc1(bs*n+bs,&b->solve_work);CHKERRQ(ierr); 53335aa4fcfSShri Abhyankar /* In b structure: Free imax, ilen, old a, old j. 53435aa4fcfSShri Abhyankar Allocate bdiag, solve_work, new a, new j */ 5353bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)fact,(bdiag[0]+1) * (sizeof(PetscInt)+bs2*sizeof(PetscScalar)));CHKERRQ(ierr); 53635aa4fcfSShri Abhyankar b->maxnz = b->nz = bdiag[0]+1; 53726fbe8dcSKarl Rupp 538ae3d28f0SHong Zhang fact->info.factor_mallocs = reallocs; 539ae3d28f0SHong Zhang fact->info.fill_ratio_given = f; 540ae3d28f0SHong Zhang fact->info.fill_ratio_needed = ((PetscReal)(bdiag[0]+1))/((PetscReal)ai[n]); 54126fbe8dcSKarl Rupp 5424dd39f65SShri Abhyankar ierr = MatSeqBAIJSetNumericFactorization(fact,both_identity);CHKERRQ(ierr); 54335aa4fcfSShri Abhyankar PetscFunctionReturn(0); 54435aa4fcfSShri Abhyankar } 54535aa4fcfSShri Abhyankar 5464e2b4712SSatish Balay /* 5474e2b4712SSatish Balay This code is virtually identical to MatILUFactorSymbolic_SeqAIJ 5484e2b4712SSatish Balay except that the data structure of Mat_SeqAIJ is slightly different. 5494e2b4712SSatish Balay Not a good example of code reuse. 5504e2b4712SSatish Balay */ 55106e38f1dSHong Zhang PetscErrorCode MatILUFactorSymbolic_SeqBAIJ_inplace(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info) 5524e2b4712SSatish Balay { 5534e2b4712SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b; 5544e2b4712SSatish Balay IS isicol; 5556849ba73SBarry Smith PetscErrorCode ierr; 5565d0c19d7SBarry Smith const PetscInt *r,*ic,*ai = a->i,*aj = a->j,*xi; 5575d0c19d7SBarry Smith PetscInt prow,n = a->mbs,*ainew,*ajnew,jmax,*fill,nz,*im,*ajfill,*flev,*xitmp; 558a96a251dSBarry Smith PetscInt *dloc,idx,row,m,fm,nzf,nzi,reallocate = 0,dcount = 0; 559d0f46423SBarry Smith PetscInt incrlev,nnz,i,bs = A->rmap->bs,bs2 = a->bs2,levels,diagonal_fill,dd; 560ace3abfcSBarry Smith PetscBool col_identity,row_identity,both_identity,flg; 561329f5518SBarry Smith PetscReal f; 5624e2b4712SSatish Balay 5634e2b4712SSatish Balay PetscFunctionBegin; 5646bce7ff8SHong Zhang ierr = MatMissingDiagonal_SeqBAIJ(A,&flg,&dd);CHKERRQ(ierr); 565e32f2f54SBarry Smith if (flg) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix A is missing diagonal entry in row %D",dd); 5666bce7ff8SHong Zhang 567435faa5fSBarry Smith f = info->fill; 568690b6cddSBarry Smith levels = (PetscInt)info->levels; 569690b6cddSBarry Smith diagonal_fill = (PetscInt)info->diagonal_fill; 57026fbe8dcSKarl Rupp 5714c49b128SBarry Smith ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 57216a2bf60SHong Zhang 573667159a5SBarry Smith ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 574667159a5SBarry Smith ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr); 575ace3abfcSBarry Smith both_identity = (PetscBool) (row_identity && col_identity); 576309c388cSBarry Smith 57741df41f0SMatthew Knepley if (!levels && both_identity) { /* special case copy the nonzero structure */ 57816a2bf60SHong Zhang ierr = MatDuplicateNoCreate_SeqBAIJ(fact,A,MAT_DO_NOT_COPY_VALUES,PETSC_TRUE);CHKERRQ(ierr); 5798b1456e3SHong Zhang ierr = MatSeqBAIJSetNumericFactorization_inplace(fact,both_identity);CHKERRQ(ierr); 5806bce7ff8SHong Zhang 581d5f3da31SBarry Smith fact->factortype = MAT_FACTOR_ILU; 582ae3d28f0SHong Zhang b = (Mat_SeqBAIJ*)fact->data; 583bb3d539aSBarry Smith b->row = isrow; 584bb3d539aSBarry Smith b->col = iscol; 585bb3d539aSBarry Smith ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 586bb3d539aSBarry Smith ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 587bb3d539aSBarry Smith b->icol = isicol; 588bcd9e38bSBarry Smith b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE; 58926fbe8dcSKarl Rupp 590785e854fSJed Brown ierr = PetscMalloc1((n+1)*bs,&b->solve_work);CHKERRQ(ierr); 5916bce7ff8SHong Zhang PetscFunctionReturn(0); 5926bce7ff8SHong Zhang } 5936bce7ff8SHong Zhang 5946bce7ff8SHong Zhang /* general case perform the symbolic factorization */ 5954e2b4712SSatish Balay ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 5964e2b4712SSatish Balay ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 5974e2b4712SSatish Balay 5984e2b4712SSatish Balay /* get new row pointers */ 599854ce69bSBarry Smith ierr = PetscMalloc1(n+1,&ainew);CHKERRQ(ierr); 6004e2b4712SSatish Balay ainew[0] = 0; 6014e2b4712SSatish Balay /* don't know how many column pointers are needed so estimate */ 602690b6cddSBarry Smith jmax = (PetscInt)(f*ai[n] + 1); 603854ce69bSBarry Smith ierr = PetscMalloc1(jmax,&ajnew);CHKERRQ(ierr); 6044e2b4712SSatish Balay /* ajfill is level of fill for each fill entry */ 605854ce69bSBarry Smith ierr = PetscMalloc1(jmax,&ajfill);CHKERRQ(ierr); 6064e2b4712SSatish Balay /* fill is a linked list of nonzeros in active row */ 607854ce69bSBarry Smith ierr = PetscMalloc1(n+1,&fill);CHKERRQ(ierr); 6084e2b4712SSatish Balay /* im is level for each filled value */ 609854ce69bSBarry Smith ierr = PetscMalloc1(n+1,&im);CHKERRQ(ierr); 6104e2b4712SSatish Balay /* dloc is location of diagonal in factor */ 611854ce69bSBarry Smith ierr = PetscMalloc1(n+1,&dloc);CHKERRQ(ierr); 6124e2b4712SSatish Balay dloc[0] = 0; 6134e2b4712SSatish Balay for (prow=0; prow<n; prow++) { 614435faa5fSBarry Smith 615435faa5fSBarry Smith /* copy prow into linked list */ 6164e2b4712SSatish Balay nzf = nz = ai[r[prow]+1] - ai[r[prow]]; 617e32f2f54SBarry 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); 6184e2b4712SSatish Balay xi = aj + ai[r[prow]]; 6194e2b4712SSatish Balay fill[n] = n; 620435faa5fSBarry Smith fill[prow] = -1; /* marker for diagonal entry */ 6214e2b4712SSatish Balay while (nz--) { 6224e2b4712SSatish Balay fm = n; 6234e2b4712SSatish Balay idx = ic[*xi++]; 6244e2b4712SSatish Balay do { 6254e2b4712SSatish Balay m = fm; 6264e2b4712SSatish Balay fm = fill[m]; 6274e2b4712SSatish Balay } while (fm < idx); 6284e2b4712SSatish Balay fill[m] = idx; 6294e2b4712SSatish Balay fill[idx] = fm; 6304e2b4712SSatish Balay im[idx] = 0; 6314e2b4712SSatish Balay } 632435faa5fSBarry Smith 633435faa5fSBarry Smith /* make sure diagonal entry is included */ 634435faa5fSBarry Smith if (diagonal_fill && fill[prow] == -1) { 635435faa5fSBarry Smith fm = n; 636435faa5fSBarry Smith while (fill[fm] < prow) fm = fill[fm]; 637435faa5fSBarry Smith fill[prow] = fill[fm]; /* insert diagonal into linked list */ 638435faa5fSBarry Smith fill[fm] = prow; 639435faa5fSBarry Smith im[prow] = 0; 640435faa5fSBarry Smith nzf++; 641335d9088SBarry Smith dcount++; 642435faa5fSBarry Smith } 643435faa5fSBarry Smith 6444e2b4712SSatish Balay nzi = 0; 6454e2b4712SSatish Balay row = fill[n]; 6464e2b4712SSatish Balay while (row < prow) { 6474e2b4712SSatish Balay incrlev = im[row] + 1; 6484e2b4712SSatish Balay nz = dloc[row]; 649435faa5fSBarry Smith xi = ajnew + ainew[row] + nz + 1; 6504e2b4712SSatish Balay flev = ajfill + ainew[row] + nz + 1; 6514e2b4712SSatish Balay nnz = ainew[row+1] - ainew[row] - nz - 1; 6524e2b4712SSatish Balay fm = row; 6534e2b4712SSatish Balay while (nnz-- > 0) { 6544e2b4712SSatish Balay idx = *xi++; 6554e2b4712SSatish Balay if (*flev + incrlev > levels) { 6564e2b4712SSatish Balay flev++; 6574e2b4712SSatish Balay continue; 6584e2b4712SSatish Balay } 6594e2b4712SSatish Balay do { 6604e2b4712SSatish Balay m = fm; 6614e2b4712SSatish Balay fm = fill[m]; 6624e2b4712SSatish Balay } while (fm < idx); 6634e2b4712SSatish Balay if (fm != idx) { 6644e2b4712SSatish Balay im[idx] = *flev + incrlev; 6654e2b4712SSatish Balay fill[m] = idx; 6664e2b4712SSatish Balay fill[idx] = fm; 6674e2b4712SSatish Balay fm = idx; 6684e2b4712SSatish Balay nzf++; 66926fbe8dcSKarl Rupp } else if (im[idx] > *flev + incrlev) im[idx] = *flev+incrlev; 6704e2b4712SSatish Balay flev++; 6714e2b4712SSatish Balay } 6724e2b4712SSatish Balay row = fill[row]; 6734e2b4712SSatish Balay nzi++; 6744e2b4712SSatish Balay } 6754e2b4712SSatish Balay /* copy new filled row into permanent storage */ 6764e2b4712SSatish Balay ainew[prow+1] = ainew[prow] + nzf; 6774e2b4712SSatish Balay if (ainew[prow+1] > jmax) { 678ecf371e4SBarry Smith 679ecf371e4SBarry Smith /* estimate how much additional space we will need */ 680ecf371e4SBarry Smith /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */ 681ecf371e4SBarry Smith /* just double the memory each time */ 682690b6cddSBarry Smith PetscInt maxadd = jmax; 683ecf371e4SBarry Smith /* maxadd = (int)(((f*ai[n]+1)*(n-prow+5))/n); */ 6844e2b4712SSatish Balay if (maxadd < nzf) maxadd = (n-prow)*(nzf+1); 6854e2b4712SSatish Balay jmax += maxadd; 686ecf371e4SBarry Smith 687ecf371e4SBarry Smith /* allocate a longer ajnew and ajfill */ 688785e854fSJed Brown ierr = PetscMalloc1(jmax,&xitmp);CHKERRQ(ierr); 6895d0c19d7SBarry Smith ierr = PetscMemcpy(xitmp,ajnew,ainew[prow]*sizeof(PetscInt));CHKERRQ(ierr); 690606d414cSSatish Balay ierr = PetscFree(ajnew);CHKERRQ(ierr); 6915d0c19d7SBarry Smith ajnew = xitmp; 692785e854fSJed Brown ierr = PetscMalloc1(jmax,&xitmp);CHKERRQ(ierr); 6935d0c19d7SBarry Smith ierr = PetscMemcpy(xitmp,ajfill,ainew[prow]*sizeof(PetscInt));CHKERRQ(ierr); 694606d414cSSatish Balay ierr = PetscFree(ajfill);CHKERRQ(ierr); 6955d0c19d7SBarry Smith ajfill = xitmp; 696eb150c5cSKris Buschelman reallocate++; /* count how many reallocations are needed */ 6974e2b4712SSatish Balay } 6985d0c19d7SBarry Smith xitmp = ajnew + ainew[prow]; 6994e2b4712SSatish Balay flev = ajfill + ainew[prow]; 7004e2b4712SSatish Balay dloc[prow] = nzi; 7014e2b4712SSatish Balay fm = fill[n]; 7024e2b4712SSatish Balay while (nzf--) { 7035d0c19d7SBarry Smith *xitmp++ = fm; 7044e2b4712SSatish Balay *flev++ = im[fm]; 7054e2b4712SSatish Balay fm = fill[fm]; 7064e2b4712SSatish Balay } 707435faa5fSBarry Smith /* make sure row has diagonal entry */ 708f23aa3ddSBarry 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\ 7092401956bSBarry Smith try running with -pc_factor_nonzeros_along_diagonal or -pc_factor_diagonal_fill",prow); 710435faa5fSBarry Smith } 711606d414cSSatish Balay ierr = PetscFree(ajfill);CHKERRQ(ierr); 7124e2b4712SSatish Balay ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 7134e2b4712SSatish Balay ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 714606d414cSSatish Balay ierr = PetscFree(fill);CHKERRQ(ierr); 715606d414cSSatish Balay ierr = PetscFree(im);CHKERRQ(ierr); 7164e2b4712SSatish Balay 7176cf91177SBarry Smith #if defined(PETSC_USE_INFO) 7184e2b4712SSatish Balay { 719329f5518SBarry Smith PetscReal af = ((PetscReal)ainew[n])/((PetscReal)ai[n]); 72057622a8eSBarry Smith ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocate,(double)f,(double)af);CHKERRQ(ierr); 72157622a8eSBarry Smith ierr = PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);CHKERRQ(ierr); 72257622a8eSBarry Smith ierr = PetscInfo1(A,"PCFactorSetFill(pc,%g);\n",(double)af);CHKERRQ(ierr); 723ae15b995SBarry Smith ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr); 724335d9088SBarry Smith if (diagonal_fill) { 725ae15b995SBarry Smith ierr = PetscInfo1(A,"Detected and replaced %D missing diagonals\n",dcount);CHKERRQ(ierr); 726335d9088SBarry Smith } 7274e2b4712SSatish Balay } 72863ba0a88SBarry Smith #endif 7294e2b4712SSatish Balay 7304e2b4712SSatish Balay /* put together the new matrix */ 731367daffbSBarry Smith ierr = MatSeqBAIJSetPreallocation(fact,bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr); 7323bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)fact,(PetscObject)isicol);CHKERRQ(ierr); 733ae3d28f0SHong Zhang b = (Mat_SeqBAIJ*)fact->data; 73426fbe8dcSKarl Rupp 735e6b907acSBarry Smith b->free_a = PETSC_TRUE; 736e6b907acSBarry Smith b->free_ij = PETSC_TRUE; 7377c922b88SBarry Smith b->singlemalloc = PETSC_FALSE; 73826fbe8dcSKarl Rupp 739785e854fSJed Brown ierr = PetscMalloc1(bs2*ainew[n],&b->a);CHKERRQ(ierr); 74026fbe8dcSKarl Rupp 7414e2b4712SSatish Balay b->j = ajnew; 7424e2b4712SSatish Balay b->i = ainew; 7434e2b4712SSatish Balay for (i=0; i<n; i++) dloc[i] += ainew[i]; 7444e2b4712SSatish Balay b->diag = dloc; 7457f53bb6cSHong Zhang b->free_diag = PETSC_TRUE; 7464e2b4712SSatish Balay b->ilen = 0; 7474e2b4712SSatish Balay b->imax = 0; 7484e2b4712SSatish Balay b->row = isrow; 7494e2b4712SSatish Balay b->col = iscol; 750bcd9e38bSBarry Smith b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE; 75126fbe8dcSKarl Rupp 752c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 753c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 754e51c0b9cSSatish Balay b->icol = isicol; 755854ce69bSBarry Smith ierr = PetscMalloc1(bs*n+bs,&b->solve_work);CHKERRQ(ierr); 7564e2b4712SSatish Balay /* In b structure: Free imax, ilen, old a, old j. 7574e2b4712SSatish Balay Allocate dloc, solve_work, new a, new j */ 7583bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)fact,(ainew[n]-n)*(sizeof(PetscInt))+bs2*ainew[n]*sizeof(PetscScalar));CHKERRQ(ierr); 7594e2b4712SSatish Balay b->maxnz = b->nz = ainew[n]; 7604e2b4712SSatish Balay 761ae3d28f0SHong Zhang fact->info.factor_mallocs = reallocate; 762ae3d28f0SHong Zhang fact->info.fill_ratio_given = f; 763ae3d28f0SHong Zhang fact->info.fill_ratio_needed = ((PetscReal)ainew[n])/((PetscReal)ai[prow]); 7646bce7ff8SHong Zhang 7658b1456e3SHong Zhang ierr = MatSeqBAIJSetNumericFactorization_inplace(fact,both_identity);CHKERRQ(ierr); 7668661488fSKris Buschelman PetscFunctionReturn(0); 7678661488fSKris Buschelman } 7688661488fSKris Buschelman 769dfbe8321SBarry Smith PetscErrorCode MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE(Mat A) 7707e7071cdSKris Buschelman { 77112272027SHong Zhang /* Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; */ 77212272027SHong Zhang /* int i,*AJ=a->j,nz=a->nz; */ 7735fd66863SKarl Rupp 7745a9542e3SKris Buschelman PetscFunctionBegin; 7757cf1b8d3SKris Buschelman /* Undo Column scaling */ 7767cf1b8d3SKris Buschelman /* while (nz--) { */ 7777cf1b8d3SKris Buschelman /* AJ[i] = AJ[i]/4; */ 7787cf1b8d3SKris Buschelman /* } */ 779c115a38dSKris Buschelman /* This should really invoke a push/pop logic, but we don't have that yet. */ 7800298fd71SBarry Smith A->ops->setunfactored = NULL; 7817cf1b8d3SKris Buschelman PetscFunctionReturn(0); 7827cf1b8d3SKris Buschelman } 7837cf1b8d3SKris Buschelman 784dfbe8321SBarry Smith PetscErrorCode MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE_usj(Mat A) 7857cf1b8d3SKris Buschelman { 7867cf1b8d3SKris Buschelman Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 787b24ad042SBarry Smith PetscInt *AJ=a->j,nz=a->nz; 7882aa5897fSKris Buschelman unsigned short *aj=(unsigned short*)AJ; 7895fd66863SKarl Rupp 7905a9542e3SKris Buschelman PetscFunctionBegin; 7910b9da03eSKris Buschelman /* Is this really necessary? */ 79220235379SKris Buschelman while (nz--) { 7930b9da03eSKris Buschelman AJ[nz] = (int)((unsigned int)aj[nz]); /* First extend, then convert to signed. */ 7947e7071cdSKris Buschelman } 7950298fd71SBarry Smith A->ops->setunfactored = NULL; 7967e7071cdSKris Buschelman PetscFunctionReturn(0); 7977e7071cdSKris Buschelman } 7987e7071cdSKris Buschelman 799732ee342SKris Buschelman 800