1be1d678aSKris Buschelman 24e2b4712SSatish Balay /* 34e2b4712SSatish Balay Factorization code for BAIJ format. 44e2b4712SSatish Balay */ 54e2b4712SSatish Balay 6c6db04a5SJed Brown #include <../src/mat/impls/baij/seq/baij.h> 7af0996ceSBarry Smith #include <petsc/private/kernels/blockinvert.h> 8c6db04a5SJed Brown #include <petscbt.h> 9c6db04a5SJed Brown #include <../src/mat/utils/freespace.h> 104e2b4712SSatish Balay 114e2b4712SSatish Balay /* ----------------------------------------------------------------*/ 1209573ac7SBarry Smith extern PetscErrorCode MatDuplicateNoCreate_SeqBAIJ(Mat,Mat,MatDuplicateOption,PetscBool); 136bce7ff8SHong Zhang 142b0b2ea7SShri Abhyankar #undef __FUNCT__ 1529a97285SShri Abhyankar #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_15_NaturalOrdering" 16766f9fbaSBarry Smith /* 17766f9fbaSBarry Smith This is not much faster than MatLUFactorNumeric_SeqBAIJ_N() but the solve is faster at least sometimes 18766f9fbaSBarry Smith */ 1929a97285SShri Abhyankar PetscErrorCode MatLUFactorNumeric_SeqBAIJ_15_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info) 202b0b2ea7SShri Abhyankar { 212b0b2ea7SShri Abhyankar Mat C =B; 222b0b2ea7SShri Abhyankar Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data; 232b0b2ea7SShri Abhyankar PetscErrorCode ierr; 24766f9fbaSBarry Smith PetscInt i,j,k,ipvt[15]; 25766f9fbaSBarry Smith const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*ajtmp,*bjtmp,*bdiag=b->diag,*pj; 26766f9fbaSBarry Smith PetscInt nz,nzL,row; 27766f9fbaSBarry Smith MatScalar *rtmp,*pc,*mwork,*pv,*vv,work[225]; 28766f9fbaSBarry Smith const MatScalar *v,*aa=a->a; 292b0b2ea7SShri Abhyankar PetscInt bs2 = a->bs2,bs=A->rmap->bs,flg; 300fa040f9SShri Abhyankar PetscInt sol_ver; 31a455e926SHong Zhang PetscBool allowzeropivot,zeropivotdetected; 322b0b2ea7SShri Abhyankar 332b0b2ea7SShri Abhyankar PetscFunctionBegin; 340164db54SHong Zhang allowzeropivot = PetscNot(A->erroriffailure); 35c5929fdfSBarry Smith ierr = PetscOptionsGetInt(NULL,((PetscObject)A)->prefix,"-sol_ver",&sol_ver,NULL);CHKERRQ(ierr); 360fa040f9SShri Abhyankar 372b0b2ea7SShri Abhyankar /* generate work space needed by the factorization */ 38dcca6d9dSJed Brown ierr = PetscMalloc2(bs2*n,&rtmp,bs2,&mwork);CHKERRQ(ierr); 392b0b2ea7SShri Abhyankar ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr); 402b0b2ea7SShri Abhyankar 412b0b2ea7SShri Abhyankar for (i=0; i<n; i++) { 422b0b2ea7SShri Abhyankar /* zero rtmp */ 432b0b2ea7SShri Abhyankar /* L part */ 442b0b2ea7SShri Abhyankar nz = bi[i+1] - bi[i]; 452b0b2ea7SShri Abhyankar bjtmp = bj + bi[i]; 462b0b2ea7SShri Abhyankar for (j=0; j<nz; j++) { 472b0b2ea7SShri Abhyankar ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 482b0b2ea7SShri Abhyankar } 492b0b2ea7SShri Abhyankar 502b0b2ea7SShri Abhyankar /* U part */ 512b0b2ea7SShri Abhyankar nz = bdiag[i] - bdiag[i+1]; 522b0b2ea7SShri Abhyankar bjtmp = bj + bdiag[i+1]+1; 532b0b2ea7SShri Abhyankar for (j=0; j<nz; j++) { 542b0b2ea7SShri Abhyankar ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 552b0b2ea7SShri Abhyankar } 562b0b2ea7SShri Abhyankar 572b0b2ea7SShri Abhyankar /* load in initial (unfactored row) */ 5829a97285SShri Abhyankar nz = ai[i+1] - ai[i]; 5929a97285SShri Abhyankar ajtmp = aj + ai[i]; 6029a97285SShri Abhyankar v = aa + bs2*ai[i]; 612b0b2ea7SShri Abhyankar for (j=0; j<nz; j++) { 6229a97285SShri Abhyankar ierr = PetscMemcpy(rtmp+bs2*ajtmp[j],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr); 632b0b2ea7SShri Abhyankar } 642b0b2ea7SShri Abhyankar 652b0b2ea7SShri Abhyankar /* elimination */ 662b0b2ea7SShri Abhyankar bjtmp = bj + bi[i]; 672b0b2ea7SShri Abhyankar nzL = bi[i+1] - bi[i]; 682b0b2ea7SShri Abhyankar for (k=0; k < nzL; k++) { 692b0b2ea7SShri Abhyankar row = bjtmp[k]; 702b0b2ea7SShri Abhyankar pc = rtmp + bs2*row; 71c35f09e5SBarry Smith for (flg=0,j=0; j<bs2; j++) { 72c35f09e5SBarry Smith if (pc[j]!=0.0) { 73c35f09e5SBarry Smith flg = 1; 74c35f09e5SBarry Smith break; 75c35f09e5SBarry Smith } 76c35f09e5SBarry Smith } 772b0b2ea7SShri Abhyankar if (flg) { 782b0b2ea7SShri Abhyankar pv = b->a + bs2*bdiag[row]; 7996b95a6bSBarry Smith PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); 8096b95a6bSBarry Smith /*ierr = PetscKernel_A_gets_A_times_B_15(pc,pv,mwork);CHKERRQ(ierr);*/ 812b0b2ea7SShri Abhyankar pj = b->j + bdiag[row+1]+1; /* begining of U(row,:) */ 822b0b2ea7SShri Abhyankar pv = b->a + bs2*(bdiag[row+1]+1); 832b0b2ea7SShri Abhyankar nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries inU(row,:), excluding diag */ 842b0b2ea7SShri Abhyankar for (j=0; j<nz; j++) { 85766f9fbaSBarry Smith vv = rtmp + bs2*pj[j]; 8696b95a6bSBarry Smith PetscKernel_A_gets_A_minus_B_times_C(bs,vv,pc,pv); 8796b95a6bSBarry Smith /* ierr = PetscKernel_A_gets_A_minus_B_times_C_15(vv,pc,pv);CHKERRQ(ierr); */ 882b0b2ea7SShri Abhyankar pv += bs2; 892b0b2ea7SShri Abhyankar } 90766f9fbaSBarry Smith ierr = PetscLogFlops(2*bs2*bs*(nz+1)-bs2);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */ 912b0b2ea7SShri Abhyankar } 922b0b2ea7SShri Abhyankar } 932b0b2ea7SShri Abhyankar 942b0b2ea7SShri Abhyankar /* finished row so stick it into b->a */ 952b0b2ea7SShri Abhyankar /* L part */ 962b0b2ea7SShri Abhyankar pv = b->a + bs2*bi[i]; 972b0b2ea7SShri Abhyankar pj = b->j + bi[i]; 982b0b2ea7SShri Abhyankar nz = bi[i+1] - bi[i]; 992b0b2ea7SShri Abhyankar for (j=0; j<nz; j++) { 1002b0b2ea7SShri Abhyankar ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 1012b0b2ea7SShri Abhyankar } 1022b0b2ea7SShri Abhyankar 1032b0b2ea7SShri Abhyankar /* Mark diagonal and invert diagonal for simplier triangular solves */ 1042b0b2ea7SShri Abhyankar pv = b->a + bs2*bdiag[i]; 1052b0b2ea7SShri Abhyankar pj = b->j + bdiag[i]; 1062b0b2ea7SShri Abhyankar ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr); 107a455e926SHong Zhang ierr = PetscKernel_A_gets_inverse_A_15(pv,ipvt,work,info->shiftamount,allowzeropivot,&zeropivotdetected);CHKERRQ(ierr); 108*7b6c816cSBarry Smith if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 1092b0b2ea7SShri Abhyankar 1102b0b2ea7SShri Abhyankar /* U part */ 1112b0b2ea7SShri Abhyankar pv = b->a + bs2*(bdiag[i+1]+1); 1122b0b2ea7SShri Abhyankar pj = b->j + bdiag[i+1]+1; 1132b0b2ea7SShri Abhyankar nz = bdiag[i] - bdiag[i+1] - 1; 1142b0b2ea7SShri Abhyankar for (j=0; j<nz; j++) { 1152b0b2ea7SShri Abhyankar ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 1162b0b2ea7SShri Abhyankar } 1172b0b2ea7SShri Abhyankar } 1182b0b2ea7SShri Abhyankar 1192b0b2ea7SShri Abhyankar ierr = PetscFree2(rtmp,mwork);CHKERRQ(ierr); 12026fbe8dcSKarl Rupp 121832cc040SShri Abhyankar C->ops->solve = MatSolve_SeqBAIJ_15_NaturalOrdering_ver1; 122766f9fbaSBarry Smith C->ops->solvetranspose = MatSolve_SeqBAIJ_N_NaturalOrdering; 1232b0b2ea7SShri Abhyankar C->assembled = PETSC_TRUE; 12426fbe8dcSKarl Rupp 125766f9fbaSBarry Smith ierr = PetscLogFlops(1.333333333333*bs*bs2*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 1262b0b2ea7SShri Abhyankar PetscFunctionReturn(0); 1272b0b2ea7SShri Abhyankar } 1282b0b2ea7SShri Abhyankar 1296bce7ff8SHong Zhang #undef __FUNCT__ 1304dd39f65SShri Abhyankar #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_N" 1314dd39f65SShri Abhyankar PetscErrorCode MatLUFactorNumeric_SeqBAIJ_N(Mat B,Mat A,const MatFactorInfo *info) 1326bce7ff8SHong Zhang { 1336bce7ff8SHong Zhang Mat C =B; 1346bce7ff8SHong Zhang Mat_SeqBAIJ *a =(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data; 1356bce7ff8SHong Zhang IS isrow = b->row,isicol = b->icol; 1366bce7ff8SHong Zhang PetscErrorCode ierr; 1375a586d82SBarry Smith const PetscInt *r,*ic; 1386bce7ff8SHong Zhang PetscInt i,j,k,n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 1396bce7ff8SHong Zhang PetscInt *ajtmp,*bjtmp,nz,nzL,row,*bdiag=b->diag,*pj; 140b588c5a2SHong Zhang MatScalar *rtmp,*pc,*mwork,*v,*pv,*aa=a->a; 141914a18a2SHong Zhang PetscInt bs=A->rmap->bs,bs2 = a->bs2,*v_pivots,flg; 142914a18a2SHong Zhang MatScalar *v_work; 143ace3abfcSBarry Smith PetscBool col_identity,row_identity,both_identity; 1445f8bbccaSHong Zhang PetscBool allowzeropivot,zeropivotdetected; 1456bce7ff8SHong Zhang 1466bce7ff8SHong Zhang PetscFunctionBegin; 1476bce7ff8SHong Zhang ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 1486bce7ff8SHong Zhang ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 1495f8bbccaSHong Zhang allowzeropivot = PetscNot(A->erroriffailure); 150ae3d28f0SHong Zhang 151785e854fSJed Brown ierr = PetscMalloc1(bs2*n,&rtmp);CHKERRQ(ierr); 152fca92195SBarry Smith ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr); 1536bce7ff8SHong Zhang 154914a18a2SHong Zhang /* generate work space needed by dense LU factorization */ 155dcca6d9dSJed Brown ierr = PetscMalloc3(bs,&v_work,bs2,&mwork,bs,&v_pivots);CHKERRQ(ierr); 156914a18a2SHong Zhang 1576bce7ff8SHong Zhang for (i=0; i<n; i++) { 1586bce7ff8SHong Zhang /* zero rtmp */ 1596bce7ff8SHong Zhang /* L part */ 1606bce7ff8SHong Zhang nz = bi[i+1] - bi[i]; 1616bce7ff8SHong Zhang bjtmp = bj + bi[i]; 162914a18a2SHong Zhang for (j=0; j<nz; j++) { 163914a18a2SHong Zhang ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 164914a18a2SHong Zhang } 1656bce7ff8SHong Zhang 1666bce7ff8SHong Zhang /* U part */ 1671a83e813SShri Abhyankar nz = bdiag[i] - bdiag[i+1]; 1681a83e813SShri Abhyankar bjtmp = bj + bdiag[i+1]+1; 1691a83e813SShri Abhyankar for (j=0; j<nz; j++) { 1701a83e813SShri Abhyankar ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 1711a83e813SShri Abhyankar } 1721a83e813SShri Abhyankar 1731a83e813SShri Abhyankar /* load in initial (unfactored row) */ 1741a83e813SShri Abhyankar nz = ai[r[i]+1] - ai[r[i]]; 1751a83e813SShri Abhyankar ajtmp = aj + ai[r[i]]; 1761a83e813SShri Abhyankar v = aa + bs2*ai[r[i]]; 1771a83e813SShri Abhyankar for (j=0; j<nz; j++) { 1781a83e813SShri Abhyankar ierr = PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1791a83e813SShri Abhyankar } 1801a83e813SShri Abhyankar 1811a83e813SShri Abhyankar /* elimination */ 1821a83e813SShri Abhyankar bjtmp = bj + bi[i]; 1831a83e813SShri Abhyankar nzL = bi[i+1] - bi[i]; 1841a83e813SShri Abhyankar for (k=0; k < nzL; k++) { 1851a83e813SShri Abhyankar row = bjtmp[k]; 1861a83e813SShri Abhyankar pc = rtmp + bs2*row; 187c35f09e5SBarry Smith for (flg=0,j=0; j<bs2; j++) { 188c35f09e5SBarry Smith if (pc[j]!=0.0) { 189c35f09e5SBarry Smith flg = 1; 190c35f09e5SBarry Smith break; 191c35f09e5SBarry Smith } 192c35f09e5SBarry Smith } 1931a83e813SShri Abhyankar if (flg) { 1941a83e813SShri Abhyankar pv = b->a + bs2*bdiag[row]; 19596b95a6bSBarry Smith PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); /* *pc = *pc * (*pv); */ 1961a83e813SShri Abhyankar pj = b->j + bdiag[row+1]+1; /* begining of U(row,:) */ 1971a83e813SShri Abhyankar pv = b->a + bs2*(bdiag[row+1]+1); 1981a83e813SShri Abhyankar nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries inU(row,:), excluding diag */ 1991a83e813SShri Abhyankar for (j=0; j<nz; j++) { 20096b95a6bSBarry Smith PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); 2011a83e813SShri Abhyankar } 2021a83e813SShri Abhyankar ierr = PetscLogFlops(2*bs2*bs*(nz+1)-bs2);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */ 2031a83e813SShri Abhyankar } 2041a83e813SShri Abhyankar } 2051a83e813SShri Abhyankar 2061a83e813SShri Abhyankar /* finished row so stick it into b->a */ 2071a83e813SShri Abhyankar /* L part */ 2081a83e813SShri Abhyankar pv = b->a + bs2*bi[i]; 2091a83e813SShri Abhyankar pj = b->j + bi[i]; 2101a83e813SShri Abhyankar nz = bi[i+1] - bi[i]; 2111a83e813SShri Abhyankar for (j=0; j<nz; j++) { 2121a83e813SShri Abhyankar ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 2131a83e813SShri Abhyankar } 2141a83e813SShri Abhyankar 2151a83e813SShri Abhyankar /* Mark diagonal and invert diagonal for simplier triangular solves */ 2161a83e813SShri Abhyankar pv = b->a + bs2*bdiag[i]; 2171a83e813SShri Abhyankar pj = b->j + bdiag[i]; 218e32f2f54SBarry Smith /* if (*pj != i)SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"row %d != *pj %d",i,*pj); */ 2191a83e813SShri Abhyankar ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr); 2205f8bbccaSHong Zhang 2215f8bbccaSHong Zhang ierr = PetscKernel_A_gets_inverse_A(bs,pv,v_pivots,v_work,allowzeropivot,&zeropivotdetected);CHKERRQ(ierr); 222*7b6c816cSBarry Smith if (zeropivotdetected) B->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 2231a83e813SShri Abhyankar 2241a83e813SShri Abhyankar /* U part */ 2251a83e813SShri Abhyankar pv = b->a + bs2*(bdiag[i+1]+1); 2261a83e813SShri Abhyankar pj = b->j + bdiag[i+1]+1; 2271a83e813SShri Abhyankar nz = bdiag[i] - bdiag[i+1] - 1; 2281a83e813SShri Abhyankar for (j=0; j<nz; j++) { 2291a83e813SShri Abhyankar ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 2301a83e813SShri Abhyankar } 2311a83e813SShri Abhyankar } 2321a83e813SShri Abhyankar 2331a83e813SShri Abhyankar ierr = PetscFree(rtmp);CHKERRQ(ierr); 234fca92195SBarry Smith ierr = PetscFree3(v_work,mwork,v_pivots);CHKERRQ(ierr); 2351a83e813SShri Abhyankar ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 2361a83e813SShri Abhyankar ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 2371a83e813SShri Abhyankar 238ae3d28f0SHong Zhang ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 239ae3d28f0SHong Zhang ierr = ISIdentity(isicol,&col_identity);CHKERRQ(ierr); 24026fbe8dcSKarl Rupp 241ace3abfcSBarry Smith both_identity = (PetscBool) (row_identity && col_identity); 242ae3d28f0SHong Zhang if (both_identity) { 2434dd39f65SShri Abhyankar C->ops->solve = MatSolve_SeqBAIJ_N_NaturalOrdering; 244ae3d28f0SHong Zhang } else { 2454dd39f65SShri Abhyankar C->ops->solve = MatSolve_SeqBAIJ_N; 246ae3d28f0SHong Zhang } 2474dd39f65SShri Abhyankar C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_N; 248ae3d28f0SHong Zhang 2491a83e813SShri Abhyankar C->assembled = PETSC_TRUE; 25026fbe8dcSKarl Rupp 251766f9fbaSBarry Smith ierr = PetscLogFlops(1.333333333333*bs*bs2*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 2521a83e813SShri Abhyankar PetscFunctionReturn(0); 2531a83e813SShri Abhyankar } 2541a83e813SShri Abhyankar 2556bce7ff8SHong Zhang /* 2566bce7ff8SHong Zhang ilu(0) with natural ordering under new data structure. 2574dd39f65SShri Abhyankar See MatILUFactorSymbolic_SeqAIJ_ilu0() for detailed description 2584dd39f65SShri Abhyankar because this code is almost identical to MatILUFactorSymbolic_SeqAIJ_ilu0_inplace(). 2596bce7ff8SHong Zhang */ 260c0c7eb62SShri Abhyankar 2616bce7ff8SHong Zhang #undef __FUNCT__ 2624dd39f65SShri Abhyankar #define __FUNCT__ "MatILUFactorSymbolic_SeqBAIJ_ilu0" 2634dd39f65SShri Abhyankar PetscErrorCode MatILUFactorSymbolic_SeqBAIJ_ilu0(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info) 2646bce7ff8SHong Zhang { 2656bce7ff8SHong Zhang 2666bce7ff8SHong Zhang Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b; 2676bce7ff8SHong Zhang PetscErrorCode ierr; 26816a2bf60SHong Zhang PetscInt n=a->mbs,*ai=a->i,*aj,*adiag=a->diag,bs2 = a->bs2; 26935aa4fcfSShri Abhyankar PetscInt i,j,nz,*bi,*bj,*bdiag,bi_temp; 27035aa4fcfSShri Abhyankar 27135aa4fcfSShri Abhyankar PetscFunctionBegin; 27235aa4fcfSShri Abhyankar ierr = MatDuplicateNoCreate_SeqBAIJ(fact,A,MAT_DO_NOT_COPY_VALUES,PETSC_FALSE);CHKERRQ(ierr); 27335aa4fcfSShri Abhyankar b = (Mat_SeqBAIJ*)(fact)->data; 27435aa4fcfSShri Abhyankar 27535aa4fcfSShri Abhyankar /* allocate matrix arrays for new data structure */ 276dcca6d9dSJed Brown ierr = PetscMalloc3(bs2*ai[n]+1,&b->a,ai[n]+1,&b->j,n+1,&b->i);CHKERRQ(ierr); 2773bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)fact,ai[n]*(bs2*sizeof(PetscScalar)+sizeof(PetscInt))+(n+1)*sizeof(PetscInt));CHKERRQ(ierr); 27826fbe8dcSKarl Rupp 27935aa4fcfSShri Abhyankar b->singlemalloc = PETSC_TRUE; 280379be0ddSLisandro Dalcin b->free_a = PETSC_TRUE; 281379be0ddSLisandro Dalcin b->free_ij = PETSC_TRUE; 2821e40a84eSLisandro Dalcin fact->preallocated = PETSC_TRUE; 2831e40a84eSLisandro Dalcin fact->assembled = PETSC_TRUE; 28435aa4fcfSShri Abhyankar if (!b->diag) { 285854ce69bSBarry Smith ierr = PetscMalloc1(n+1,&b->diag);CHKERRQ(ierr); 2863bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)fact,(n+1)*sizeof(PetscInt));CHKERRQ(ierr); 28735aa4fcfSShri Abhyankar } 28835aa4fcfSShri Abhyankar bdiag = b->diag; 28935aa4fcfSShri Abhyankar 29035aa4fcfSShri Abhyankar if (n > 0) { 29135aa4fcfSShri Abhyankar ierr = PetscMemzero(b->a,bs2*ai[n]*sizeof(MatScalar));CHKERRQ(ierr); 29235aa4fcfSShri Abhyankar } 29335aa4fcfSShri Abhyankar 29435aa4fcfSShri Abhyankar /* set bi and bj with new data structure */ 29535aa4fcfSShri Abhyankar bi = b->i; 29635aa4fcfSShri Abhyankar bj = b->j; 29735aa4fcfSShri Abhyankar 29835aa4fcfSShri Abhyankar /* L part */ 29935aa4fcfSShri Abhyankar bi[0] = 0; 30035aa4fcfSShri Abhyankar for (i=0; i<n; i++) { 30135aa4fcfSShri Abhyankar nz = adiag[i] - ai[i]; 30235aa4fcfSShri Abhyankar bi[i+1] = bi[i] + nz; 30335aa4fcfSShri Abhyankar aj = a->j + ai[i]; 30435aa4fcfSShri Abhyankar for (j=0; j<nz; j++) { 30535aa4fcfSShri Abhyankar *bj = aj[j]; bj++; 30635aa4fcfSShri Abhyankar } 30735aa4fcfSShri Abhyankar } 30835aa4fcfSShri Abhyankar 30935aa4fcfSShri Abhyankar /* U part */ 31035aa4fcfSShri Abhyankar bi_temp = bi[n]; 31135aa4fcfSShri Abhyankar bdiag[n] = bi[n]-1; 31235aa4fcfSShri Abhyankar for (i=n-1; i>=0; i--) { 31335aa4fcfSShri Abhyankar nz = ai[i+1] - adiag[i] - 1; 31435aa4fcfSShri Abhyankar bi_temp = bi_temp + nz + 1; 31535aa4fcfSShri Abhyankar aj = a->j + adiag[i] + 1; 31635aa4fcfSShri Abhyankar for (j=0; j<nz; j++) { 31735aa4fcfSShri Abhyankar *bj = aj[j]; bj++; 31835aa4fcfSShri Abhyankar } 31935aa4fcfSShri Abhyankar /* diag[i] */ 32035aa4fcfSShri Abhyankar *bj = i; bj++; 32135aa4fcfSShri Abhyankar bdiag[i] = bi_temp - 1; 32235aa4fcfSShri Abhyankar } 32335aa4fcfSShri Abhyankar PetscFunctionReturn(0); 32435aa4fcfSShri Abhyankar } 32535aa4fcfSShri Abhyankar 32635aa4fcfSShri Abhyankar #undef __FUNCT__ 3274dd39f65SShri Abhyankar #define __FUNCT__ "MatILUFactorSymbolic_SeqBAIJ" 3284dd39f65SShri Abhyankar PetscErrorCode MatILUFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info) 32916a2bf60SHong Zhang { 33016a2bf60SHong Zhang Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b; 33116a2bf60SHong Zhang IS isicol; 33216a2bf60SHong Zhang PetscErrorCode ierr; 33316a2bf60SHong Zhang const PetscInt *r,*ic; 3347fa3a6a0SHong Zhang PetscInt n=a->mbs,*ai=a->i,*aj=a->j,d; 33516a2bf60SHong Zhang PetscInt *bi,*cols,nnz,*cols_lvl; 33616a2bf60SHong Zhang PetscInt *bdiag,prow,fm,nzbd,reallocs=0,dcount=0; 33716a2bf60SHong Zhang PetscInt i,levels,diagonal_fill; 338ace3abfcSBarry Smith PetscBool col_identity,row_identity,both_identity; 33916a2bf60SHong Zhang PetscReal f; 3400298fd71SBarry Smith PetscInt nlnk,*lnk,*lnk_lvl=NULL; 34116a2bf60SHong Zhang PetscBT lnkbt; 34216a2bf60SHong Zhang PetscInt nzi,*bj,**bj_ptr,**bjlvl_ptr; 3430298fd71SBarry Smith PetscFreeSpaceList free_space =NULL,current_space=NULL; 3440298fd71SBarry Smith PetscFreeSpaceList free_space_lvl=NULL,current_space_lvl=NULL; 345ace3abfcSBarry Smith PetscBool missing; 3467fa3a6a0SHong Zhang PetscInt bs=A->rmap->bs,bs2=a->bs2; 34716a2bf60SHong Zhang 34816a2bf60SHong Zhang PetscFunctionBegin; 349e32f2f54SBarry 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); 3506ba06ab7SHong Zhang if (bs>1) { /* check shifttype */ 3516c4ed002SBarry 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"); 3526ba06ab7SHong Zhang } 3536ba06ab7SHong Zhang 35416a2bf60SHong Zhang ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr); 355e32f2f54SBarry Smith if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d); 35616a2bf60SHong Zhang 35716a2bf60SHong Zhang f = info->fill; 35816a2bf60SHong Zhang levels = (PetscInt)info->levels; 35916a2bf60SHong Zhang diagonal_fill = (PetscInt)info->diagonal_fill; 36026fbe8dcSKarl Rupp 36116a2bf60SHong Zhang ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 36216a2bf60SHong Zhang 36316a2bf60SHong Zhang ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 36416a2bf60SHong Zhang ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr); 36526fbe8dcSKarl Rupp 366ace3abfcSBarry Smith both_identity = (PetscBool) (row_identity && col_identity); 36716a2bf60SHong Zhang 3687fa3a6a0SHong Zhang if (!levels && both_identity) { 36916a2bf60SHong Zhang /* special case: ilu(0) with natural ordering */ 3704dd39f65SShri Abhyankar ierr = MatILUFactorSymbolic_SeqBAIJ_ilu0(fact,A,isrow,iscol,info);CHKERRQ(ierr); 3714dd39f65SShri Abhyankar ierr = MatSeqBAIJSetNumericFactorization(fact,both_identity);CHKERRQ(ierr); 37235aa4fcfSShri Abhyankar 373d5f3da31SBarry Smith fact->factortype = MAT_FACTOR_ILU; 37435aa4fcfSShri Abhyankar (fact)->info.factor_mallocs = 0; 37535aa4fcfSShri Abhyankar (fact)->info.fill_ratio_given = info->fill; 37635aa4fcfSShri Abhyankar (fact)->info.fill_ratio_needed = 1.0; 37726fbe8dcSKarl Rupp 37835aa4fcfSShri Abhyankar b = (Mat_SeqBAIJ*)(fact)->data; 37935aa4fcfSShri Abhyankar b->row = isrow; 38035aa4fcfSShri Abhyankar b->col = iscol; 38135aa4fcfSShri Abhyankar b->icol = isicol; 38235aa4fcfSShri Abhyankar ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 38335aa4fcfSShri Abhyankar ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 38435aa4fcfSShri Abhyankar b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE; 38526fbe8dcSKarl Rupp 386785e854fSJed Brown ierr = PetscMalloc1((n+1)*bs,&b->solve_work);CHKERRQ(ierr); 38735aa4fcfSShri Abhyankar PetscFunctionReturn(0); 38835aa4fcfSShri Abhyankar } 38935aa4fcfSShri Abhyankar 39035aa4fcfSShri Abhyankar ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 39135aa4fcfSShri Abhyankar ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 39235aa4fcfSShri Abhyankar 39335aa4fcfSShri Abhyankar /* get new row pointers */ 394854ce69bSBarry Smith ierr = PetscMalloc1(n+1,&bi);CHKERRQ(ierr); 39535aa4fcfSShri Abhyankar bi[0] = 0; 39635aa4fcfSShri Abhyankar /* bdiag is location of diagonal in factor */ 397854ce69bSBarry Smith ierr = PetscMalloc1(n+1,&bdiag);CHKERRQ(ierr); 39835aa4fcfSShri Abhyankar bdiag[0] = 0; 39935aa4fcfSShri Abhyankar 400dcca6d9dSJed Brown ierr = PetscMalloc2(n,&bj_ptr,n,&bjlvl_ptr);CHKERRQ(ierr); 40135aa4fcfSShri Abhyankar 40235aa4fcfSShri Abhyankar /* create a linked list for storing column indices of the active row */ 40335aa4fcfSShri Abhyankar nlnk = n + 1; 40435aa4fcfSShri Abhyankar ierr = PetscIncompleteLLCreate(n,n,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 40535aa4fcfSShri Abhyankar 40635aa4fcfSShri Abhyankar /* initial FreeSpace size is f*(ai[n]+1) */ 407f91af8c7SBarry Smith ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(f,ai[n]+1),&free_space);CHKERRQ(ierr); 40835aa4fcfSShri Abhyankar current_space = free_space; 409f91af8c7SBarry Smith ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(f,ai[n]+1),&free_space_lvl);CHKERRQ(ierr); 41035aa4fcfSShri Abhyankar current_space_lvl = free_space_lvl; 41135aa4fcfSShri Abhyankar 41235aa4fcfSShri Abhyankar for (i=0; i<n; i++) { 41335aa4fcfSShri Abhyankar nzi = 0; 41435aa4fcfSShri Abhyankar /* copy current row into linked list */ 41535aa4fcfSShri Abhyankar nnz = ai[r[i]+1] - ai[r[i]]; 416e32f2f54SBarry 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); 41735aa4fcfSShri Abhyankar cols = aj + ai[r[i]]; 41835aa4fcfSShri Abhyankar lnk[i] = -1; /* marker to indicate if diagonal exists */ 41935aa4fcfSShri Abhyankar ierr = PetscIncompleteLLInit(nnz,cols,n,ic,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 42035aa4fcfSShri Abhyankar nzi += nlnk; 42135aa4fcfSShri Abhyankar 42235aa4fcfSShri Abhyankar /* make sure diagonal entry is included */ 42335aa4fcfSShri Abhyankar if (diagonal_fill && lnk[i] == -1) { 42435aa4fcfSShri Abhyankar fm = n; 42535aa4fcfSShri Abhyankar while (lnk[fm] < i) fm = lnk[fm]; 42635aa4fcfSShri Abhyankar lnk[i] = lnk[fm]; /* insert diagonal into linked list */ 42735aa4fcfSShri Abhyankar lnk[fm] = i; 42835aa4fcfSShri Abhyankar lnk_lvl[i] = 0; 42935aa4fcfSShri Abhyankar nzi++; dcount++; 43035aa4fcfSShri Abhyankar } 43135aa4fcfSShri Abhyankar 43235aa4fcfSShri Abhyankar /* add pivot rows into the active row */ 43335aa4fcfSShri Abhyankar nzbd = 0; 43435aa4fcfSShri Abhyankar prow = lnk[n]; 43535aa4fcfSShri Abhyankar while (prow < i) { 43635aa4fcfSShri Abhyankar nnz = bdiag[prow]; 43735aa4fcfSShri Abhyankar cols = bj_ptr[prow] + nnz + 1; 43835aa4fcfSShri Abhyankar cols_lvl = bjlvl_ptr[prow] + nnz + 1; 43935aa4fcfSShri Abhyankar nnz = bi[prow+1] - bi[prow] - nnz - 1; 44026fbe8dcSKarl Rupp 44135aa4fcfSShri Abhyankar ierr = PetscILULLAddSorted(nnz,cols,levels,cols_lvl,prow,nlnk,lnk,lnk_lvl,lnkbt,prow);CHKERRQ(ierr); 44235aa4fcfSShri Abhyankar nzi += nlnk; 44335aa4fcfSShri Abhyankar prow = lnk[prow]; 44435aa4fcfSShri Abhyankar nzbd++; 44535aa4fcfSShri Abhyankar } 44635aa4fcfSShri Abhyankar bdiag[i] = nzbd; 44735aa4fcfSShri Abhyankar bi[i+1] = bi[i] + nzi; 44835aa4fcfSShri Abhyankar 44935aa4fcfSShri Abhyankar /* if free space is not available, make more free space */ 45035aa4fcfSShri Abhyankar if (current_space->local_remaining<nzi) { 451f91af8c7SBarry Smith nnz = PetscIntMultTruncate(2,PetscIntMultTruncate(nzi,(n - i))); /* estimated and max additional space needed */ 45235aa4fcfSShri Abhyankar ierr = PetscFreeSpaceGet(nnz,¤t_space);CHKERRQ(ierr); 45335aa4fcfSShri Abhyankar ierr = PetscFreeSpaceGet(nnz,¤t_space_lvl);CHKERRQ(ierr); 45435aa4fcfSShri Abhyankar reallocs++; 45535aa4fcfSShri Abhyankar } 45635aa4fcfSShri Abhyankar 45735aa4fcfSShri Abhyankar /* copy data into free_space and free_space_lvl, then initialize lnk */ 45835aa4fcfSShri Abhyankar ierr = PetscIncompleteLLClean(n,n,nzi,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr); 45926fbe8dcSKarl Rupp 46035aa4fcfSShri Abhyankar bj_ptr[i] = current_space->array; 46135aa4fcfSShri Abhyankar bjlvl_ptr[i] = current_space_lvl->array; 46235aa4fcfSShri Abhyankar 46335aa4fcfSShri Abhyankar /* make sure the active row i has diagonal entry */ 46465e19b50SBarry 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); 46535aa4fcfSShri Abhyankar 46635aa4fcfSShri Abhyankar current_space->array += nzi; 46735aa4fcfSShri Abhyankar current_space->local_used += nzi; 46835aa4fcfSShri Abhyankar current_space->local_remaining -= nzi; 46926fbe8dcSKarl Rupp 47035aa4fcfSShri Abhyankar current_space_lvl->array += nzi; 47135aa4fcfSShri Abhyankar current_space_lvl->local_used += nzi; 47235aa4fcfSShri Abhyankar current_space_lvl->local_remaining -= nzi; 47335aa4fcfSShri Abhyankar } 47435aa4fcfSShri Abhyankar 47535aa4fcfSShri Abhyankar ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 47635aa4fcfSShri Abhyankar ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 47735aa4fcfSShri Abhyankar 47835aa4fcfSShri Abhyankar /* copy free_space into bj and free free_space; set bi, bj, bdiag in new datastructure; */ 479854ce69bSBarry Smith ierr = PetscMalloc1(bi[n]+1,&bj);CHKERRQ(ierr); 4802ce24eb6SHong Zhang ierr = PetscFreeSpaceContiguous_LU(&free_space,bj,n,bi,bdiag);CHKERRQ(ierr); 48135aa4fcfSShri Abhyankar 48235aa4fcfSShri Abhyankar ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 48335aa4fcfSShri Abhyankar ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr); 484fca92195SBarry Smith ierr = PetscFree2(bj_ptr,bjlvl_ptr);CHKERRQ(ierr); 48535aa4fcfSShri Abhyankar 48635aa4fcfSShri Abhyankar #if defined(PETSC_USE_INFO) 48735aa4fcfSShri Abhyankar { 488aef85c9fSShri Abhyankar PetscReal af = ((PetscReal)(bdiag[0]+1))/((PetscReal)ai[n]); 48957622a8eSBarry Smith ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)f,(double)af);CHKERRQ(ierr); 49057622a8eSBarry Smith ierr = PetscInfo1(A,"Run with -[sub_]pc_factor_fill %g or use \n",(double)af);CHKERRQ(ierr); 49157622a8eSBarry Smith ierr = PetscInfo1(A,"PCFactorSetFill([sub]pc,%g);\n",(double)af);CHKERRQ(ierr); 49235aa4fcfSShri Abhyankar ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr); 49335aa4fcfSShri Abhyankar if (diagonal_fill) { 494955c1f14SBarry Smith ierr = PetscInfo1(A,"Detected and replaced %D missing diagonals\n",dcount);CHKERRQ(ierr); 49535aa4fcfSShri Abhyankar } 49635aa4fcfSShri Abhyankar } 49735aa4fcfSShri Abhyankar #endif 49835aa4fcfSShri Abhyankar 49935aa4fcfSShri Abhyankar /* put together the new matrix */ 500367daffbSBarry Smith ierr = MatSeqBAIJSetPreallocation(fact,bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr); 5013bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)fact,(PetscObject)isicol);CHKERRQ(ierr); 50226fbe8dcSKarl Rupp 50335aa4fcfSShri Abhyankar b = (Mat_SeqBAIJ*)(fact)->data; 50435aa4fcfSShri Abhyankar b->free_a = PETSC_TRUE; 50535aa4fcfSShri Abhyankar b->free_ij = PETSC_TRUE; 50635aa4fcfSShri Abhyankar b->singlemalloc = PETSC_FALSE; 50726fbe8dcSKarl Rupp 508854ce69bSBarry Smith ierr = PetscMalloc1(bs2*(bdiag[0]+1),&b->a);CHKERRQ(ierr); 50926fbe8dcSKarl Rupp 51035aa4fcfSShri Abhyankar b->j = bj; 51135aa4fcfSShri Abhyankar b->i = bi; 51235aa4fcfSShri Abhyankar b->diag = bdiag; 51335aa4fcfSShri Abhyankar b->free_diag = PETSC_TRUE; 51435aa4fcfSShri Abhyankar b->ilen = 0; 51535aa4fcfSShri Abhyankar b->imax = 0; 51635aa4fcfSShri Abhyankar b->row = isrow; 51735aa4fcfSShri Abhyankar b->col = iscol; 51835aa4fcfSShri Abhyankar ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 51935aa4fcfSShri Abhyankar ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 52035aa4fcfSShri Abhyankar b->icol = isicol; 52126fbe8dcSKarl Rupp 522854ce69bSBarry Smith ierr = PetscMalloc1(bs*n+bs,&b->solve_work);CHKERRQ(ierr); 52335aa4fcfSShri Abhyankar /* In b structure: Free imax, ilen, old a, old j. 52435aa4fcfSShri Abhyankar Allocate bdiag, solve_work, new a, new j */ 5253bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)fact,(bdiag[0]+1) * (sizeof(PetscInt)+bs2*sizeof(PetscScalar)));CHKERRQ(ierr); 52635aa4fcfSShri Abhyankar b->maxnz = b->nz = bdiag[0]+1; 52726fbe8dcSKarl Rupp 528ae3d28f0SHong Zhang fact->info.factor_mallocs = reallocs; 529ae3d28f0SHong Zhang fact->info.fill_ratio_given = f; 530ae3d28f0SHong Zhang fact->info.fill_ratio_needed = ((PetscReal)(bdiag[0]+1))/((PetscReal)ai[n]); 53126fbe8dcSKarl Rupp 5324dd39f65SShri Abhyankar ierr = MatSeqBAIJSetNumericFactorization(fact,both_identity);CHKERRQ(ierr); 53335aa4fcfSShri Abhyankar PetscFunctionReturn(0); 53435aa4fcfSShri Abhyankar } 53535aa4fcfSShri Abhyankar 5364e2b4712SSatish Balay /* 5374e2b4712SSatish Balay This code is virtually identical to MatILUFactorSymbolic_SeqAIJ 5384e2b4712SSatish Balay except that the data structure of Mat_SeqAIJ is slightly different. 5394e2b4712SSatish Balay Not a good example of code reuse. 5404e2b4712SSatish Balay */ 5414a2ae208SSatish Balay #undef __FUNCT__ 54206e38f1dSHong Zhang #define __FUNCT__ "MatILUFactorSymbolic_SeqBAIJ_inplace" 54306e38f1dSHong Zhang PetscErrorCode MatILUFactorSymbolic_SeqBAIJ_inplace(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info) 5444e2b4712SSatish Balay { 5454e2b4712SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b; 5464e2b4712SSatish Balay IS isicol; 5476849ba73SBarry Smith PetscErrorCode ierr; 5485d0c19d7SBarry Smith const PetscInt *r,*ic,*ai = a->i,*aj = a->j,*xi; 5495d0c19d7SBarry Smith PetscInt prow,n = a->mbs,*ainew,*ajnew,jmax,*fill,nz,*im,*ajfill,*flev,*xitmp; 550a96a251dSBarry Smith PetscInt *dloc,idx,row,m,fm,nzf,nzi,reallocate = 0,dcount = 0; 551d0f46423SBarry Smith PetscInt incrlev,nnz,i,bs = A->rmap->bs,bs2 = a->bs2,levels,diagonal_fill,dd; 552ace3abfcSBarry Smith PetscBool col_identity,row_identity,both_identity,flg; 553329f5518SBarry Smith PetscReal f; 5544e2b4712SSatish Balay 5554e2b4712SSatish Balay PetscFunctionBegin; 5566bce7ff8SHong Zhang ierr = MatMissingDiagonal_SeqBAIJ(A,&flg,&dd);CHKERRQ(ierr); 557e32f2f54SBarry Smith if (flg) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix A is missing diagonal entry in row %D",dd); 5586bce7ff8SHong Zhang 559435faa5fSBarry Smith f = info->fill; 560690b6cddSBarry Smith levels = (PetscInt)info->levels; 561690b6cddSBarry Smith diagonal_fill = (PetscInt)info->diagonal_fill; 56226fbe8dcSKarl Rupp 5634c49b128SBarry Smith ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 56416a2bf60SHong Zhang 565667159a5SBarry Smith ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 566667159a5SBarry Smith ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr); 567ace3abfcSBarry Smith both_identity = (PetscBool) (row_identity && col_identity); 568309c388cSBarry Smith 56941df41f0SMatthew Knepley if (!levels && both_identity) { /* special case copy the nonzero structure */ 57016a2bf60SHong Zhang ierr = MatDuplicateNoCreate_SeqBAIJ(fact,A,MAT_DO_NOT_COPY_VALUES,PETSC_TRUE);CHKERRQ(ierr); 5718b1456e3SHong Zhang ierr = MatSeqBAIJSetNumericFactorization_inplace(fact,both_identity);CHKERRQ(ierr); 5726bce7ff8SHong Zhang 573d5f3da31SBarry Smith fact->factortype = MAT_FACTOR_ILU; 574ae3d28f0SHong Zhang b = (Mat_SeqBAIJ*)fact->data; 575bb3d539aSBarry Smith b->row = isrow; 576bb3d539aSBarry Smith b->col = iscol; 577bb3d539aSBarry Smith ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 578bb3d539aSBarry Smith ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 579bb3d539aSBarry Smith b->icol = isicol; 580bcd9e38bSBarry Smith b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE; 58126fbe8dcSKarl Rupp 582785e854fSJed Brown ierr = PetscMalloc1((n+1)*bs,&b->solve_work);CHKERRQ(ierr); 5836bce7ff8SHong Zhang PetscFunctionReturn(0); 5846bce7ff8SHong Zhang } 5856bce7ff8SHong Zhang 5866bce7ff8SHong Zhang /* general case perform the symbolic factorization */ 5874e2b4712SSatish Balay ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 5884e2b4712SSatish Balay ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 5894e2b4712SSatish Balay 5904e2b4712SSatish Balay /* get new row pointers */ 591854ce69bSBarry Smith ierr = PetscMalloc1(n+1,&ainew);CHKERRQ(ierr); 5924e2b4712SSatish Balay ainew[0] = 0; 5934e2b4712SSatish Balay /* don't know how many column pointers are needed so estimate */ 594690b6cddSBarry Smith jmax = (PetscInt)(f*ai[n] + 1); 595854ce69bSBarry Smith ierr = PetscMalloc1(jmax,&ajnew);CHKERRQ(ierr); 5964e2b4712SSatish Balay /* ajfill is level of fill for each fill entry */ 597854ce69bSBarry Smith ierr = PetscMalloc1(jmax,&ajfill);CHKERRQ(ierr); 5984e2b4712SSatish Balay /* fill is a linked list of nonzeros in active row */ 599854ce69bSBarry Smith ierr = PetscMalloc1(n+1,&fill);CHKERRQ(ierr); 6004e2b4712SSatish Balay /* im is level for each filled value */ 601854ce69bSBarry Smith ierr = PetscMalloc1(n+1,&im);CHKERRQ(ierr); 6024e2b4712SSatish Balay /* dloc is location of diagonal in factor */ 603854ce69bSBarry Smith ierr = PetscMalloc1(n+1,&dloc);CHKERRQ(ierr); 6044e2b4712SSatish Balay dloc[0] = 0; 6054e2b4712SSatish Balay for (prow=0; prow<n; prow++) { 606435faa5fSBarry Smith 607435faa5fSBarry Smith /* copy prow into linked list */ 6084e2b4712SSatish Balay nzf = nz = ai[r[prow]+1] - ai[r[prow]]; 609e32f2f54SBarry 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); 6104e2b4712SSatish Balay xi = aj + ai[r[prow]]; 6114e2b4712SSatish Balay fill[n] = n; 612435faa5fSBarry Smith fill[prow] = -1; /* marker for diagonal entry */ 6134e2b4712SSatish Balay while (nz--) { 6144e2b4712SSatish Balay fm = n; 6154e2b4712SSatish Balay idx = ic[*xi++]; 6164e2b4712SSatish Balay do { 6174e2b4712SSatish Balay m = fm; 6184e2b4712SSatish Balay fm = fill[m]; 6194e2b4712SSatish Balay } while (fm < idx); 6204e2b4712SSatish Balay fill[m] = idx; 6214e2b4712SSatish Balay fill[idx] = fm; 6224e2b4712SSatish Balay im[idx] = 0; 6234e2b4712SSatish Balay } 624435faa5fSBarry Smith 625435faa5fSBarry Smith /* make sure diagonal entry is included */ 626435faa5fSBarry Smith if (diagonal_fill && fill[prow] == -1) { 627435faa5fSBarry Smith fm = n; 628435faa5fSBarry Smith while (fill[fm] < prow) fm = fill[fm]; 629435faa5fSBarry Smith fill[prow] = fill[fm]; /* insert diagonal into linked list */ 630435faa5fSBarry Smith fill[fm] = prow; 631435faa5fSBarry Smith im[prow] = 0; 632435faa5fSBarry Smith nzf++; 633335d9088SBarry Smith dcount++; 634435faa5fSBarry Smith } 635435faa5fSBarry Smith 6364e2b4712SSatish Balay nzi = 0; 6374e2b4712SSatish Balay row = fill[n]; 6384e2b4712SSatish Balay while (row < prow) { 6394e2b4712SSatish Balay incrlev = im[row] + 1; 6404e2b4712SSatish Balay nz = dloc[row]; 641435faa5fSBarry Smith xi = ajnew + ainew[row] + nz + 1; 6424e2b4712SSatish Balay flev = ajfill + ainew[row] + nz + 1; 6434e2b4712SSatish Balay nnz = ainew[row+1] - ainew[row] - nz - 1; 6444e2b4712SSatish Balay fm = row; 6454e2b4712SSatish Balay while (nnz-- > 0) { 6464e2b4712SSatish Balay idx = *xi++; 6474e2b4712SSatish Balay if (*flev + incrlev > levels) { 6484e2b4712SSatish Balay flev++; 6494e2b4712SSatish Balay continue; 6504e2b4712SSatish Balay } 6514e2b4712SSatish Balay do { 6524e2b4712SSatish Balay m = fm; 6534e2b4712SSatish Balay fm = fill[m]; 6544e2b4712SSatish Balay } while (fm < idx); 6554e2b4712SSatish Balay if (fm != idx) { 6564e2b4712SSatish Balay im[idx] = *flev + incrlev; 6574e2b4712SSatish Balay fill[m] = idx; 6584e2b4712SSatish Balay fill[idx] = fm; 6594e2b4712SSatish Balay fm = idx; 6604e2b4712SSatish Balay nzf++; 66126fbe8dcSKarl Rupp } else if (im[idx] > *flev + incrlev) im[idx] = *flev+incrlev; 6624e2b4712SSatish Balay flev++; 6634e2b4712SSatish Balay } 6644e2b4712SSatish Balay row = fill[row]; 6654e2b4712SSatish Balay nzi++; 6664e2b4712SSatish Balay } 6674e2b4712SSatish Balay /* copy new filled row into permanent storage */ 6684e2b4712SSatish Balay ainew[prow+1] = ainew[prow] + nzf; 6694e2b4712SSatish Balay if (ainew[prow+1] > jmax) { 670ecf371e4SBarry Smith 671ecf371e4SBarry Smith /* estimate how much additional space we will need */ 672ecf371e4SBarry Smith /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */ 673ecf371e4SBarry Smith /* just double the memory each time */ 674690b6cddSBarry Smith PetscInt maxadd = jmax; 675ecf371e4SBarry Smith /* maxadd = (int)(((f*ai[n]+1)*(n-prow+5))/n); */ 6764e2b4712SSatish Balay if (maxadd < nzf) maxadd = (n-prow)*(nzf+1); 6774e2b4712SSatish Balay jmax += maxadd; 678ecf371e4SBarry Smith 679ecf371e4SBarry Smith /* allocate a longer ajnew and ajfill */ 680785e854fSJed Brown ierr = PetscMalloc1(jmax,&xitmp);CHKERRQ(ierr); 6815d0c19d7SBarry Smith ierr = PetscMemcpy(xitmp,ajnew,ainew[prow]*sizeof(PetscInt));CHKERRQ(ierr); 682606d414cSSatish Balay ierr = PetscFree(ajnew);CHKERRQ(ierr); 6835d0c19d7SBarry Smith ajnew = xitmp; 684785e854fSJed Brown ierr = PetscMalloc1(jmax,&xitmp);CHKERRQ(ierr); 6855d0c19d7SBarry Smith ierr = PetscMemcpy(xitmp,ajfill,ainew[prow]*sizeof(PetscInt));CHKERRQ(ierr); 686606d414cSSatish Balay ierr = PetscFree(ajfill);CHKERRQ(ierr); 6875d0c19d7SBarry Smith ajfill = xitmp; 688eb150c5cSKris Buschelman reallocate++; /* count how many reallocations are needed */ 6894e2b4712SSatish Balay } 6905d0c19d7SBarry Smith xitmp = ajnew + ainew[prow]; 6914e2b4712SSatish Balay flev = ajfill + ainew[prow]; 6924e2b4712SSatish Balay dloc[prow] = nzi; 6934e2b4712SSatish Balay fm = fill[n]; 6944e2b4712SSatish Balay while (nzf--) { 6955d0c19d7SBarry Smith *xitmp++ = fm; 6964e2b4712SSatish Balay *flev++ = im[fm]; 6974e2b4712SSatish Balay fm = fill[fm]; 6984e2b4712SSatish Balay } 699435faa5fSBarry Smith /* make sure row has diagonal entry */ 700f23aa3ddSBarry 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\ 7012401956bSBarry Smith try running with -pc_factor_nonzeros_along_diagonal or -pc_factor_diagonal_fill",prow); 702435faa5fSBarry Smith } 703606d414cSSatish Balay ierr = PetscFree(ajfill);CHKERRQ(ierr); 7044e2b4712SSatish Balay ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 7054e2b4712SSatish Balay ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 706606d414cSSatish Balay ierr = PetscFree(fill);CHKERRQ(ierr); 707606d414cSSatish Balay ierr = PetscFree(im);CHKERRQ(ierr); 7084e2b4712SSatish Balay 7096cf91177SBarry Smith #if defined(PETSC_USE_INFO) 7104e2b4712SSatish Balay { 711329f5518SBarry Smith PetscReal af = ((PetscReal)ainew[n])/((PetscReal)ai[n]); 71257622a8eSBarry Smith ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocate,(double)f,(double)af);CHKERRQ(ierr); 71357622a8eSBarry Smith ierr = PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);CHKERRQ(ierr); 71457622a8eSBarry Smith ierr = PetscInfo1(A,"PCFactorSetFill(pc,%g);\n",(double)af);CHKERRQ(ierr); 715ae15b995SBarry Smith ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr); 716335d9088SBarry Smith if (diagonal_fill) { 717ae15b995SBarry Smith ierr = PetscInfo1(A,"Detected and replaced %D missing diagonals\n",dcount);CHKERRQ(ierr); 718335d9088SBarry Smith } 7194e2b4712SSatish Balay } 72063ba0a88SBarry Smith #endif 7214e2b4712SSatish Balay 7224e2b4712SSatish Balay /* put together the new matrix */ 723367daffbSBarry Smith ierr = MatSeqBAIJSetPreallocation(fact,bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr); 7243bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)fact,(PetscObject)isicol);CHKERRQ(ierr); 725ae3d28f0SHong Zhang b = (Mat_SeqBAIJ*)fact->data; 72626fbe8dcSKarl Rupp 727e6b907acSBarry Smith b->free_a = PETSC_TRUE; 728e6b907acSBarry Smith b->free_ij = PETSC_TRUE; 7297c922b88SBarry Smith b->singlemalloc = PETSC_FALSE; 73026fbe8dcSKarl Rupp 731785e854fSJed Brown ierr = PetscMalloc1(bs2*ainew[n],&b->a);CHKERRQ(ierr); 73226fbe8dcSKarl Rupp 7334e2b4712SSatish Balay b->j = ajnew; 7344e2b4712SSatish Balay b->i = ainew; 7354e2b4712SSatish Balay for (i=0; i<n; i++) dloc[i] += ainew[i]; 7364e2b4712SSatish Balay b->diag = dloc; 7377f53bb6cSHong Zhang b->free_diag = PETSC_TRUE; 7384e2b4712SSatish Balay b->ilen = 0; 7394e2b4712SSatish Balay b->imax = 0; 7404e2b4712SSatish Balay b->row = isrow; 7414e2b4712SSatish Balay b->col = iscol; 742bcd9e38bSBarry Smith b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE; 74326fbe8dcSKarl Rupp 744c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 745c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 746e51c0b9cSSatish Balay b->icol = isicol; 747854ce69bSBarry Smith ierr = PetscMalloc1(bs*n+bs,&b->solve_work);CHKERRQ(ierr); 7484e2b4712SSatish Balay /* In b structure: Free imax, ilen, old a, old j. 7494e2b4712SSatish Balay Allocate dloc, solve_work, new a, new j */ 7503bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)fact,(ainew[n]-n)*(sizeof(PetscInt))+bs2*ainew[n]*sizeof(PetscScalar));CHKERRQ(ierr); 7514e2b4712SSatish Balay b->maxnz = b->nz = ainew[n]; 7524e2b4712SSatish Balay 753ae3d28f0SHong Zhang fact->info.factor_mallocs = reallocate; 754ae3d28f0SHong Zhang fact->info.fill_ratio_given = f; 755ae3d28f0SHong Zhang fact->info.fill_ratio_needed = ((PetscReal)ainew[n])/((PetscReal)ai[prow]); 7566bce7ff8SHong Zhang 7578b1456e3SHong Zhang ierr = MatSeqBAIJSetNumericFactorization_inplace(fact,both_identity);CHKERRQ(ierr); 7588661488fSKris Buschelman PetscFunctionReturn(0); 7598661488fSKris Buschelman } 7608661488fSKris Buschelman 761732ee342SKris Buschelman #undef __FUNCT__ 7627e7071cdSKris Buschelman #define __FUNCT__ "MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE" 763dfbe8321SBarry Smith PetscErrorCode MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE(Mat A) 7647e7071cdSKris Buschelman { 76512272027SHong Zhang /* Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; */ 76612272027SHong Zhang /* int i,*AJ=a->j,nz=a->nz; */ 7675fd66863SKarl Rupp 7685a9542e3SKris Buschelman PetscFunctionBegin; 7697cf1b8d3SKris Buschelman /* Undo Column scaling */ 7707cf1b8d3SKris Buschelman /* while (nz--) { */ 7717cf1b8d3SKris Buschelman /* AJ[i] = AJ[i]/4; */ 7727cf1b8d3SKris Buschelman /* } */ 773c115a38dSKris Buschelman /* This should really invoke a push/pop logic, but we don't have that yet. */ 7740298fd71SBarry Smith A->ops->setunfactored = NULL; 7757cf1b8d3SKris Buschelman PetscFunctionReturn(0); 7767cf1b8d3SKris Buschelman } 7777cf1b8d3SKris Buschelman 7787cf1b8d3SKris Buschelman #undef __FUNCT__ 7797cf1b8d3SKris Buschelman #define __FUNCT__ "MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE_usj" 780dfbe8321SBarry Smith PetscErrorCode MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE_usj(Mat A) 7817cf1b8d3SKris Buschelman { 7827cf1b8d3SKris Buschelman Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 783b24ad042SBarry Smith PetscInt *AJ=a->j,nz=a->nz; 7842aa5897fSKris Buschelman unsigned short *aj=(unsigned short*)AJ; 7855fd66863SKarl Rupp 7865a9542e3SKris Buschelman PetscFunctionBegin; 7870b9da03eSKris Buschelman /* Is this really necessary? */ 78820235379SKris Buschelman while (nz--) { 7890b9da03eSKris Buschelman AJ[nz] = (int)((unsigned int)aj[nz]); /* First extend, then convert to signed. */ 7907e7071cdSKris Buschelman } 7910298fd71SBarry Smith A->ops->setunfactored = NULL; 7927e7071cdSKris Buschelman PetscFunctionReturn(0); 7937e7071cdSKris Buschelman } 7947e7071cdSKris Buschelman 795732ee342SKris Buschelman 796