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; 21766f9fbaSBarry Smith PetscInt i,j,k,ipvt[15]; 22766f9fbaSBarry Smith const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*ajtmp,*bjtmp,*bdiag=b->diag,*pj; 23766f9fbaSBarry Smith PetscInt nz,nzL,row; 24766f9fbaSBarry Smith MatScalar *rtmp,*pc,*mwork,*pv,*vv,work[225]; 25766f9fbaSBarry Smith const MatScalar *v,*aa=a->a; 262b0b2ea7SShri Abhyankar PetscInt bs2 = a->bs2,bs=A->rmap->bs,flg; 270fa040f9SShri Abhyankar PetscInt sol_ver; 28a455e926SHong Zhang PetscBool allowzeropivot,zeropivotdetected; 292b0b2ea7SShri Abhyankar 302b0b2ea7SShri Abhyankar PetscFunctionBegin; 310164db54SHong Zhang allowzeropivot = PetscNot(A->erroriffailure); 325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,((PetscObject)A)->prefix,"-sol_ver",&sol_ver,NULL)); 330fa040f9SShri Abhyankar 342b0b2ea7SShri Abhyankar /* generate work space needed by the factorization */ 355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(bs2*n,&rtmp,bs2,&mwork)); 365f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArrayzero(rtmp,bs2*n)); 372b0b2ea7SShri Abhyankar 382b0b2ea7SShri Abhyankar for (i=0; i<n; i++) { 392b0b2ea7SShri Abhyankar /* zero rtmp */ 402b0b2ea7SShri Abhyankar /* L part */ 412b0b2ea7SShri Abhyankar nz = bi[i+1] - bi[i]; 422b0b2ea7SShri Abhyankar bjtmp = bj + bi[i]; 432b0b2ea7SShri Abhyankar for (j=0; j<nz; j++) { 445f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArrayzero(rtmp+bs2*bjtmp[j],bs2)); 452b0b2ea7SShri Abhyankar } 462b0b2ea7SShri Abhyankar 472b0b2ea7SShri Abhyankar /* U part */ 482b0b2ea7SShri Abhyankar nz = bdiag[i] - bdiag[i+1]; 492b0b2ea7SShri Abhyankar bjtmp = bj + bdiag[i+1]+1; 502b0b2ea7SShri Abhyankar for (j=0; j<nz; j++) { 515f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArrayzero(rtmp+bs2*bjtmp[j],bs2)); 522b0b2ea7SShri Abhyankar } 532b0b2ea7SShri Abhyankar 542b0b2ea7SShri Abhyankar /* load in initial (unfactored row) */ 5529a97285SShri Abhyankar nz = ai[i+1] - ai[i]; 5629a97285SShri Abhyankar ajtmp = aj + ai[i]; 5729a97285SShri Abhyankar v = aa + bs2*ai[i]; 582b0b2ea7SShri Abhyankar for (j=0; j<nz; j++) { 595f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(rtmp+bs2*ajtmp[j],v+bs2*j,bs2)); 602b0b2ea7SShri Abhyankar } 612b0b2ea7SShri Abhyankar 622b0b2ea7SShri Abhyankar /* elimination */ 632b0b2ea7SShri Abhyankar bjtmp = bj + bi[i]; 642b0b2ea7SShri Abhyankar nzL = bi[i+1] - bi[i]; 652b0b2ea7SShri Abhyankar for (k=0; k < nzL; k++) { 662b0b2ea7SShri Abhyankar row = bjtmp[k]; 672b0b2ea7SShri Abhyankar pc = rtmp + bs2*row; 68c35f09e5SBarry Smith for (flg=0,j=0; j<bs2; j++) { 69c35f09e5SBarry Smith if (pc[j]!=0.0) { 70c35f09e5SBarry Smith flg = 1; 71c35f09e5SBarry Smith break; 72c35f09e5SBarry Smith } 73c35f09e5SBarry Smith } 742b0b2ea7SShri Abhyankar if (flg) { 752b0b2ea7SShri Abhyankar pv = b->a + bs2*bdiag[row]; 7696b95a6bSBarry Smith PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); 775f80ce2aSJacob Faibussowitsch /*CHKERRQ(PetscKernel_A_gets_A_times_B_15(pc,pv,mwork));*/ 78a5b23f4aSJose E. Roman pj = b->j + bdiag[row+1]+1; /* beginning of U(row,:) */ 792b0b2ea7SShri Abhyankar pv = b->a + bs2*(bdiag[row+1]+1); 802b0b2ea7SShri Abhyankar nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries inU(row,:), excluding diag */ 812b0b2ea7SShri Abhyankar for (j=0; j<nz; j++) { 82766f9fbaSBarry Smith vv = rtmp + bs2*pj[j]; 8396b95a6bSBarry Smith PetscKernel_A_gets_A_minus_B_times_C(bs,vv,pc,pv); 845f80ce2aSJacob Faibussowitsch /* CHKERRQ(PetscKernel_A_gets_A_minus_B_times_C_15(vv,pc,pv)); */ 852b0b2ea7SShri Abhyankar pv += bs2; 862b0b2ea7SShri Abhyankar } 875f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogFlops(2.0*bs2*bs*(nz+1)-bs2)); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */ 882b0b2ea7SShri Abhyankar } 892b0b2ea7SShri Abhyankar } 902b0b2ea7SShri Abhyankar 912b0b2ea7SShri Abhyankar /* finished row so stick it into b->a */ 922b0b2ea7SShri Abhyankar /* L part */ 932b0b2ea7SShri Abhyankar pv = b->a + bs2*bi[i]; 942b0b2ea7SShri Abhyankar pj = b->j + bi[i]; 952b0b2ea7SShri Abhyankar nz = bi[i+1] - bi[i]; 962b0b2ea7SShri Abhyankar for (j=0; j<nz; j++) { 975f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(pv+bs2*j,rtmp+bs2*pj[j],bs2)); 982b0b2ea7SShri Abhyankar } 992b0b2ea7SShri Abhyankar 100a5b23f4aSJose E. Roman /* Mark diagonal and invert diagonal for simpler triangular solves */ 1012b0b2ea7SShri Abhyankar pv = b->a + bs2*bdiag[i]; 1022b0b2ea7SShri Abhyankar pj = b->j + bdiag[i]; 1035f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(pv,rtmp+bs2*pj[0],bs2)); 1045f80ce2aSJacob Faibussowitsch CHKERRQ(PetscKernel_A_gets_inverse_A_15(pv,ipvt,work,info->shiftamount,allowzeropivot,&zeropivotdetected)); 1057b6c816cSBarry Smith if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 1062b0b2ea7SShri Abhyankar 1072b0b2ea7SShri Abhyankar /* U part */ 1082b0b2ea7SShri Abhyankar pv = b->a + bs2*(bdiag[i+1]+1); 1092b0b2ea7SShri Abhyankar pj = b->j + bdiag[i+1]+1; 1102b0b2ea7SShri Abhyankar nz = bdiag[i] - bdiag[i+1] - 1; 1112b0b2ea7SShri Abhyankar for (j=0; j<nz; j++) { 1125f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(pv+bs2*j,rtmp+bs2*pj[j],bs2)); 1132b0b2ea7SShri Abhyankar } 1142b0b2ea7SShri Abhyankar } 1152b0b2ea7SShri Abhyankar 1165f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(rtmp,mwork)); 11726fbe8dcSKarl Rupp 118832cc040SShri Abhyankar C->ops->solve = MatSolve_SeqBAIJ_15_NaturalOrdering_ver1; 119766f9fbaSBarry Smith C->ops->solvetranspose = MatSolve_SeqBAIJ_N_NaturalOrdering; 1202b0b2ea7SShri Abhyankar C->assembled = PETSC_TRUE; 12126fbe8dcSKarl Rupp 1225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogFlops(1.333333333333*bs*bs2*b->mbs)); /* from inverting diagonal blocks */ 1232b0b2ea7SShri Abhyankar PetscFunctionReturn(0); 1242b0b2ea7SShri Abhyankar } 1252b0b2ea7SShri Abhyankar 1264dd39f65SShri Abhyankar PetscErrorCode MatLUFactorNumeric_SeqBAIJ_N(Mat B,Mat A,const MatFactorInfo *info) 1276bce7ff8SHong Zhang { 1286bce7ff8SHong Zhang Mat C =B; 1296bce7ff8SHong Zhang Mat_SeqBAIJ *a =(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data; 1306bce7ff8SHong Zhang IS isrow = b->row,isicol = b->icol; 1315a586d82SBarry Smith const PetscInt *r,*ic; 1326bce7ff8SHong Zhang PetscInt i,j,k,n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 1336bce7ff8SHong Zhang PetscInt *ajtmp,*bjtmp,nz,nzL,row,*bdiag=b->diag,*pj; 134b588c5a2SHong Zhang MatScalar *rtmp,*pc,*mwork,*v,*pv,*aa=a->a; 135914a18a2SHong Zhang PetscInt bs=A->rmap->bs,bs2 = a->bs2,*v_pivots,flg; 136914a18a2SHong Zhang MatScalar *v_work; 137ace3abfcSBarry Smith PetscBool col_identity,row_identity,both_identity; 1385f8bbccaSHong Zhang PetscBool allowzeropivot,zeropivotdetected; 1396bce7ff8SHong Zhang 1406bce7ff8SHong Zhang PetscFunctionBegin; 1415f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(isrow,&r)); 1425f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(isicol,&ic)); 1435f8bbccaSHong Zhang allowzeropivot = PetscNot(A->erroriffailure); 144ae3d28f0SHong Zhang 1455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscCalloc1(bs2*n,&rtmp)); 1466bce7ff8SHong Zhang 147914a18a2SHong Zhang /* generate work space needed by dense LU factorization */ 1485f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc3(bs,&v_work,bs2,&mwork,bs,&v_pivots)); 149914a18a2SHong Zhang 1506bce7ff8SHong Zhang for (i=0; i<n; i++) { 1516bce7ff8SHong Zhang /* zero rtmp */ 1526bce7ff8SHong Zhang /* L part */ 1536bce7ff8SHong Zhang nz = bi[i+1] - bi[i]; 1546bce7ff8SHong Zhang bjtmp = bj + bi[i]; 155914a18a2SHong Zhang for (j=0; j<nz; j++) { 1565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArrayzero(rtmp+bs2*bjtmp[j],bs2)); 157914a18a2SHong Zhang } 1586bce7ff8SHong Zhang 1596bce7ff8SHong Zhang /* U part */ 1601a83e813SShri Abhyankar nz = bdiag[i] - bdiag[i+1]; 1611a83e813SShri Abhyankar bjtmp = bj + bdiag[i+1]+1; 1621a83e813SShri Abhyankar for (j=0; j<nz; j++) { 1635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArrayzero(rtmp+bs2*bjtmp[j],bs2)); 1641a83e813SShri Abhyankar } 1651a83e813SShri Abhyankar 1661a83e813SShri Abhyankar /* load in initial (unfactored row) */ 1671a83e813SShri Abhyankar nz = ai[r[i]+1] - ai[r[i]]; 1681a83e813SShri Abhyankar ajtmp = aj + ai[r[i]]; 1691a83e813SShri Abhyankar v = aa + bs2*ai[r[i]]; 1701a83e813SShri Abhyankar for (j=0; j<nz; j++) { 1715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(rtmp+bs2*ic[ajtmp[j]],v+bs2*j,bs2)); 1721a83e813SShri Abhyankar } 1731a83e813SShri Abhyankar 1741a83e813SShri Abhyankar /* elimination */ 1751a83e813SShri Abhyankar bjtmp = bj + bi[i]; 1761a83e813SShri Abhyankar nzL = bi[i+1] - bi[i]; 1771a83e813SShri Abhyankar for (k=0; k < nzL; k++) { 1781a83e813SShri Abhyankar row = bjtmp[k]; 1791a83e813SShri Abhyankar pc = rtmp + bs2*row; 180c35f09e5SBarry Smith for (flg=0,j=0; j<bs2; j++) { 181c35f09e5SBarry Smith if (pc[j]!=0.0) { 182c35f09e5SBarry Smith flg = 1; 183c35f09e5SBarry Smith break; 184c35f09e5SBarry Smith } 185c35f09e5SBarry Smith } 1861a83e813SShri Abhyankar if (flg) { 1871a83e813SShri Abhyankar pv = b->a + bs2*bdiag[row]; 18896b95a6bSBarry Smith PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); /* *pc = *pc * (*pv); */ 189a5b23f4aSJose E. Roman pj = b->j + bdiag[row+1]+1; /* beginning of U(row,:) */ 1901a83e813SShri Abhyankar pv = b->a + bs2*(bdiag[row+1]+1); 1911a83e813SShri Abhyankar nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries inU(row,:), excluding diag */ 1921a83e813SShri Abhyankar for (j=0; j<nz; j++) { 19396b95a6bSBarry Smith PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); 1941a83e813SShri Abhyankar } 1955f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogFlops(2.0*bs2*bs*(nz+1)-bs2)); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */ 1961a83e813SShri Abhyankar } 1971a83e813SShri Abhyankar } 1981a83e813SShri Abhyankar 1991a83e813SShri Abhyankar /* finished row so stick it into b->a */ 2001a83e813SShri Abhyankar /* L part */ 2011a83e813SShri Abhyankar pv = b->a + bs2*bi[i]; 2021a83e813SShri Abhyankar pj = b->j + bi[i]; 2031a83e813SShri Abhyankar nz = bi[i+1] - bi[i]; 2041a83e813SShri Abhyankar for (j=0; j<nz; j++) { 2055f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(pv+bs2*j,rtmp+bs2*pj[j],bs2)); 2061a83e813SShri Abhyankar } 2071a83e813SShri Abhyankar 208a5b23f4aSJose E. Roman /* Mark diagonal and invert diagonal for simpler triangular solves */ 2091a83e813SShri Abhyankar pv = b->a + bs2*bdiag[i]; 2101a83e813SShri Abhyankar pj = b->j + bdiag[i]; 2115f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(pv,rtmp+bs2*pj[0],bs2)); 2125f8bbccaSHong Zhang 2135f80ce2aSJacob Faibussowitsch CHKERRQ(PetscKernel_A_gets_inverse_A(bs,pv,v_pivots,v_work,allowzeropivot,&zeropivotdetected)); 2147b6c816cSBarry Smith if (zeropivotdetected) B->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 2151a83e813SShri Abhyankar 2161a83e813SShri Abhyankar /* U part */ 2171a83e813SShri Abhyankar pv = b->a + bs2*(bdiag[i+1]+1); 2181a83e813SShri Abhyankar pj = b->j + bdiag[i+1]+1; 2191a83e813SShri Abhyankar nz = bdiag[i] - bdiag[i+1] - 1; 2201a83e813SShri Abhyankar for (j=0; j<nz; j++) { 2215f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(pv+bs2*j,rtmp+bs2*pj[j],bs2)); 2221a83e813SShri Abhyankar } 2231a83e813SShri Abhyankar } 2241a83e813SShri Abhyankar 2255f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(rtmp)); 2265f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree3(v_work,mwork,v_pivots)); 2275f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(isicol,&ic)); 2285f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(isrow,&r)); 2291a83e813SShri Abhyankar 2305f80ce2aSJacob Faibussowitsch CHKERRQ(ISIdentity(isrow,&row_identity)); 2315f80ce2aSJacob Faibussowitsch CHKERRQ(ISIdentity(isicol,&col_identity)); 23226fbe8dcSKarl Rupp 233ace3abfcSBarry Smith both_identity = (PetscBool) (row_identity && col_identity); 234ae3d28f0SHong Zhang if (both_identity) { 235ba7f0461SHong Zhang switch (bs) { 23696e086a2SDaniel Kokron case 9: 2375f70456aSHong Zhang #if defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX2__) && defined(__FMA__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES) 23896e086a2SDaniel Kokron C->ops->solve = MatSolve_SeqBAIJ_9_NaturalOrdering; 239aee5c371SBarry Smith #else 240aee5c371SBarry Smith C->ops->solve = MatSolve_SeqBAIJ_N_NaturalOrdering; 241aee5c371SBarry Smith #endif 24296e086a2SDaniel Kokron break; 243ba7f0461SHong Zhang case 11: 244ba7f0461SHong Zhang C->ops->solve = MatSolve_SeqBAIJ_11_NaturalOrdering; 245ba7f0461SHong Zhang break; 246ba7f0461SHong Zhang case 12: 247ba7f0461SHong Zhang C->ops->solve = MatSolve_SeqBAIJ_12_NaturalOrdering; 248ba7f0461SHong Zhang break; 249ba7f0461SHong Zhang case 13: 250ba7f0461SHong Zhang C->ops->solve = MatSolve_SeqBAIJ_13_NaturalOrdering; 251ba7f0461SHong Zhang break; 252ba7f0461SHong Zhang case 14: 253ba7f0461SHong Zhang C->ops->solve = MatSolve_SeqBAIJ_14_NaturalOrdering; 254ba7f0461SHong Zhang break; 255ba7f0461SHong Zhang default: 2564dd39f65SShri Abhyankar C->ops->solve = MatSolve_SeqBAIJ_N_NaturalOrdering; 257ba7f0461SHong Zhang break; 258ba7f0461SHong Zhang } 259ae3d28f0SHong Zhang } else { 2604dd39f65SShri Abhyankar C->ops->solve = MatSolve_SeqBAIJ_N; 261ae3d28f0SHong Zhang } 2624dd39f65SShri Abhyankar C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_N; 263ae3d28f0SHong Zhang 2641a83e813SShri Abhyankar C->assembled = PETSC_TRUE; 26526fbe8dcSKarl Rupp 2665f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogFlops(1.333333333333*bs*bs2*b->mbs)); /* from inverting diagonal blocks */ 2671a83e813SShri Abhyankar PetscFunctionReturn(0); 2681a83e813SShri Abhyankar } 2691a83e813SShri Abhyankar 2706bce7ff8SHong Zhang /* 2716bce7ff8SHong Zhang ilu(0) with natural ordering under new data structure. 2724dd39f65SShri Abhyankar See MatILUFactorSymbolic_SeqAIJ_ilu0() for detailed description 2734dd39f65SShri Abhyankar because this code is almost identical to MatILUFactorSymbolic_SeqAIJ_ilu0_inplace(). 2746bce7ff8SHong Zhang */ 275c0c7eb62SShri Abhyankar 2764dd39f65SShri Abhyankar PetscErrorCode MatILUFactorSymbolic_SeqBAIJ_ilu0(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info) 2776bce7ff8SHong Zhang { 2786bce7ff8SHong Zhang 2796bce7ff8SHong Zhang Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b; 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; 2845f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicateNoCreate_SeqBAIJ(fact,A,MAT_DO_NOT_COPY_VALUES,PETSC_FALSE)); 28535aa4fcfSShri Abhyankar b = (Mat_SeqBAIJ*)(fact)->data; 28635aa4fcfSShri Abhyankar 28735aa4fcfSShri Abhyankar /* allocate matrix arrays for new data structure */ 2885f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc3(bs2*ai[n]+1,&b->a,ai[n]+1,&b->j,n+1,&b->i)); 2895f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogObjectMemory((PetscObject)fact,ai[n]*(bs2*sizeof(PetscScalar)+sizeof(PetscInt))+(n+1)*sizeof(PetscInt))); 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) { 2975f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(n+1,&b->diag)); 2985f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogObjectMemory((PetscObject)fact,(n+1)*sizeof(PetscInt))); 29935aa4fcfSShri Abhyankar } 30035aa4fcfSShri Abhyankar bdiag = b->diag; 30135aa4fcfSShri Abhyankar 30235aa4fcfSShri Abhyankar if (n > 0) { 3035f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArrayzero(b->a,bs2*ai[n])); 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 const PetscInt *r,*ic; 3437fa3a6a0SHong Zhang PetscInt n=a->mbs,*ai=a->i,*aj=a->j,d; 34416a2bf60SHong Zhang PetscInt *bi,*cols,nnz,*cols_lvl; 34516a2bf60SHong Zhang PetscInt *bdiag,prow,fm,nzbd,reallocs=0,dcount=0; 34616a2bf60SHong Zhang PetscInt i,levels,diagonal_fill; 347ace3abfcSBarry Smith PetscBool col_identity,row_identity,both_identity; 34816a2bf60SHong Zhang PetscReal f; 3490298fd71SBarry Smith PetscInt nlnk,*lnk,*lnk_lvl=NULL; 35016a2bf60SHong Zhang PetscBT lnkbt; 35116a2bf60SHong Zhang PetscInt nzi,*bj,**bj_ptr,**bjlvl_ptr; 3520298fd71SBarry Smith PetscFreeSpaceList free_space =NULL,current_space=NULL; 3530298fd71SBarry Smith PetscFreeSpaceList free_space_lvl=NULL,current_space_lvl=NULL; 354ace3abfcSBarry Smith PetscBool missing; 3557fa3a6a0SHong Zhang PetscInt bs=A->rmap->bs,bs2=a->bs2; 35616a2bf60SHong Zhang 35716a2bf60SHong Zhang PetscFunctionBegin; 3582c71b3e2SJacob Faibussowitsch PetscCheckFalse(A->rmap->n != A->cmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Must be square matrix, rows %" PetscInt_FMT " columns %" PetscInt_FMT,A->rmap->n,A->cmap->n); 3596ba06ab7SHong Zhang if (bs>1) { /* check shifttype */ 3602c71b3e2SJacob Faibussowitsch PetscCheckFalse(info->shifttype == MAT_SHIFT_NONZERO || info->shifttype == MAT_SHIFT_POSITIVE_DEFINITE,PETSC_COMM_SELF,PETSC_ERR_SUP,"Only MAT_SHIFT_NONE and MAT_SHIFT_INBLOCKS are supported for BAIJ matrix"); 3616ba06ab7SHong Zhang } 3626ba06ab7SHong Zhang 3635f80ce2aSJacob Faibussowitsch CHKERRQ(MatMissingDiagonal(A,&missing,&d)); 364*28b400f6SJacob Faibussowitsch PetscCheck(!missing,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %" PetscInt_FMT,d); 36516a2bf60SHong Zhang 36616a2bf60SHong Zhang f = info->fill; 36716a2bf60SHong Zhang levels = (PetscInt)info->levels; 36816a2bf60SHong Zhang diagonal_fill = (PetscInt)info->diagonal_fill; 36926fbe8dcSKarl Rupp 3705f80ce2aSJacob Faibussowitsch CHKERRQ(ISInvertPermutation(iscol,PETSC_DECIDE,&isicol)); 37116a2bf60SHong Zhang 3725f80ce2aSJacob Faibussowitsch CHKERRQ(ISIdentity(isrow,&row_identity)); 3735f80ce2aSJacob Faibussowitsch CHKERRQ(ISIdentity(iscol,&col_identity)); 37426fbe8dcSKarl Rupp 375ace3abfcSBarry Smith both_identity = (PetscBool) (row_identity && col_identity); 37616a2bf60SHong Zhang 3777fa3a6a0SHong Zhang if (!levels && both_identity) { 37816a2bf60SHong Zhang /* special case: ilu(0) with natural ordering */ 3795f80ce2aSJacob Faibussowitsch CHKERRQ(MatILUFactorSymbolic_SeqBAIJ_ilu0(fact,A,isrow,iscol,info)); 3805f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqBAIJSetNumericFactorization(fact,both_identity)); 38135aa4fcfSShri Abhyankar 382d5f3da31SBarry Smith fact->factortype = MAT_FACTOR_ILU; 38335aa4fcfSShri Abhyankar (fact)->info.factor_mallocs = 0; 38435aa4fcfSShri Abhyankar (fact)->info.fill_ratio_given = info->fill; 38535aa4fcfSShri Abhyankar (fact)->info.fill_ratio_needed = 1.0; 38626fbe8dcSKarl Rupp 38735aa4fcfSShri Abhyankar b = (Mat_SeqBAIJ*)(fact)->data; 38835aa4fcfSShri Abhyankar b->row = isrow; 38935aa4fcfSShri Abhyankar b->col = iscol; 39035aa4fcfSShri Abhyankar b->icol = isicol; 3915f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)isrow)); 3925f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)iscol)); 39335aa4fcfSShri Abhyankar b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE; 39426fbe8dcSKarl Rupp 3955f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1((n+1)*bs,&b->solve_work)); 39635aa4fcfSShri Abhyankar PetscFunctionReturn(0); 39735aa4fcfSShri Abhyankar } 39835aa4fcfSShri Abhyankar 3995f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(isrow,&r)); 4005f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(isicol,&ic)); 40135aa4fcfSShri Abhyankar 40235aa4fcfSShri Abhyankar /* get new row pointers */ 4035f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(n+1,&bi)); 40435aa4fcfSShri Abhyankar bi[0] = 0; 40535aa4fcfSShri Abhyankar /* bdiag is location of diagonal in factor */ 4065f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(n+1,&bdiag)); 40735aa4fcfSShri Abhyankar bdiag[0] = 0; 40835aa4fcfSShri Abhyankar 4095f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(n,&bj_ptr,n,&bjlvl_ptr)); 41035aa4fcfSShri Abhyankar 41135aa4fcfSShri Abhyankar /* create a linked list for storing column indices of the active row */ 41235aa4fcfSShri Abhyankar nlnk = n + 1; 4135f80ce2aSJacob Faibussowitsch CHKERRQ(PetscIncompleteLLCreate(n,n,nlnk,lnk,lnk_lvl,lnkbt)); 41435aa4fcfSShri Abhyankar 41535aa4fcfSShri Abhyankar /* initial FreeSpace size is f*(ai[n]+1) */ 4165f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFreeSpaceGet(PetscRealIntMultTruncate(f,ai[n]+1),&free_space)); 41735aa4fcfSShri Abhyankar current_space = free_space; 4185f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFreeSpaceGet(PetscRealIntMultTruncate(f,ai[n]+1),&free_space_lvl)); 41935aa4fcfSShri Abhyankar current_space_lvl = free_space_lvl; 42035aa4fcfSShri Abhyankar 42135aa4fcfSShri Abhyankar for (i=0; i<n; i++) { 42235aa4fcfSShri Abhyankar nzi = 0; 42335aa4fcfSShri Abhyankar /* copy current row into linked list */ 42435aa4fcfSShri Abhyankar nnz = ai[r[i]+1] - ai[r[i]]; 42594bad497SJacob Faibussowitsch PetscCheck(nnz,PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %" PetscInt_FMT " in permuted ordering %" PetscInt_FMT,r[i],i); 42635aa4fcfSShri Abhyankar cols = aj + ai[r[i]]; 42735aa4fcfSShri Abhyankar lnk[i] = -1; /* marker to indicate if diagonal exists */ 42894bad497SJacob Faibussowitsch CHKERRQ(PetscIncompleteLLInit(nnz,cols,n,ic,&nlnk,lnk,lnk_lvl,lnkbt)); 42935aa4fcfSShri Abhyankar nzi += nlnk; 43035aa4fcfSShri Abhyankar 43135aa4fcfSShri Abhyankar /* make sure diagonal entry is included */ 43235aa4fcfSShri Abhyankar if (diagonal_fill && lnk[i] == -1) { 43335aa4fcfSShri Abhyankar fm = n; 43435aa4fcfSShri Abhyankar while (lnk[fm] < i) fm = lnk[fm]; 43535aa4fcfSShri Abhyankar lnk[i] = lnk[fm]; /* insert diagonal into linked list */ 43635aa4fcfSShri Abhyankar lnk[fm] = i; 43735aa4fcfSShri Abhyankar lnk_lvl[i] = 0; 43835aa4fcfSShri Abhyankar nzi++; dcount++; 43935aa4fcfSShri Abhyankar } 44035aa4fcfSShri Abhyankar 44135aa4fcfSShri Abhyankar /* add pivot rows into the active row */ 44235aa4fcfSShri Abhyankar nzbd = 0; 44335aa4fcfSShri Abhyankar prow = lnk[n]; 44435aa4fcfSShri Abhyankar while (prow < i) { 44535aa4fcfSShri Abhyankar nnz = bdiag[prow]; 44635aa4fcfSShri Abhyankar cols = bj_ptr[prow] + nnz + 1; 44735aa4fcfSShri Abhyankar cols_lvl = bjlvl_ptr[prow] + nnz + 1; 44835aa4fcfSShri Abhyankar nnz = bi[prow+1] - bi[prow] - nnz - 1; 44926fbe8dcSKarl Rupp 45094bad497SJacob Faibussowitsch CHKERRQ(PetscILULLAddSorted(nnz,cols,levels,cols_lvl,prow,&nlnk,lnk,lnk_lvl,lnkbt,prow)); 45135aa4fcfSShri Abhyankar nzi += nlnk; 45235aa4fcfSShri Abhyankar prow = lnk[prow]; 45335aa4fcfSShri Abhyankar nzbd++; 45435aa4fcfSShri Abhyankar } 45535aa4fcfSShri Abhyankar bdiag[i] = nzbd; 45635aa4fcfSShri Abhyankar bi[i+1] = bi[i] + nzi; 45735aa4fcfSShri Abhyankar 45835aa4fcfSShri Abhyankar /* if free space is not available, make more free space */ 45935aa4fcfSShri Abhyankar if (current_space->local_remaining<nzi) { 460f91af8c7SBarry Smith nnz = PetscIntMultTruncate(2,PetscIntMultTruncate(nzi,(n - i))); /* estimated and max additional space needed */ 4615f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFreeSpaceGet(nnz,¤t_space)); 4625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFreeSpaceGet(nnz,¤t_space_lvl)); 46335aa4fcfSShri Abhyankar reallocs++; 46435aa4fcfSShri Abhyankar } 46535aa4fcfSShri Abhyankar 46635aa4fcfSShri Abhyankar /* copy data into free_space and free_space_lvl, then initialize lnk */ 4675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscIncompleteLLClean(n,n,nzi,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt)); 46826fbe8dcSKarl Rupp 46935aa4fcfSShri Abhyankar bj_ptr[i] = current_space->array; 47035aa4fcfSShri Abhyankar bjlvl_ptr[i] = current_space_lvl->array; 47135aa4fcfSShri Abhyankar 47235aa4fcfSShri Abhyankar /* make sure the active row i has diagonal entry */ 47394bad497SJacob Faibussowitsch PetscCheck(*(bj_ptr[i]+bdiag[i]) == i,PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Row %" PetscInt_FMT " has missing diagonal in factored matrix\ntry running with -pc_factor_nonzeros_along_diagonal or -pc_factor_diagonal_fill",i); 47435aa4fcfSShri Abhyankar 47535aa4fcfSShri Abhyankar current_space->array += nzi; 47635aa4fcfSShri Abhyankar current_space->local_used += nzi; 47735aa4fcfSShri Abhyankar current_space->local_remaining -= nzi; 47826fbe8dcSKarl Rupp 47935aa4fcfSShri Abhyankar current_space_lvl->array += nzi; 48035aa4fcfSShri Abhyankar current_space_lvl->local_used += nzi; 48135aa4fcfSShri Abhyankar current_space_lvl->local_remaining -= nzi; 48235aa4fcfSShri Abhyankar } 48335aa4fcfSShri Abhyankar 4845f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(isrow,&r)); 4855f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(isicol,&ic)); 48635aa4fcfSShri Abhyankar 48735aa4fcfSShri Abhyankar /* copy free_space into bj and free free_space; set bi, bj, bdiag in new datastructure; */ 4885f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(bi[n]+1,&bj)); 4895f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFreeSpaceContiguous_LU(&free_space,bj,n,bi,bdiag)); 49035aa4fcfSShri Abhyankar 4915f80ce2aSJacob Faibussowitsch CHKERRQ(PetscIncompleteLLDestroy(lnk,lnkbt)); 4925f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFreeSpaceDestroy(free_space_lvl)); 4935f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(bj_ptr,bjlvl_ptr)); 49435aa4fcfSShri Abhyankar 49535aa4fcfSShri Abhyankar #if defined(PETSC_USE_INFO) 49635aa4fcfSShri Abhyankar { 497aef85c9fSShri Abhyankar PetscReal af = ((PetscReal)(bdiag[0]+1))/((PetscReal)ai[n]); 4985f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(A,"Reallocs %" PetscInt_FMT " Fill ratio:given %g needed %g\n",reallocs,(double)f,(double)af)); 4995f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(A,"Run with -[sub_]pc_factor_fill %g or use \n",(double)af)); 5005f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(A,"PCFactorSetFill([sub]pc,%g);\n",(double)af)); 5015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(A,"for best performance.\n")); 50235aa4fcfSShri Abhyankar if (diagonal_fill) { 5035f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(A,"Detected and replaced %" PetscInt_FMT " missing diagonals\n",dcount)); 50435aa4fcfSShri Abhyankar } 50535aa4fcfSShri Abhyankar } 50635aa4fcfSShri Abhyankar #endif 50735aa4fcfSShri Abhyankar 50835aa4fcfSShri Abhyankar /* put together the new matrix */ 5095f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqBAIJSetPreallocation(fact,bs,MAT_SKIP_ALLOCATION,NULL)); 5105f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogObjectParent((PetscObject)fact,(PetscObject)isicol)); 51126fbe8dcSKarl Rupp 51235aa4fcfSShri Abhyankar b = (Mat_SeqBAIJ*)(fact)->data; 51335aa4fcfSShri Abhyankar b->free_a = PETSC_TRUE; 51435aa4fcfSShri Abhyankar b->free_ij = PETSC_TRUE; 51535aa4fcfSShri Abhyankar b->singlemalloc = PETSC_FALSE; 51626fbe8dcSKarl Rupp 5175f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(bs2*(bdiag[0]+1),&b->a)); 51826fbe8dcSKarl Rupp 51935aa4fcfSShri Abhyankar b->j = bj; 52035aa4fcfSShri Abhyankar b->i = bi; 52135aa4fcfSShri Abhyankar b->diag = bdiag; 52235aa4fcfSShri Abhyankar b->free_diag = PETSC_TRUE; 523f4259b30SLisandro Dalcin b->ilen = NULL; 524f4259b30SLisandro Dalcin b->imax = NULL; 52535aa4fcfSShri Abhyankar b->row = isrow; 52635aa4fcfSShri Abhyankar b->col = iscol; 5275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)isrow)); 5285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)iscol)); 52935aa4fcfSShri Abhyankar b->icol = isicol; 53026fbe8dcSKarl Rupp 5315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(bs*n+bs,&b->solve_work)); 53235aa4fcfSShri Abhyankar /* In b structure: Free imax, ilen, old a, old j. 53335aa4fcfSShri Abhyankar Allocate bdiag, solve_work, new a, new j */ 5345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogObjectMemory((PetscObject)fact,(bdiag[0]+1) * (sizeof(PetscInt)+bs2*sizeof(PetscScalar)))); 53535aa4fcfSShri Abhyankar b->maxnz = b->nz = bdiag[0]+1; 53626fbe8dcSKarl Rupp 537ae3d28f0SHong Zhang fact->info.factor_mallocs = reallocs; 538ae3d28f0SHong Zhang fact->info.fill_ratio_given = f; 539ae3d28f0SHong Zhang fact->info.fill_ratio_needed = ((PetscReal)(bdiag[0]+1))/((PetscReal)ai[n]); 54026fbe8dcSKarl Rupp 5415f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqBAIJSetNumericFactorization(fact,both_identity)); 54235aa4fcfSShri Abhyankar PetscFunctionReturn(0); 54335aa4fcfSShri Abhyankar } 54435aa4fcfSShri Abhyankar 5454e2b4712SSatish Balay /* 5464e2b4712SSatish Balay This code is virtually identical to MatILUFactorSymbolic_SeqAIJ 5474e2b4712SSatish Balay except that the data structure of Mat_SeqAIJ is slightly different. 5484e2b4712SSatish Balay Not a good example of code reuse. 5494e2b4712SSatish Balay */ 55006e38f1dSHong Zhang PetscErrorCode MatILUFactorSymbolic_SeqBAIJ_inplace(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info) 5514e2b4712SSatish Balay { 5524e2b4712SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b; 5534e2b4712SSatish Balay IS isicol; 5545d0c19d7SBarry Smith const PetscInt *r,*ic,*ai = a->i,*aj = a->j,*xi; 5555d0c19d7SBarry Smith PetscInt prow,n = a->mbs,*ainew,*ajnew,jmax,*fill,nz,*im,*ajfill,*flev,*xitmp; 556a96a251dSBarry Smith PetscInt *dloc,idx,row,m,fm,nzf,nzi,reallocate = 0,dcount = 0; 557d0f46423SBarry Smith PetscInt incrlev,nnz,i,bs = A->rmap->bs,bs2 = a->bs2,levels,diagonal_fill,dd; 558ace3abfcSBarry Smith PetscBool col_identity,row_identity,both_identity,flg; 559329f5518SBarry Smith PetscReal f; 5604e2b4712SSatish Balay 5614e2b4712SSatish Balay PetscFunctionBegin; 5625f80ce2aSJacob Faibussowitsch CHKERRQ(MatMissingDiagonal_SeqBAIJ(A,&flg,&dd)); 563*28b400f6SJacob Faibussowitsch PetscCheck(!flg,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix A is missing diagonal entry in row %" PetscInt_FMT,dd); 5646bce7ff8SHong Zhang 565435faa5fSBarry Smith f = info->fill; 566690b6cddSBarry Smith levels = (PetscInt)info->levels; 567690b6cddSBarry Smith diagonal_fill = (PetscInt)info->diagonal_fill; 56826fbe8dcSKarl Rupp 5695f80ce2aSJacob Faibussowitsch CHKERRQ(ISInvertPermutation(iscol,PETSC_DECIDE,&isicol)); 57016a2bf60SHong Zhang 5715f80ce2aSJacob Faibussowitsch CHKERRQ(ISIdentity(isrow,&row_identity)); 5725f80ce2aSJacob Faibussowitsch CHKERRQ(ISIdentity(iscol,&col_identity)); 573ace3abfcSBarry Smith both_identity = (PetscBool) (row_identity && col_identity); 574309c388cSBarry Smith 57541df41f0SMatthew Knepley if (!levels && both_identity) { /* special case copy the nonzero structure */ 5765f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicateNoCreate_SeqBAIJ(fact,A,MAT_DO_NOT_COPY_VALUES,PETSC_TRUE)); 5775f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqBAIJSetNumericFactorization_inplace(fact,both_identity)); 5786bce7ff8SHong Zhang 579d5f3da31SBarry Smith fact->factortype = MAT_FACTOR_ILU; 580ae3d28f0SHong Zhang b = (Mat_SeqBAIJ*)fact->data; 581bb3d539aSBarry Smith b->row = isrow; 582bb3d539aSBarry Smith b->col = iscol; 5835f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)isrow)); 5845f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)iscol)); 585bb3d539aSBarry Smith b->icol = isicol; 586bcd9e38bSBarry Smith b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE; 58726fbe8dcSKarl Rupp 5885f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1((n+1)*bs,&b->solve_work)); 5896bce7ff8SHong Zhang PetscFunctionReturn(0); 5906bce7ff8SHong Zhang } 5916bce7ff8SHong Zhang 5926bce7ff8SHong Zhang /* general case perform the symbolic factorization */ 5935f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(isrow,&r)); 5945f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(isicol,&ic)); 5954e2b4712SSatish Balay 5964e2b4712SSatish Balay /* get new row pointers */ 5975f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(n+1,&ainew)); 5984e2b4712SSatish Balay ainew[0] = 0; 5994e2b4712SSatish Balay /* don't know how many column pointers are needed so estimate */ 600690b6cddSBarry Smith jmax = (PetscInt)(f*ai[n] + 1); 6015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(jmax,&ajnew)); 6024e2b4712SSatish Balay /* ajfill is level of fill for each fill entry */ 6035f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(jmax,&ajfill)); 6044e2b4712SSatish Balay /* fill is a linked list of nonzeros in active row */ 6055f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(n+1,&fill)); 6064e2b4712SSatish Balay /* im is level for each filled value */ 6075f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(n+1,&im)); 6084e2b4712SSatish Balay /* dloc is location of diagonal in factor */ 6095f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(n+1,&dloc)); 6104e2b4712SSatish Balay dloc[0] = 0; 6114e2b4712SSatish Balay for (prow=0; prow<n; prow++) { 612435faa5fSBarry Smith 613435faa5fSBarry Smith /* copy prow into linked list */ 6144e2b4712SSatish Balay nzf = nz = ai[r[prow]+1] - ai[r[prow]]; 615*28b400f6SJacob Faibussowitsch PetscCheck(nz,PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %" PetscInt_FMT " in permuted ordering %" PetscInt_FMT,r[prow],prow); 6164e2b4712SSatish Balay xi = aj + ai[r[prow]]; 6174e2b4712SSatish Balay fill[n] = n; 618435faa5fSBarry Smith fill[prow] = -1; /* marker for diagonal entry */ 6194e2b4712SSatish Balay while (nz--) { 6204e2b4712SSatish Balay fm = n; 6214e2b4712SSatish Balay idx = ic[*xi++]; 6224e2b4712SSatish Balay do { 6234e2b4712SSatish Balay m = fm; 6244e2b4712SSatish Balay fm = fill[m]; 6254e2b4712SSatish Balay } while (fm < idx); 6264e2b4712SSatish Balay fill[m] = idx; 6274e2b4712SSatish Balay fill[idx] = fm; 6284e2b4712SSatish Balay im[idx] = 0; 6294e2b4712SSatish Balay } 630435faa5fSBarry Smith 631435faa5fSBarry Smith /* make sure diagonal entry is included */ 632435faa5fSBarry Smith if (diagonal_fill && fill[prow] == -1) { 633435faa5fSBarry Smith fm = n; 634435faa5fSBarry Smith while (fill[fm] < prow) fm = fill[fm]; 635435faa5fSBarry Smith fill[prow] = fill[fm]; /* insert diagonal into linked list */ 636435faa5fSBarry Smith fill[fm] = prow; 637435faa5fSBarry Smith im[prow] = 0; 638435faa5fSBarry Smith nzf++; 639335d9088SBarry Smith dcount++; 640435faa5fSBarry Smith } 641435faa5fSBarry Smith 6424e2b4712SSatish Balay nzi = 0; 6434e2b4712SSatish Balay row = fill[n]; 6444e2b4712SSatish Balay while (row < prow) { 6454e2b4712SSatish Balay incrlev = im[row] + 1; 6464e2b4712SSatish Balay nz = dloc[row]; 647435faa5fSBarry Smith xi = ajnew + ainew[row] + nz + 1; 6484e2b4712SSatish Balay flev = ajfill + ainew[row] + nz + 1; 6494e2b4712SSatish Balay nnz = ainew[row+1] - ainew[row] - nz - 1; 6504e2b4712SSatish Balay fm = row; 6514e2b4712SSatish Balay while (nnz-- > 0) { 6524e2b4712SSatish Balay idx = *xi++; 6534e2b4712SSatish Balay if (*flev + incrlev > levels) { 6544e2b4712SSatish Balay flev++; 6554e2b4712SSatish Balay continue; 6564e2b4712SSatish Balay } 6574e2b4712SSatish Balay do { 6584e2b4712SSatish Balay m = fm; 6594e2b4712SSatish Balay fm = fill[m]; 6604e2b4712SSatish Balay } while (fm < idx); 6614e2b4712SSatish Balay if (fm != idx) { 6624e2b4712SSatish Balay im[idx] = *flev + incrlev; 6634e2b4712SSatish Balay fill[m] = idx; 6644e2b4712SSatish Balay fill[idx] = fm; 6654e2b4712SSatish Balay fm = idx; 6664e2b4712SSatish Balay nzf++; 66726fbe8dcSKarl Rupp } else if (im[idx] > *flev + incrlev) im[idx] = *flev+incrlev; 6684e2b4712SSatish Balay flev++; 6694e2b4712SSatish Balay } 6704e2b4712SSatish Balay row = fill[row]; 6714e2b4712SSatish Balay nzi++; 6724e2b4712SSatish Balay } 6734e2b4712SSatish Balay /* copy new filled row into permanent storage */ 6744e2b4712SSatish Balay ainew[prow+1] = ainew[prow] + nzf; 6754e2b4712SSatish Balay if (ainew[prow+1] > jmax) { 676ecf371e4SBarry Smith 677ecf371e4SBarry Smith /* estimate how much additional space we will need */ 678ecf371e4SBarry Smith /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */ 679ecf371e4SBarry Smith /* just double the memory each time */ 680690b6cddSBarry Smith PetscInt maxadd = jmax; 681ecf371e4SBarry Smith /* maxadd = (int)(((f*ai[n]+1)*(n-prow+5))/n); */ 6824e2b4712SSatish Balay if (maxadd < nzf) maxadd = (n-prow)*(nzf+1); 6834e2b4712SSatish Balay jmax += maxadd; 684ecf371e4SBarry Smith 685ecf371e4SBarry Smith /* allocate a longer ajnew and ajfill */ 6865f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(jmax,&xitmp)); 6875f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(xitmp,ajnew,ainew[prow])); 6885f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(ajnew)); 6895d0c19d7SBarry Smith ajnew = xitmp; 6905f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(jmax,&xitmp)); 6915f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(xitmp,ajfill,ainew[prow])); 6925f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(ajfill)); 6935d0c19d7SBarry Smith ajfill = xitmp; 694eb150c5cSKris Buschelman reallocate++; /* count how many reallocations are needed */ 6954e2b4712SSatish Balay } 6965d0c19d7SBarry Smith xitmp = ajnew + ainew[prow]; 6974e2b4712SSatish Balay flev = ajfill + ainew[prow]; 6984e2b4712SSatish Balay dloc[prow] = nzi; 6994e2b4712SSatish Balay fm = fill[n]; 7004e2b4712SSatish Balay while (nzf--) { 7015d0c19d7SBarry Smith *xitmp++ = fm; 7024e2b4712SSatish Balay *flev++ = im[fm]; 7034e2b4712SSatish Balay fm = fill[fm]; 7044e2b4712SSatish Balay } 705435faa5fSBarry Smith /* make sure row has diagonal entry */ 7062c71b3e2SJacob Faibussowitsch PetscCheckFalse(ajnew[ainew[prow]+dloc[prow]] != prow,PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Row %" PetscInt_FMT " has missing diagonal in factored matrix\n\ 7072401956bSBarry Smith try running with -pc_factor_nonzeros_along_diagonal or -pc_factor_diagonal_fill",prow); 708435faa5fSBarry Smith } 7095f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(ajfill)); 7105f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(isrow,&r)); 7115f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(isicol,&ic)); 7125f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(fill)); 7135f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(im)); 7144e2b4712SSatish Balay 7156cf91177SBarry Smith #if defined(PETSC_USE_INFO) 7164e2b4712SSatish Balay { 717329f5518SBarry Smith PetscReal af = ((PetscReal)ainew[n])/((PetscReal)ai[n]); 7185f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(A,"Reallocs %" PetscInt_FMT " Fill ratio:given %g needed %g\n",reallocate,(double)f,(double)af)); 7195f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(A,"Run with -pc_factor_fill %g or use \n",(double)af)); 7205f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(A,"PCFactorSetFill(pc,%g);\n",(double)af)); 7215f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(A,"for best performance.\n")); 722335d9088SBarry Smith if (diagonal_fill) { 7235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(A,"Detected and replaced %" PetscInt_FMT " missing diagonals\n",dcount)); 724335d9088SBarry Smith } 7254e2b4712SSatish Balay } 72663ba0a88SBarry Smith #endif 7274e2b4712SSatish Balay 7284e2b4712SSatish Balay /* put together the new matrix */ 7295f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqBAIJSetPreallocation(fact,bs,MAT_SKIP_ALLOCATION,NULL)); 7305f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogObjectParent((PetscObject)fact,(PetscObject)isicol)); 731ae3d28f0SHong Zhang b = (Mat_SeqBAIJ*)fact->data; 73226fbe8dcSKarl Rupp 733e6b907acSBarry Smith b->free_a = PETSC_TRUE; 734e6b907acSBarry Smith b->free_ij = PETSC_TRUE; 7357c922b88SBarry Smith b->singlemalloc = PETSC_FALSE; 73626fbe8dcSKarl Rupp 7375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(bs2*ainew[n],&b->a)); 73826fbe8dcSKarl Rupp 7394e2b4712SSatish Balay b->j = ajnew; 7404e2b4712SSatish Balay b->i = ainew; 7414e2b4712SSatish Balay for (i=0; i<n; i++) dloc[i] += ainew[i]; 7424e2b4712SSatish Balay b->diag = dloc; 7437f53bb6cSHong Zhang b->free_diag = PETSC_TRUE; 744f4259b30SLisandro Dalcin b->ilen = NULL; 745f4259b30SLisandro Dalcin b->imax = NULL; 7464e2b4712SSatish Balay b->row = isrow; 7474e2b4712SSatish Balay b->col = iscol; 748bcd9e38bSBarry Smith b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE; 74926fbe8dcSKarl Rupp 7505f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)isrow)); 7515f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)iscol)); 752e51c0b9cSSatish Balay b->icol = isicol; 7535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(bs*n+bs,&b->solve_work)); 7544e2b4712SSatish Balay /* In b structure: Free imax, ilen, old a, old j. 7554e2b4712SSatish Balay Allocate dloc, solve_work, new a, new j */ 7565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogObjectMemory((PetscObject)fact,(ainew[n]-n)*(sizeof(PetscInt))+bs2*ainew[n]*sizeof(PetscScalar))); 7574e2b4712SSatish Balay b->maxnz = b->nz = ainew[n]; 7584e2b4712SSatish Balay 759ae3d28f0SHong Zhang fact->info.factor_mallocs = reallocate; 760ae3d28f0SHong Zhang fact->info.fill_ratio_given = f; 761ae3d28f0SHong Zhang fact->info.fill_ratio_needed = ((PetscReal)ainew[n])/((PetscReal)ai[prow]); 7626bce7ff8SHong Zhang 7635f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqBAIJSetNumericFactorization_inplace(fact,both_identity)); 7648661488fSKris Buschelman PetscFunctionReturn(0); 7658661488fSKris Buschelman } 7668661488fSKris Buschelman 767dfbe8321SBarry Smith PetscErrorCode MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE(Mat A) 7687e7071cdSKris Buschelman { 76912272027SHong Zhang /* Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; */ 77012272027SHong Zhang /* int i,*AJ=a->j,nz=a->nz; */ 7715fd66863SKarl Rupp 7725a9542e3SKris Buschelman PetscFunctionBegin; 7737cf1b8d3SKris Buschelman /* Undo Column scaling */ 7747cf1b8d3SKris Buschelman /* while (nz--) { */ 7757cf1b8d3SKris Buschelman /* AJ[i] = AJ[i]/4; */ 7767cf1b8d3SKris Buschelman /* } */ 777c115a38dSKris Buschelman /* This should really invoke a push/pop logic, but we don't have that yet. */ 7780298fd71SBarry Smith A->ops->setunfactored = NULL; 7797cf1b8d3SKris Buschelman PetscFunctionReturn(0); 7807cf1b8d3SKris Buschelman } 7817cf1b8d3SKris Buschelman 782dfbe8321SBarry Smith PetscErrorCode MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE_usj(Mat A) 7837cf1b8d3SKris Buschelman { 7847cf1b8d3SKris Buschelman Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 785b24ad042SBarry Smith PetscInt *AJ=a->j,nz=a->nz; 7862aa5897fSKris Buschelman unsigned short *aj=(unsigned short*)AJ; 7875fd66863SKarl Rupp 7885a9542e3SKris Buschelman PetscFunctionBegin; 7890b9da03eSKris Buschelman /* Is this really necessary? */ 79020235379SKris Buschelman while (nz--) { 7910b9da03eSKris Buschelman AJ[nz] = (int)((unsigned int)aj[nz]); /* First extend, then convert to signed. */ 7927e7071cdSKris Buschelman } 7930298fd71SBarry Smith A->ops->setunfactored = NULL; 7947e7071cdSKris Buschelman PetscFunctionReturn(0); 7957e7071cdSKris Buschelman } 796