1*be1d678aSKris Buschelman #define PETSCMAT_DLL 2*be1d678aSKris Buschelman 35d517e7eSBarry Smith /* 4ec3a8b7bSBarry Smith Factorization code for BAIJ format. 55d517e7eSBarry Smith */ 670f55243SBarry Smith #include "src/mat/impls/baij/seq/baij.h" 74078e994SBarry Smith #include "src/inline/ilu.h" 8ec3a8b7bSBarry Smith 96d84be18SBarry Smith /* ------------------------------------------------------------*/ 106d84be18SBarry Smith /* 114eeb42bcSBarry Smith Version for when blocks are 2 by 2 124eeb42bcSBarry Smith */ 134a2ae208SSatish Balay #undef __FUNCT__ 144a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2" 15af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2(Mat A,MatFactorInfo *info,Mat *B) 164eeb42bcSBarry Smith { 174eeb42bcSBarry Smith Mat C = *B; 184eeb42bcSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; 197cdf2dbbSSatish Balay IS isrow = b->row,isicol = b->icol; 206849ba73SBarry Smith PetscErrorCode ierr; 21690b6cddSBarry Smith PetscInt *r,*ic,i,j,n = a->mbs,*bi = b->i,*bj = b->j; 22690b6cddSBarry Smith PetscInt *ajtmpold,*ajtmp,nz,row; 23690b6cddSBarry Smith PetscInt *diag_offset=b->diag,idx,*ai=a->i,*aj=a->j,*pj; 24329f5518SBarry Smith MatScalar *pv,*v,*rtmp,m1,m2,m3,m4,*pc,*w,*x,x1,x2,x3,x4; 252515f8d2SSatish Balay MatScalar p1,p2,p3,p4; 263f1db9ecSBarry Smith MatScalar *ba = b->a,*aa = a->a; 274eeb42bcSBarry Smith 283a40ed3dSBarry Smith PetscFunctionBegin; 294eeb42bcSBarry Smith ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 304eeb42bcSBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 31b0a32e0cSBarry Smith ierr = PetscMalloc(4*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 324eeb42bcSBarry Smith 334eeb42bcSBarry Smith for (i=0; i<n; i++) { 344078e994SBarry Smith nz = bi[i+1] - bi[i]; 354078e994SBarry Smith ajtmp = bj + bi[i]; 364eeb42bcSBarry Smith for (j=0; j<nz; j++) { 374eeb42bcSBarry Smith x = rtmp+4*ajtmp[j]; x[0] = x[1] = x[2] = x[3] = 0.0; 384eeb42bcSBarry Smith } 394eeb42bcSBarry Smith /* load in initial (unfactored row) */ 404eeb42bcSBarry Smith idx = r[i]; 414078e994SBarry Smith nz = ai[idx+1] - ai[idx]; 424078e994SBarry Smith ajtmpold = aj + ai[idx]; 434078e994SBarry Smith v = aa + 4*ai[idx]; 444eeb42bcSBarry Smith for (j=0; j<nz; j++) { 454eeb42bcSBarry Smith x = rtmp+4*ic[ajtmpold[j]]; 464eeb42bcSBarry Smith x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3]; 474eeb42bcSBarry Smith v += 4; 484eeb42bcSBarry Smith } 494eeb42bcSBarry Smith row = *ajtmp++; 504eeb42bcSBarry Smith while (row < i) { 514eeb42bcSBarry Smith pc = rtmp + 4*row; 524eeb42bcSBarry Smith p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3]; 5388685aaeSLois Curfman McInnes if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0) { 544078e994SBarry Smith pv = ba + 4*diag_offset[row]; 554078e994SBarry Smith pj = bj + diag_offset[row] + 1; 564eeb42bcSBarry Smith x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 574eeb42bcSBarry Smith pc[0] = m1 = p1*x1 + p3*x2; 584eeb42bcSBarry Smith pc[1] = m2 = p2*x1 + p4*x2; 594eeb42bcSBarry Smith pc[2] = m3 = p1*x3 + p3*x4; 604eeb42bcSBarry Smith pc[3] = m4 = p2*x3 + p4*x4; 614078e994SBarry Smith nz = bi[row+1] - diag_offset[row] - 1; 624eeb42bcSBarry Smith pv += 4; 634eeb42bcSBarry Smith for (j=0; j<nz; j++) { 644eeb42bcSBarry Smith x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 654eeb42bcSBarry Smith x = rtmp + 4*pj[j]; 664eeb42bcSBarry Smith x[0] -= m1*x1 + m3*x2; 674eeb42bcSBarry Smith x[1] -= m2*x1 + m4*x2; 684eeb42bcSBarry Smith x[2] -= m1*x3 + m3*x4; 694eeb42bcSBarry Smith x[3] -= m2*x3 + m4*x4; 704eeb42bcSBarry Smith pv += 4; 714eeb42bcSBarry Smith } 72efee365bSSatish Balay ierr = PetscLogFlops(16*nz+12);CHKERRQ(ierr); 734eeb42bcSBarry Smith } 744eeb42bcSBarry Smith row = *ajtmp++; 754eeb42bcSBarry Smith } 764eeb42bcSBarry Smith /* finished row so stick it into b->a */ 774078e994SBarry Smith pv = ba + 4*bi[i]; 784078e994SBarry Smith pj = bj + bi[i]; 794078e994SBarry Smith nz = bi[i+1] - bi[i]; 804eeb42bcSBarry Smith for (j=0; j<nz; j++) { 814eeb42bcSBarry Smith x = rtmp+4*pj[j]; 824eeb42bcSBarry Smith pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3]; 834eeb42bcSBarry Smith pv += 4; 844eeb42bcSBarry Smith } 854eeb42bcSBarry Smith /* invert diagonal block */ 864078e994SBarry Smith w = ba + 4*diag_offset[i]; 874cd38560SBarry Smith ierr = Kernel_A_gets_inverse_A_2(w);CHKERRQ(ierr); 884eeb42bcSBarry Smith } 894eeb42bcSBarry Smith 90606d414cSSatish Balay ierr = PetscFree(rtmp);CHKERRQ(ierr); 914eeb42bcSBarry Smith ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 924eeb42bcSBarry Smith ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 934eeb42bcSBarry Smith C->factor = FACTOR_LU; 944eeb42bcSBarry Smith C->assembled = PETSC_TRUE; 95efee365bSSatish Balay ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 963a40ed3dSBarry Smith PetscFunctionReturn(0); 974eeb42bcSBarry Smith } 984cd38560SBarry Smith /* 994cd38560SBarry Smith Version for when blocks are 2 by 2 Using natural ordering 1004cd38560SBarry Smith */ 1014a2ae208SSatish Balay #undef __FUNCT__ 1024a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering" 103af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering(Mat A,MatFactorInfo *info,Mat *B) 1044cd38560SBarry Smith { 1054cd38560SBarry Smith Mat C = *B; 1064cd38560SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; 107dfbe8321SBarry Smith PetscErrorCode ierr; 108690b6cddSBarry Smith PetscInt i,j,n = a->mbs,*bi = b->i,*bj = b->j; 109690b6cddSBarry Smith PetscInt *ajtmpold,*ajtmp,nz,row; 110690b6cddSBarry Smith PetscInt *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj; 111329f5518SBarry Smith MatScalar *pv,*v,*rtmp,*pc,*w,*x; 1122515f8d2SSatish Balay MatScalar p1,p2,p3,p4,m1,m2,m3,m4,x1,x2,x3,x4; 1134cd38560SBarry Smith MatScalar *ba = b->a,*aa = a->a; 1144cd38560SBarry Smith 1154cd38560SBarry Smith PetscFunctionBegin; 116b0a32e0cSBarry Smith ierr = PetscMalloc(4*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 1174cd38560SBarry Smith 1184cd38560SBarry Smith for (i=0; i<n; i++) { 1194cd38560SBarry Smith nz = bi[i+1] - bi[i]; 1204cd38560SBarry Smith ajtmp = bj + bi[i]; 1214cd38560SBarry Smith for (j=0; j<nz; j++) { 1224cd38560SBarry Smith x = rtmp+4*ajtmp[j]; 1234cd38560SBarry Smith x[0] = x[1] = x[2] = x[3] = 0.0; 1244cd38560SBarry Smith } 1254cd38560SBarry Smith /* load in initial (unfactored row) */ 1264cd38560SBarry Smith nz = ai[i+1] - ai[i]; 1274cd38560SBarry Smith ajtmpold = aj + ai[i]; 1284cd38560SBarry Smith v = aa + 4*ai[i]; 1294cd38560SBarry Smith for (j=0; j<nz; j++) { 1304cd38560SBarry Smith x = rtmp+4*ajtmpold[j]; 1314cd38560SBarry Smith x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3]; 1324cd38560SBarry Smith v += 4; 1334cd38560SBarry Smith } 1344cd38560SBarry Smith row = *ajtmp++; 1354cd38560SBarry Smith while (row < i) { 1364cd38560SBarry Smith pc = rtmp + 4*row; 1374cd38560SBarry Smith p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3]; 1384cd38560SBarry Smith if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0) { 1394cd38560SBarry Smith pv = ba + 4*diag_offset[row]; 1404cd38560SBarry Smith pj = bj + diag_offset[row] + 1; 1414cd38560SBarry Smith x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 1424cd38560SBarry Smith pc[0] = m1 = p1*x1 + p3*x2; 1434cd38560SBarry Smith pc[1] = m2 = p2*x1 + p4*x2; 1444cd38560SBarry Smith pc[2] = m3 = p1*x3 + p3*x4; 1454cd38560SBarry Smith pc[3] = m4 = p2*x3 + p4*x4; 1464cd38560SBarry Smith nz = bi[row+1] - diag_offset[row] - 1; 1474cd38560SBarry Smith pv += 4; 1484cd38560SBarry Smith for (j=0; j<nz; j++) { 1494cd38560SBarry Smith x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 1504cd38560SBarry Smith x = rtmp + 4*pj[j]; 1514cd38560SBarry Smith x[0] -= m1*x1 + m3*x2; 1524cd38560SBarry Smith x[1] -= m2*x1 + m4*x2; 1534cd38560SBarry Smith x[2] -= m1*x3 + m3*x4; 1544cd38560SBarry Smith x[3] -= m2*x3 + m4*x4; 1554cd38560SBarry Smith pv += 4; 1564cd38560SBarry Smith } 157efee365bSSatish Balay ierr = PetscLogFlops(16*nz+12);CHKERRQ(ierr); 1584cd38560SBarry Smith } 1594cd38560SBarry Smith row = *ajtmp++; 1604cd38560SBarry Smith } 1614cd38560SBarry Smith /* finished row so stick it into b->a */ 1624cd38560SBarry Smith pv = ba + 4*bi[i]; 1634cd38560SBarry Smith pj = bj + bi[i]; 1644cd38560SBarry Smith nz = bi[i+1] - bi[i]; 1654cd38560SBarry Smith for (j=0; j<nz; j++) { 1664cd38560SBarry Smith x = rtmp+4*pj[j]; 1674cd38560SBarry Smith pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3]; 1684cd38560SBarry Smith pv += 4; 1694cd38560SBarry Smith } 1704cd38560SBarry Smith /* invert diagonal block */ 1714cd38560SBarry Smith w = ba + 4*diag_offset[i]; 1724cd38560SBarry Smith ierr = Kernel_A_gets_inverse_A_2(w);CHKERRQ(ierr); 1734cd38560SBarry Smith /*Kernel_A_gets_inverse_A(bs,w,v_pivots,v_work);*/ 1744cd38560SBarry Smith } 1754cd38560SBarry Smith 176606d414cSSatish Balay ierr = PetscFree(rtmp);CHKERRQ(ierr); 1774cd38560SBarry Smith C->factor = FACTOR_LU; 1784cd38560SBarry Smith C->assembled = PETSC_TRUE; 179efee365bSSatish Balay ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 1804cd38560SBarry Smith PetscFunctionReturn(0); 1814cd38560SBarry Smith } 1827fc0212eSBarry Smith 1837fc0212eSBarry Smith /* ----------------------------------------------------------- */ 1847fc0212eSBarry Smith /* 1857fc0212eSBarry Smith Version for when blocks are 1 by 1. 1867fc0212eSBarry Smith */ 1874a2ae208SSatish Balay #undef __FUNCT__ 1884a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_1" 189af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqBAIJ_1(Mat A,MatFactorInfo *info,Mat *B) 1907fc0212eSBarry Smith { 1917fc0212eSBarry Smith Mat C = *B; 1927fc0212eSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; 1937cdf2dbbSSatish Balay IS isrow = b->row,isicol = b->icol; 1946849ba73SBarry Smith PetscErrorCode ierr; 195690b6cddSBarry Smith PetscInt *r,*ic,i,j,n = a->mbs,*bi = b->i,*bj = b->j; 196690b6cddSBarry Smith PetscInt *ajtmpold,*ajtmp,nz,row,*ai = a->i,*aj = a->j; 197690b6cddSBarry Smith PetscInt *diag_offset = b->diag,diag,*pj; 198329f5518SBarry Smith MatScalar *pv,*v,*rtmp,multiplier,*pc; 1993f1db9ecSBarry Smith MatScalar *ba = b->a,*aa = a->a; 2007fc0212eSBarry Smith 2013a40ed3dSBarry Smith PetscFunctionBegin; 2027fc0212eSBarry Smith ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 2037fc0212eSBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 204b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 2057fc0212eSBarry Smith 2067fc0212eSBarry Smith for (i=0; i<n; i++) { 2074078e994SBarry Smith nz = bi[i+1] - bi[i]; 2084078e994SBarry Smith ajtmp = bj + bi[i]; 2097fc0212eSBarry Smith for (j=0; j<nz; j++) rtmp[ajtmp[j]] = 0.0; 2107fc0212eSBarry Smith 2117fc0212eSBarry Smith /* load in initial (unfactored row) */ 2124078e994SBarry Smith nz = ai[r[i]+1] - ai[r[i]]; 2134078e994SBarry Smith ajtmpold = aj + ai[r[i]]; 2144078e994SBarry Smith v = aa + ai[r[i]]; 2157fc0212eSBarry Smith for (j=0; j<nz; j++) rtmp[ic[ajtmpold[j]]] = v[j]; 2167fc0212eSBarry Smith 2177fc0212eSBarry Smith row = *ajtmp++; 2187fc0212eSBarry Smith while (row < i) { 2197fc0212eSBarry Smith pc = rtmp + row; 2207fc0212eSBarry Smith if (*pc != 0.0) { 2214078e994SBarry Smith pv = ba + diag_offset[row]; 2224078e994SBarry Smith pj = bj + diag_offset[row] + 1; 2237fc0212eSBarry Smith multiplier = *pc * *pv++; 2247fc0212eSBarry Smith *pc = multiplier; 2254078e994SBarry Smith nz = bi[row+1] - diag_offset[row] - 1; 2267fc0212eSBarry Smith for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j]; 227efee365bSSatish Balay ierr = PetscLogFlops(1+2*nz);CHKERRQ(ierr); 2287fc0212eSBarry Smith } 2297fc0212eSBarry Smith row = *ajtmp++; 2307fc0212eSBarry Smith } 2317fc0212eSBarry Smith /* finished row so stick it into b->a */ 2324078e994SBarry Smith pv = ba + bi[i]; 2334078e994SBarry Smith pj = bj + bi[i]; 2344078e994SBarry Smith nz = bi[i+1] - bi[i]; 2357fc0212eSBarry Smith for (j=0; j<nz; j++) {pv[j] = rtmp[pj[j]];} 2364078e994SBarry Smith diag = diag_offset[i] - bi[i]; 2377fc0212eSBarry Smith /* check pivot entry for current row */ 238a73cf429SBarry Smith if (pv[diag] == 0.0) { 23929bbc08cSBarry Smith SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot"); 2407fc0212eSBarry Smith } 2417fc0212eSBarry Smith pv[diag] = 1.0/pv[diag]; 2427fc0212eSBarry Smith } 2437fc0212eSBarry Smith 244606d414cSSatish Balay ierr = PetscFree(rtmp);CHKERRQ(ierr); 2457fc0212eSBarry Smith ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 2467fc0212eSBarry Smith ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 2477fc0212eSBarry Smith C->factor = FACTOR_LU; 2487fc0212eSBarry Smith C->assembled = PETSC_TRUE; 249efee365bSSatish Balay ierr = PetscLogFlops(C->n);CHKERRQ(ierr); 2503a40ed3dSBarry Smith PetscFunctionReturn(0); 2517fc0212eSBarry Smith } 2527fc0212eSBarry Smith 253273d9f13SBarry Smith 2545d517e7eSBarry Smith /* ----------------------------------------------------------- */ 2554a2ae208SSatish Balay #undef __FUNCT__ 2564a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqBAIJ" 257dfbe8321SBarry Smith PetscErrorCode MatLUFactor_SeqBAIJ(Mat A,IS row,IS col,MatFactorInfo *info) 2585d517e7eSBarry Smith { 259dfbe8321SBarry Smith PetscErrorCode ierr; 2605d517e7eSBarry Smith Mat C; 2615d517e7eSBarry Smith 2623a40ed3dSBarry Smith PetscFunctionBegin; 263b9b97703SBarry Smith ierr = MatLUFactorSymbolic(A,row,col,info,&C);CHKERRQ(ierr); 264af281ebdSHong Zhang ierr = MatLUFactorNumeric(A,info,&C);CHKERRQ(ierr); 265273d9f13SBarry Smith ierr = MatHeaderCopy(A,C);CHKERRQ(ierr); 26652e6d16bSBarry Smith ierr = PetscLogObjectParent(A,((Mat_SeqBAIJ*)(A->data))->icol);CHKERRQ(ierr); 2673a40ed3dSBarry Smith PetscFunctionReturn(0); 2685d517e7eSBarry Smith } 2694cd38560SBarry Smith 270c05c3958SHong Zhang #include "src/mat/impls/sbaij/seq/sbaij.h" 271c05c3958SHong Zhang #undef __FUNCT__ 272c05c3958SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N" 273af281ebdSHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N(Mat A,MatFactorInfo *info,Mat *B) 274c05c3958SHong Zhang { 275f3086b4bSHong Zhang PetscErrorCode ierr; 27678910aadSHong Zhang Mat C = *B; 27778910aadSHong Zhang Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data; 27878910aadSHong Zhang Mat_SeqSBAIJ *b=(Mat_SeqSBAIJ*)C->data; 27978910aadSHong Zhang IS ip=b->row; 28078910aadSHong Zhang PetscInt *rip,i,j,mbs=a->mbs,bs=A->bs,*bi=b->i,*bj=b->j,*bcol; 28178910aadSHong Zhang PetscInt *ai=a->i,*aj=a->j; 28278910aadSHong Zhang PetscInt k,jmin,jmax,*jl,*il,col,nexti,ili,nz; 28378910aadSHong Zhang MatScalar *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi; 2843cea4cbeSHong Zhang PetscReal zeropivot,rs,shiftnz; 2853cea4cbeSHong Zhang PetscTruth shiftpd; 2863cea4cbeSHong Zhang ChShift_Ctx sctx; 2873cea4cbeSHong Zhang PetscInt newshift; 28878910aadSHong Zhang 289c05c3958SHong Zhang PetscFunctionBegin; 2906ad2eaddSHong Zhang if (bs > 1) { 2916ad2eaddSHong Zhang if (!a->sbaijMat){ 292ceb03754SKris Buschelman ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr); 2936ad2eaddSHong Zhang } 294f3086b4bSHong Zhang ierr = (a->sbaijMat)->ops->choleskyfactornumeric(a->sbaijMat,info,B);CHKERRQ(ierr); 2956ad2eaddSHong Zhang ierr = MatDestroy(a->sbaijMat);CHKERRQ(ierr); 2966ad2eaddSHong Zhang a->sbaijMat = PETSC_NULL; 2976ad2eaddSHong Zhang PetscFunctionReturn(0); 2986ad2eaddSHong Zhang } 29978910aadSHong Zhang 30078910aadSHong Zhang /* initialization */ 3013cea4cbeSHong Zhang shiftnz = info->shiftnz; 3023cea4cbeSHong Zhang shiftpd = info->shiftpd; 3033cea4cbeSHong Zhang zeropivot = info->zeropivot; 3043cea4cbeSHong Zhang 3056ad2eaddSHong Zhang ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr); 30678910aadSHong Zhang nz = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar); 30778910aadSHong Zhang ierr = PetscMalloc(nz,&il);CHKERRQ(ierr); 30878910aadSHong Zhang jl = il + mbs; 30978910aadSHong Zhang rtmp = (MatScalar*)(jl + mbs); 31078910aadSHong Zhang 3113cea4cbeSHong Zhang sctx.shift_amount = 0; 3123cea4cbeSHong Zhang sctx.nshift = 0; 31378910aadSHong Zhang do { 3143cea4cbeSHong Zhang sctx.chshift = PETSC_FALSE; 31578910aadSHong Zhang for (i=0; i<mbs; i++) { 31678910aadSHong Zhang rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0; 31778910aadSHong Zhang } 31878910aadSHong Zhang 31978910aadSHong Zhang for (k = 0; k<mbs; k++){ 32078910aadSHong Zhang bval = ba + bi[k]; 32178910aadSHong Zhang /* initialize k-th row by the perm[k]-th row of A */ 32278910aadSHong Zhang jmin = ai[rip[k]]; jmax = ai[rip[k]+1]; 32378910aadSHong Zhang for (j = jmin; j < jmax; j++){ 32478910aadSHong Zhang col = rip[aj[j]]; 32578910aadSHong Zhang if (col >= k){ /* only take upper triangular entry */ 32678910aadSHong Zhang rtmp[col] = aa[j]; 32778910aadSHong Zhang *bval++ = 0.0; /* for in-place factorization */ 32878910aadSHong Zhang } 32978910aadSHong Zhang } 3303cea4cbeSHong Zhang 3313cea4cbeSHong Zhang /* shift the diagonal of the matrix */ 3323cea4cbeSHong Zhang if (sctx.nshift) rtmp[k] += sctx.shift_amount; 33378910aadSHong Zhang 33478910aadSHong Zhang /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 33578910aadSHong Zhang dk = rtmp[k]; 33678910aadSHong Zhang i = jl[k]; /* first row to be added to k_th row */ 33778910aadSHong Zhang 33878910aadSHong Zhang while (i < k){ 33978910aadSHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 34078910aadSHong Zhang 34178910aadSHong Zhang /* compute multiplier, update diag(k) and U(i,k) */ 34278910aadSHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 34378910aadSHong Zhang uikdi = - ba[ili]*ba[bi[i]]; /* diagonal(k) */ 34478910aadSHong Zhang dk += uikdi*ba[ili]; 34578910aadSHong Zhang ba[ili] = uikdi; /* -U(i,k) */ 34678910aadSHong Zhang 34778910aadSHong Zhang /* add multiple of row i to k-th row */ 34878910aadSHong Zhang jmin = ili + 1; jmax = bi[i+1]; 34978910aadSHong Zhang if (jmin < jmax){ 35078910aadSHong Zhang for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j]; 35178910aadSHong Zhang /* update il and jl for row i */ 35278910aadSHong Zhang il[i] = jmin; 35378910aadSHong Zhang j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; 35478910aadSHong Zhang } 35578910aadSHong Zhang i = nexti; 35678910aadSHong Zhang } 35778910aadSHong Zhang 3583cea4cbeSHong Zhang /* shift the diagonals when zero pivot is detected */ 3593cea4cbeSHong Zhang /* compute rs=sum of abs(off-diagonal) */ 3603cea4cbeSHong Zhang rs = 0.0; 36178910aadSHong Zhang jmin = bi[k]+1; 36278910aadSHong Zhang nz = bi[k+1] - jmin; 36378910aadSHong Zhang if (nz){ 36478910aadSHong Zhang bcol = bj + jmin; 36578910aadSHong Zhang while (nz--){ 36689f845c8SHong Zhang rs += PetscAbsScalar(rtmp[*bcol]); 36789f845c8SHong Zhang bcol++; 36878910aadSHong Zhang } 36978910aadSHong Zhang } 3703cea4cbeSHong Zhang 3713cea4cbeSHong Zhang sctx.rs = rs; 3723cea4cbeSHong Zhang sctx.pv = dk; 3733e8c821fSHong Zhang ierr = MatCholeskyCheckShift_inline(info,sctx,newshift);CHKERRQ(ierr); 3743cea4cbeSHong Zhang if (newshift == 1){ 3753cea4cbeSHong Zhang break; /* sctx.shift_amount is updated */ 3763cea4cbeSHong Zhang } else if (newshift == -1){ 3773cea4cbeSHong Zhang SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g * rs %g",k,PetscAbsScalar(dk),zeropivot,rs); 37878910aadSHong Zhang } 37978910aadSHong Zhang 38078910aadSHong Zhang /* copy data into U(k,:) */ 38178910aadSHong Zhang ba[bi[k]] = 1.0/dk; /* U(k,k) */ 38278910aadSHong Zhang jmin = bi[k]+1; jmax = bi[k+1]; 38378910aadSHong Zhang if (jmin < jmax) { 38478910aadSHong Zhang for (j=jmin; j<jmax; j++){ 38578910aadSHong Zhang col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0; 38678910aadSHong Zhang } 38778910aadSHong Zhang /* add the k-th row into il and jl */ 38878910aadSHong Zhang il[k] = jmin; 38978910aadSHong Zhang i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; 39078910aadSHong Zhang } 39178910aadSHong Zhang } 3923cea4cbeSHong Zhang } while (sctx.chshift); 39378910aadSHong Zhang ierr = PetscFree(il);CHKERRQ(ierr); 39478910aadSHong Zhang 39578910aadSHong Zhang ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr); 39678910aadSHong Zhang C->factor = FACTOR_CHOLESKY; 39778910aadSHong Zhang C->assembled = PETSC_TRUE; 39878910aadSHong Zhang C->preallocated = PETSC_TRUE; 399efee365bSSatish Balay ierr = PetscLogFlops(C->m);CHKERRQ(ierr); 4003cea4cbeSHong Zhang if (sctx.nshift){ 4013cea4cbeSHong Zhang if (shiftnz) { 40263ba0a88SBarry Smith ierr = PetscLogInfo((0,"MatCholeskyFactorNumeric_SeqBAIJ_1: number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount));CHKERRQ(ierr); 4033cea4cbeSHong Zhang } else if (shiftpd) { 40463ba0a88SBarry Smith ierr = PetscLogInfo((0,"MatCholeskyFactorNumeric_SeqBAIJ_1: number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount));CHKERRQ(ierr); 40578910aadSHong Zhang } 40678910aadSHong Zhang } 407c05c3958SHong Zhang PetscFunctionReturn(0); 408c05c3958SHong Zhang } 4094cd38560SBarry Smith 410c05c3958SHong Zhang #undef __FUNCT__ 411c05c3958SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering" 412af281ebdSHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering(Mat A,MatFactorInfo *info,Mat *fact) 413c05c3958SHong Zhang { 41478910aadSHong Zhang Mat C = *fact; 41578910aadSHong Zhang Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data; 41678910aadSHong Zhang Mat_SeqSBAIJ *b=(Mat_SeqSBAIJ*)C->data; 41778910aadSHong Zhang PetscErrorCode ierr; 4183cea4cbeSHong Zhang PetscInt i,j,am=a->mbs; 41978910aadSHong Zhang PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 4203cea4cbeSHong Zhang PetscInt k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz; 42178910aadSHong Zhang MatScalar *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval; 4223cea4cbeSHong Zhang PetscReal zeropivot,rs,shiftnz; 4233cea4cbeSHong Zhang PetscTruth shiftpd; 4243cea4cbeSHong Zhang ChShift_Ctx sctx; 4253cea4cbeSHong Zhang PetscInt newshift; 42678910aadSHong Zhang 427c05c3958SHong Zhang PetscFunctionBegin; 42878910aadSHong Zhang /* initialization */ 4293cea4cbeSHong Zhang shiftnz = info->shiftnz; 4303cea4cbeSHong Zhang shiftpd = info->shiftpd; 4313cea4cbeSHong Zhang zeropivot = info->zeropivot; 4323cea4cbeSHong Zhang 43378910aadSHong Zhang nz = (2*am+1)*sizeof(PetscInt)+am*sizeof(MatScalar); 43478910aadSHong Zhang ierr = PetscMalloc(nz,&il);CHKERRQ(ierr); 43578910aadSHong Zhang jl = il + am; 43678910aadSHong Zhang rtmp = (MatScalar*)(jl + am); 43778910aadSHong Zhang 4383cea4cbeSHong Zhang sctx.shift_amount = 0; 4393cea4cbeSHong Zhang sctx.nshift = 0; 44078910aadSHong Zhang do { 4413cea4cbeSHong Zhang sctx.chshift = PETSC_FALSE; 44278910aadSHong Zhang for (i=0; i<am; i++) { 44378910aadSHong Zhang rtmp[i] = 0.0; jl[i] = am; il[0] = 0; 44478910aadSHong Zhang } 44578910aadSHong Zhang 44678910aadSHong Zhang for (k = 0; k<am; k++){ 44778910aadSHong Zhang /* initialize k-th row with elements nonzero in row perm(k) of A */ 44878910aadSHong Zhang nz = ai[k+1] - ai[k]; 44978910aadSHong Zhang acol = aj + ai[k]; 45078910aadSHong Zhang aval = aa + ai[k]; 45178910aadSHong Zhang bval = ba + bi[k]; 45278910aadSHong Zhang while (nz -- ){ 45378910aadSHong Zhang if (*acol < k) { /* skip lower triangular entries */ 45478910aadSHong Zhang acol++; aval++; 45578910aadSHong Zhang } else { 45678910aadSHong Zhang rtmp[*acol++] = *aval++; 45778910aadSHong Zhang *bval++ = 0.0; /* for in-place factorization */ 45878910aadSHong Zhang } 45978910aadSHong Zhang } 4603cea4cbeSHong Zhang 4613cea4cbeSHong Zhang /* shift the diagonal of the matrix */ 4623cea4cbeSHong Zhang if (sctx.nshift) rtmp[k] += sctx.shift_amount; 46378910aadSHong Zhang 46478910aadSHong Zhang /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 46578910aadSHong Zhang dk = rtmp[k]; 46678910aadSHong Zhang i = jl[k]; /* first row to be added to k_th row */ 46778910aadSHong Zhang 46878910aadSHong Zhang while (i < k){ 46978910aadSHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 47078910aadSHong Zhang /* compute multiplier, update D(k) and U(i,k) */ 47178910aadSHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 47278910aadSHong Zhang uikdi = - ba[ili]*ba[bi[i]]; 47378910aadSHong Zhang dk += uikdi*ba[ili]; 47478910aadSHong Zhang ba[ili] = uikdi; /* -U(i,k) */ 47578910aadSHong Zhang 47678910aadSHong Zhang /* add multiple of row i to k-th row ... */ 47778910aadSHong Zhang jmin = ili + 1; 47878910aadSHong Zhang nz = bi[i+1] - jmin; 47978910aadSHong Zhang if (nz > 0){ 48078910aadSHong Zhang bcol = bj + jmin; 48178910aadSHong Zhang bval = ba + jmin; 48278910aadSHong Zhang while (nz --) rtmp[*bcol++] += uikdi*(*bval++); 48378910aadSHong Zhang /* update il and jl for i-th row */ 48478910aadSHong Zhang il[i] = jmin; 48578910aadSHong Zhang j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; 48678910aadSHong Zhang } 48778910aadSHong Zhang i = nexti; 48878910aadSHong Zhang } 48978910aadSHong Zhang 4903cea4cbeSHong Zhang /* shift the diagonals when zero pivot is detected */ 4913cea4cbeSHong Zhang /* compute rs=sum of abs(off-diagonal) */ 4923cea4cbeSHong Zhang rs = 0.0; 49378910aadSHong Zhang jmin = bi[k]+1; 49478910aadSHong Zhang nz = bi[k+1] - jmin; 49578910aadSHong Zhang if (nz){ 49678910aadSHong Zhang bcol = bj + jmin; 49778910aadSHong Zhang while (nz--){ 49889f845c8SHong Zhang rs += PetscAbsScalar(rtmp[*bcol]); 49989f845c8SHong Zhang bcol++; 50078910aadSHong Zhang } 50178910aadSHong Zhang } 5023cea4cbeSHong Zhang 5033cea4cbeSHong Zhang sctx.rs = rs; 5043cea4cbeSHong Zhang sctx.pv = dk; 5053e8c821fSHong Zhang ierr = MatCholeskyCheckShift_inline(info,sctx,newshift);CHKERRQ(ierr); 5063cea4cbeSHong Zhang if (newshift == 1){ 5073cea4cbeSHong Zhang break; /* sctx.shift_amount is updated */ 5083cea4cbeSHong Zhang } else if (newshift == -1){ 5093cea4cbeSHong Zhang SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g * rs %g",k,PetscAbsScalar(dk),zeropivot,rs); 51078910aadSHong Zhang } 51178910aadSHong Zhang 51278910aadSHong Zhang /* copy data into U(k,:) */ 51378910aadSHong Zhang ba[bi[k]] = 1.0/dk; 51478910aadSHong Zhang jmin = bi[k]+1; 51578910aadSHong Zhang nz = bi[k+1] - jmin; 51678910aadSHong Zhang if (nz){ 51778910aadSHong Zhang bcol = bj + jmin; 51878910aadSHong Zhang bval = ba + jmin; 51978910aadSHong Zhang while (nz--){ 52078910aadSHong Zhang *bval++ = rtmp[*bcol]; 52178910aadSHong Zhang rtmp[*bcol++] = 0.0; 52278910aadSHong Zhang } 52378910aadSHong Zhang /* add k-th row into il and jl */ 52478910aadSHong Zhang il[k] = jmin; 52578910aadSHong Zhang i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; 52678910aadSHong Zhang } 52778910aadSHong Zhang } 5283cea4cbeSHong Zhang } while (sctx.chshift); 52978910aadSHong Zhang ierr = PetscFree(il);CHKERRQ(ierr); 53078910aadSHong Zhang 53178910aadSHong Zhang C->factor = FACTOR_CHOLESKY; 53278910aadSHong Zhang C->assembled = PETSC_TRUE; 53378910aadSHong Zhang C->preallocated = PETSC_TRUE; 534efee365bSSatish Balay ierr = PetscLogFlops(C->m);CHKERRQ(ierr); 5353cea4cbeSHong Zhang if (sctx.nshift){ 5363cea4cbeSHong Zhang if (shiftnz) { 53763ba0a88SBarry Smith ierr = PetscLogInfo((0,"MatCholeskyFactorNumeric_SeqBAIJ_1_NaturalOrdering: number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount));CHKERRQ(ierr); 5383cea4cbeSHong Zhang } else if (shiftpd) { 53963ba0a88SBarry Smith ierr = PetscLogInfo((0,"MatCholeskyFactorNumeric_SeqBAIJ_1_NaturalOrdering: number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount));CHKERRQ(ierr); 54078910aadSHong Zhang } 54178910aadSHong Zhang } 542c05c3958SHong Zhang PetscFunctionReturn(0); 543c05c3958SHong Zhang } 544c05c3958SHong Zhang 545c05c3958SHong Zhang #include "petscbt.h" 546c05c3958SHong Zhang #include "src/mat/utils/freespace.h" 547c05c3958SHong Zhang #undef __FUNCT__ 548c05c3958SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqBAIJ" 549c05c3958SHong Zhang PetscErrorCode MatICCFactorSymbolic_SeqBAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact) 550c05c3958SHong Zhang { 55178910aadSHong Zhang Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 55278910aadSHong Zhang Mat_SeqSBAIJ *b; 55378910aadSHong Zhang Mat B; 55478910aadSHong Zhang PetscErrorCode ierr; 55578910aadSHong Zhang PetscTruth perm_identity; 55678910aadSHong Zhang PetscInt reallocs=0,*rip,i,*ai=a->i,*aj=a->j,am=a->mbs,bs=A->bs,*ui; 55778910aadSHong Zhang PetscInt jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow; 558cfdb8b8aSHong Zhang PetscInt nlnk,*lnk,*lnk_lvl=PETSC_NULL,ncols,ncols_upper,*cols,*cols_lvl,*uj,**uj_ptr,**uj_lvl_ptr; 55978910aadSHong Zhang PetscReal fill=info->fill,levels=info->levels; 56078910aadSHong Zhang FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 56178910aadSHong Zhang FreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL; 56278910aadSHong Zhang PetscBT lnkbt; 56378910aadSHong Zhang 564c05c3958SHong Zhang PetscFunctionBegin; 5656ad2eaddSHong Zhang if (bs > 1){ 5666ad2eaddSHong Zhang if (!a->sbaijMat){ 567ceb03754SKris Buschelman ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr); 5686ad2eaddSHong Zhang } 5696ad2eaddSHong Zhang ierr = MatICCFactorSymbolic(a->sbaijMat,perm,info,fact);CHKERRQ(ierr); 5706ad2eaddSHong Zhang B = *fact; 5716ad2eaddSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N; 5726ad2eaddSHong Zhang PetscFunctionReturn(0); 5736ad2eaddSHong Zhang } 5746ad2eaddSHong Zhang 57578910aadSHong Zhang ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 57678910aadSHong Zhang ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 57778910aadSHong Zhang 57878910aadSHong Zhang /* special case that simply copies fill pattern */ 57978910aadSHong Zhang if (!levels && perm_identity) { 58078910aadSHong Zhang ierr = MatMarkDiagonal_SeqBAIJ(A);CHKERRQ(ierr); 58178910aadSHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 58278910aadSHong Zhang for (i=0; i<am; i++) { 58378910aadSHong Zhang ui[i] = ai[i+1] - a->diag[i]; /* ui: rowlengths - changes when !perm_identity */ 58478910aadSHong Zhang } 58578910aadSHong Zhang ierr = MatCreate(PETSC_COMM_SELF,am,am,am,am,fact);CHKERRQ(ierr); 58678910aadSHong Zhang B = *fact; 58778910aadSHong Zhang ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr); 58878910aadSHong Zhang ierr = MatSeqSBAIJSetPreallocation(B,1,0,ui);CHKERRQ(ierr); 58978910aadSHong Zhang 59078910aadSHong Zhang b = (Mat_SeqSBAIJ*)B->data; 59178910aadSHong Zhang uj = b->j; 59278910aadSHong Zhang for (i=0; i<am; i++) { 59378910aadSHong Zhang aj = a->j + a->diag[i]; 59478910aadSHong Zhang for (j=0; j<ui[i]; j++){ 59578910aadSHong Zhang *uj++ = *aj++; 59678910aadSHong Zhang } 59778910aadSHong Zhang b->ilen[i] = ui[i]; 59878910aadSHong Zhang } 59978910aadSHong Zhang ierr = PetscFree(ui);CHKERRQ(ierr); 60078910aadSHong Zhang ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 60178910aadSHong Zhang ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 60278910aadSHong Zhang 60378910aadSHong Zhang B->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 60478910aadSHong Zhang B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 60578910aadSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering; 60678910aadSHong Zhang PetscFunctionReturn(0); 60778910aadSHong Zhang } 60878910aadSHong Zhang 60978910aadSHong Zhang /* initialization */ 61078910aadSHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 61178910aadSHong Zhang ui[0] = 0; 61278910aadSHong Zhang ierr = PetscMalloc((2*am+1)*sizeof(PetscInt),&cols_lvl);CHKERRQ(ierr); 61378910aadSHong Zhang 61478910aadSHong Zhang /* jl: linked list for storing indices of the pivot rows 61578910aadSHong Zhang il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */ 61678910aadSHong Zhang ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt*),&jl);CHKERRQ(ierr); 61778910aadSHong Zhang il = jl + am; 61878910aadSHong Zhang uj_ptr = (PetscInt**)(il + am); 61978910aadSHong Zhang uj_lvl_ptr = (PetscInt**)(uj_ptr + am); 62078910aadSHong Zhang for (i=0; i<am; i++){ 62178910aadSHong Zhang jl[i] = am; il[i] = 0; 62278910aadSHong Zhang } 62378910aadSHong Zhang 62478910aadSHong Zhang /* create and initialize a linked list for storing column indices of the active row k */ 62578910aadSHong Zhang nlnk = am + 1; 62678910aadSHong Zhang ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 62778910aadSHong Zhang 62878910aadSHong Zhang /* initial FreeSpace size is fill*(ai[am]+1) */ 62978910aadSHong Zhang ierr = GetMoreSpace((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr); 63078910aadSHong Zhang current_space = free_space; 63178910aadSHong Zhang ierr = GetMoreSpace((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr); 63278910aadSHong Zhang current_space_lvl = free_space_lvl; 63378910aadSHong Zhang 63478910aadSHong Zhang for (k=0; k<am; k++){ /* for each active row k */ 63578910aadSHong Zhang /* initialize lnk by the column indices of row rip[k] of A */ 63678910aadSHong Zhang nzk = 0; 63778910aadSHong Zhang ncols = ai[rip[k]+1] - ai[rip[k]]; 63878910aadSHong Zhang ncols_upper = 0; 63978910aadSHong Zhang cols = cols_lvl + am; 64078910aadSHong Zhang for (j=0; j<ncols; j++){ 64178910aadSHong Zhang i = rip[*(aj + ai[rip[k]] + j)]; 64278910aadSHong Zhang if (i >= k){ /* only take upper triangular entry */ 64378910aadSHong Zhang cols[ncols_upper] = i; 64478910aadSHong Zhang cols_lvl[ncols_upper] = -1; /* initialize level for nonzero entries */ 64578910aadSHong Zhang ncols_upper++; 64678910aadSHong Zhang } 64778910aadSHong Zhang } 64878910aadSHong Zhang ierr = PetscIncompleteLLAdd(ncols_upper,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 64978910aadSHong Zhang nzk += nlnk; 65078910aadSHong Zhang 65178910aadSHong Zhang /* update lnk by computing fill-in for each pivot row to be merged in */ 65278910aadSHong Zhang prow = jl[k]; /* 1st pivot row */ 65378910aadSHong Zhang 65478910aadSHong Zhang while (prow < k){ 65578910aadSHong Zhang nextprow = jl[prow]; 65678910aadSHong Zhang 65778910aadSHong Zhang /* merge prow into k-th row */ 65878910aadSHong Zhang jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:am-1) */ 65978910aadSHong Zhang jmax = ui[prow+1]; 66078910aadSHong Zhang ncols = jmax-jmin; 66178910aadSHong Zhang i = jmin - ui[prow]; 66278910aadSHong Zhang cols = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */ 66378910aadSHong Zhang for (j=0; j<ncols; j++) cols_lvl[j] = *(uj_lvl_ptr[prow] + i + j); 6645a8e39fbSHong Zhang ierr = PetscIncompleteLLAddSorted(ncols,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 66578910aadSHong Zhang nzk += nlnk; 66678910aadSHong Zhang 66778910aadSHong Zhang /* update il and jl for prow */ 66878910aadSHong Zhang if (jmin < jmax){ 66978910aadSHong Zhang il[prow] = jmin; 67078910aadSHong Zhang j = *cols; jl[prow] = jl[j]; jl[j] = prow; 67178910aadSHong Zhang } 67278910aadSHong Zhang prow = nextprow; 67378910aadSHong Zhang } 67478910aadSHong Zhang 67578910aadSHong Zhang /* if free space is not available, make more free space */ 67678910aadSHong Zhang if (current_space->local_remaining<nzk) { 67778910aadSHong Zhang i = am - k + 1; /* num of unfactored rows */ 67878910aadSHong Zhang i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 67978910aadSHong Zhang ierr = GetMoreSpace(i,¤t_space);CHKERRQ(ierr); 68078910aadSHong Zhang ierr = GetMoreSpace(i,¤t_space_lvl);CHKERRQ(ierr); 68178910aadSHong Zhang reallocs++; 68278910aadSHong Zhang } 68378910aadSHong Zhang 68478910aadSHong Zhang /* copy data into free_space and free_space_lvl, then initialize lnk */ 68578910aadSHong Zhang ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr); 68678910aadSHong Zhang 68778910aadSHong Zhang /* add the k-th row into il and jl */ 68878910aadSHong Zhang if (nzk-1 > 0){ 68978910aadSHong Zhang i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */ 69078910aadSHong Zhang jl[k] = jl[i]; jl[i] = k; 69178910aadSHong Zhang il[k] = ui[k] + 1; 69278910aadSHong Zhang } 69378910aadSHong Zhang uj_ptr[k] = current_space->array; 69478910aadSHong Zhang uj_lvl_ptr[k] = current_space_lvl->array; 69578910aadSHong Zhang 69678910aadSHong Zhang current_space->array += nzk; 69778910aadSHong Zhang current_space->local_used += nzk; 69878910aadSHong Zhang current_space->local_remaining -= nzk; 69978910aadSHong Zhang 70078910aadSHong Zhang current_space_lvl->array += nzk; 70178910aadSHong Zhang current_space_lvl->local_used += nzk; 70278910aadSHong Zhang current_space_lvl->local_remaining -= nzk; 70378910aadSHong Zhang 70478910aadSHong Zhang ui[k+1] = ui[k] + nzk; 70578910aadSHong Zhang } 70678910aadSHong Zhang 70763ba0a88SBarry Smith #if defined(PETSC_USE_DEBUG) 70878910aadSHong Zhang if (ai[am] != 0) { 70978910aadSHong Zhang PetscReal af = ((PetscReal)(2*ui[am]-am))/((PetscReal)ai[am]); 71063ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatICCFactorSymbolic_SeqBAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,fill,af));CHKERRQ(ierr); 71163ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatICCFactorSymbolic_SeqBAIJ:Run with -pc_cholesky_fill %g or use \n",af));CHKERRQ(ierr); 71263ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatICCFactorSymbolic_SeqBAIJ:PCCholeskySetFill(pc,%g) for best performance.\n",af));CHKERRQ(ierr); 71378910aadSHong Zhang } else { 71463ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatICCFactorSymbolic_SeqBAIJ:Empty matrix.\n"));CHKERRQ(ierr); 71578910aadSHong Zhang } 71663ba0a88SBarry Smith #endif 71778910aadSHong Zhang 71878910aadSHong Zhang ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 71978910aadSHong Zhang ierr = PetscFree(jl);CHKERRQ(ierr); 72078910aadSHong Zhang ierr = PetscFree(cols_lvl);CHKERRQ(ierr); 72178910aadSHong Zhang 72278910aadSHong Zhang /* destroy list of free space and other temporary array(s) */ 72378910aadSHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 72478910aadSHong Zhang ierr = MakeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 72578910aadSHong Zhang ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 72678910aadSHong Zhang ierr = DestroySpace(free_space_lvl);CHKERRQ(ierr); 72778910aadSHong Zhang 72878910aadSHong Zhang /* put together the new matrix in MATSEQSBAIJ format */ 72978910aadSHong Zhang ierr = MatCreate(PETSC_COMM_SELF,am,am,am,am,fact);CHKERRQ(ierr); 73078910aadSHong Zhang B = *fact; 73178910aadSHong Zhang ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr); 73278910aadSHong Zhang ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr); 73378910aadSHong Zhang 73478910aadSHong Zhang b = (Mat_SeqSBAIJ*)B->data; 73578910aadSHong Zhang ierr = PetscFree(b->imax);CHKERRQ(ierr); 73678910aadSHong Zhang b->singlemalloc = PETSC_FALSE; 73778910aadSHong Zhang /* the next line frees the default space generated by the Create() */ 73878910aadSHong Zhang ierr = PetscFree(b->a);CHKERRQ(ierr); 73978910aadSHong Zhang ierr = PetscFree(b->ilen);CHKERRQ(ierr); 74078910aadSHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 74178910aadSHong Zhang b->j = uj; 74278910aadSHong Zhang b->i = ui; 74378910aadSHong Zhang b->diag = 0; 74478910aadSHong Zhang b->ilen = 0; 74578910aadSHong Zhang b->imax = 0; 74678910aadSHong Zhang b->row = perm; 74778910aadSHong Zhang b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 74878910aadSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 74978910aadSHong Zhang b->icol = perm; 75078910aadSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 75178910aadSHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 75252e6d16bSBarry Smith ierr = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 75378910aadSHong Zhang b->maxnz = b->nz = ui[am]; 75478910aadSHong Zhang 75578910aadSHong Zhang B->factor = FACTOR_CHOLESKY; 75678910aadSHong Zhang B->info.factor_mallocs = reallocs; 75778910aadSHong Zhang B->info.fill_ratio_given = fill; 75878910aadSHong Zhang if (ai[am] != 0) { 75978910aadSHong Zhang B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]); 76078910aadSHong Zhang } else { 76178910aadSHong Zhang B->info.fill_ratio_needed = 0.0; 76278910aadSHong Zhang } 76378910aadSHong Zhang if (perm_identity){ 76478910aadSHong Zhang B->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 76578910aadSHong Zhang B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 76678910aadSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering; 76778910aadSHong Zhang } else { 76878910aadSHong Zhang (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N; 76978910aadSHong Zhang } 770c05c3958SHong Zhang PetscFunctionReturn(0); 771c05c3958SHong Zhang } 772c05c3958SHong Zhang 773c05c3958SHong Zhang #undef __FUNCT__ 774c05c3958SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqBAIJ" 775c05c3958SHong Zhang PetscErrorCode MatCholeskyFactorSymbolic_SeqBAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact) 776c05c3958SHong Zhang { 77778910aadSHong Zhang Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 77878910aadSHong Zhang Mat_SeqSBAIJ *b; 77978910aadSHong Zhang Mat B; 78078910aadSHong Zhang PetscErrorCode ierr; 78178910aadSHong Zhang PetscTruth perm_identity; 78278910aadSHong Zhang PetscReal fill = info->fill; 78378910aadSHong Zhang PetscInt *rip,*riip,i,mbs=a->mbs,bs=A->bs,*ai=a->i,*aj=a->j,reallocs=0,prow; 78478910aadSHong Zhang PetscInt *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow; 78578910aadSHong Zhang PetscInt nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr; 78678910aadSHong Zhang FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 78778910aadSHong Zhang PetscBT lnkbt; 78878910aadSHong Zhang IS iperm; 78978910aadSHong Zhang 790c05c3958SHong Zhang PetscFunctionBegin; 7916ad2eaddSHong Zhang if (bs > 1) { /* convert to seqsbaij */ 7926ad2eaddSHong Zhang if (!a->sbaijMat){ 793ceb03754SKris Buschelman ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr); 7946ad2eaddSHong Zhang } 7956ad2eaddSHong Zhang ierr = MatCholeskyFactorSymbolic(a->sbaijMat,perm,info,fact);CHKERRQ(ierr); 7966ad2eaddSHong Zhang B = *fact; 7976ad2eaddSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N; 7986ad2eaddSHong Zhang PetscFunctionReturn(0); 7996ad2eaddSHong Zhang } 8006ad2eaddSHong Zhang 80178910aadSHong Zhang /* check whether perm is the identity mapping */ 80278910aadSHong Zhang ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 80378910aadSHong Zhang ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 80478910aadSHong Zhang 80578910aadSHong Zhang if (!perm_identity){ 80678910aadSHong Zhang /* check if perm is symmetric! */ 80778910aadSHong Zhang ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr); 80878910aadSHong Zhang ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr); 80978910aadSHong Zhang for (i=0; i<mbs; i++) { 81078910aadSHong Zhang if (rip[i] != riip[i]) SETERRQ(PETSC_ERR_ARG_INCOMP,"Non-symmetric permutation, must use symmetric permutation"); 81178910aadSHong Zhang } 81278910aadSHong Zhang ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr); 81378910aadSHong Zhang ierr = ISDestroy(iperm);CHKERRQ(ierr); 81478910aadSHong Zhang } 81578910aadSHong Zhang 81678910aadSHong Zhang /* initialization */ 81778910aadSHong Zhang ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 81878910aadSHong Zhang ui[0] = 0; 81978910aadSHong Zhang 82078910aadSHong Zhang /* jl: linked list for storing indices of the pivot rows 82178910aadSHong Zhang il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */ 82278910aadSHong Zhang ierr = PetscMalloc((3*mbs+1)*sizeof(PetscInt)+mbs*sizeof(PetscInt*),&jl);CHKERRQ(ierr); 82378910aadSHong Zhang il = jl + mbs; 82478910aadSHong Zhang cols = il + mbs; 82578910aadSHong Zhang ui_ptr = (PetscInt**)(cols + mbs); 82678910aadSHong Zhang for (i=0; i<mbs; i++){ 82778910aadSHong Zhang jl[i] = mbs; il[i] = 0; 82878910aadSHong Zhang } 82978910aadSHong Zhang 83078910aadSHong Zhang /* create and initialize a linked list for storing column indices of the active row k */ 83178910aadSHong Zhang nlnk = mbs + 1; 83278910aadSHong Zhang ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 83378910aadSHong Zhang 83478910aadSHong Zhang /* initial FreeSpace size is fill*(ai[mbs]+1) */ 83578910aadSHong Zhang ierr = GetMoreSpace((PetscInt)(fill*(ai[mbs]+1)),&free_space);CHKERRQ(ierr); 83678910aadSHong Zhang current_space = free_space; 83778910aadSHong Zhang 83878910aadSHong Zhang for (k=0; k<mbs; k++){ /* for each active row k */ 83978910aadSHong Zhang /* initialize lnk by the column indices of row rip[k] of A */ 84078910aadSHong Zhang nzk = 0; 84178910aadSHong Zhang ncols = ai[rip[k]+1] - ai[rip[k]]; 84278910aadSHong Zhang ncols_upper = 0; 84378910aadSHong Zhang for (j=0; j<ncols; j++){ 84478910aadSHong Zhang i = rip[*(aj + ai[rip[k]] + j)]; 84578910aadSHong Zhang if (i >= k){ /* only take upper triangular entry */ 84678910aadSHong Zhang cols[ncols_upper] = i; 84778910aadSHong Zhang ncols_upper++; 84878910aadSHong Zhang } 84978910aadSHong Zhang } 85078910aadSHong Zhang ierr = PetscLLAdd(ncols_upper,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 85178910aadSHong Zhang nzk += nlnk; 85278910aadSHong Zhang 85378910aadSHong Zhang /* update lnk by computing fill-in for each pivot row to be merged in */ 85478910aadSHong Zhang prow = jl[k]; /* 1st pivot row */ 85578910aadSHong Zhang 85678910aadSHong Zhang while (prow < k){ 85778910aadSHong Zhang nextprow = jl[prow]; 85878910aadSHong Zhang /* merge prow into k-th row */ 85978910aadSHong Zhang jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:mbs-1) */ 86078910aadSHong Zhang jmax = ui[prow+1]; 86178910aadSHong Zhang ncols = jmax-jmin; 86278910aadSHong Zhang uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */ 8635a8e39fbSHong Zhang ierr = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 86478910aadSHong Zhang nzk += nlnk; 86578910aadSHong Zhang 86678910aadSHong Zhang /* update il and jl for prow */ 86778910aadSHong Zhang if (jmin < jmax){ 86878910aadSHong Zhang il[prow] = jmin; 86978910aadSHong Zhang j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow; 87078910aadSHong Zhang } 87178910aadSHong Zhang prow = nextprow; 87278910aadSHong Zhang } 87378910aadSHong Zhang 87478910aadSHong Zhang /* if free space is not available, make more free space */ 87578910aadSHong Zhang if (current_space->local_remaining<nzk) { 87678910aadSHong Zhang i = mbs - k + 1; /* num of unfactored rows */ 87778910aadSHong Zhang i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 87878910aadSHong Zhang ierr = GetMoreSpace(i,¤t_space);CHKERRQ(ierr); 87978910aadSHong Zhang reallocs++; 88078910aadSHong Zhang } 88178910aadSHong Zhang 88278910aadSHong Zhang /* copy data into free space, then initialize lnk */ 88378910aadSHong Zhang ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 88478910aadSHong Zhang 88578910aadSHong Zhang /* add the k-th row into il and jl */ 88678910aadSHong Zhang if (nzk-1 > 0){ 88778910aadSHong Zhang i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */ 88878910aadSHong Zhang jl[k] = jl[i]; jl[i] = k; 88978910aadSHong Zhang il[k] = ui[k] + 1; 89078910aadSHong Zhang } 89178910aadSHong Zhang ui_ptr[k] = current_space->array; 89278910aadSHong Zhang current_space->array += nzk; 89378910aadSHong Zhang current_space->local_used += nzk; 89478910aadSHong Zhang current_space->local_remaining -= nzk; 89578910aadSHong Zhang 89678910aadSHong Zhang ui[k+1] = ui[k] + nzk; 89778910aadSHong Zhang } 89878910aadSHong Zhang 89963ba0a88SBarry Smith #if defined(PETSC_USE_DEBUG) 90078910aadSHong Zhang if (ai[mbs] != 0) { 90178910aadSHong Zhang PetscReal af = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]); 90263ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,fill,af));CHKERRQ(ierr); 90363ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Run with -pc_cholesky_fill %g or use \n",af));CHKERRQ(ierr); 90463ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqSBAIJ:PCCholeskySetFill(pc,%g) for best performance.\n",af));CHKERRQ(ierr); 90578910aadSHong Zhang } else { 90663ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Empty matrix.\n"));CHKERRQ(ierr); 90778910aadSHong Zhang } 90863ba0a88SBarry Smith #endif 90978910aadSHong Zhang 91078910aadSHong Zhang ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 91178910aadSHong Zhang ierr = PetscFree(jl);CHKERRQ(ierr); 91278910aadSHong Zhang 91378910aadSHong Zhang /* destroy list of free space and other temporary array(s) */ 91478910aadSHong Zhang ierr = PetscMalloc((ui[mbs]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 91578910aadSHong Zhang ierr = MakeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 91678910aadSHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 91778910aadSHong Zhang 91878910aadSHong Zhang /* put together the new matrix in MATSEQSBAIJ format */ 91978910aadSHong Zhang ierr = MatCreate(PETSC_COMM_SELF,mbs,mbs,mbs,mbs,fact);CHKERRQ(ierr); 92078910aadSHong Zhang B = *fact; 92178910aadSHong Zhang ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr); 92278910aadSHong Zhang ierr = MatSeqSBAIJSetPreallocation(B,bs,0,PETSC_NULL);CHKERRQ(ierr); 92378910aadSHong Zhang 92478910aadSHong Zhang b = (Mat_SeqSBAIJ*)B->data; 92578910aadSHong Zhang ierr = PetscFree(b->imax);CHKERRQ(ierr); 92678910aadSHong Zhang b->singlemalloc = PETSC_FALSE; 92778910aadSHong Zhang /* the next line frees the default space generated by the Create() */ 92878910aadSHong Zhang ierr = PetscFree(b->a);CHKERRQ(ierr); 92978910aadSHong Zhang ierr = PetscFree(b->ilen);CHKERRQ(ierr); 93078910aadSHong Zhang ierr = PetscMalloc((ui[mbs]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 93178910aadSHong Zhang b->j = uj; 93278910aadSHong Zhang b->i = ui; 93378910aadSHong Zhang b->diag = 0; 93478910aadSHong Zhang b->ilen = 0; 93578910aadSHong Zhang b->imax = 0; 93678910aadSHong Zhang b->row = perm; 93778910aadSHong Zhang b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 93878910aadSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 93978910aadSHong Zhang b->icol = perm; 94078910aadSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 94178910aadSHong Zhang ierr = PetscMalloc((mbs+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 94252e6d16bSBarry Smith ierr = PetscLogObjectMemory(B,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 94378910aadSHong Zhang b->maxnz = b->nz = ui[mbs]; 94478910aadSHong Zhang 94578910aadSHong Zhang B->factor = FACTOR_CHOLESKY; 94678910aadSHong Zhang B->info.factor_mallocs = reallocs; 94778910aadSHong Zhang B->info.fill_ratio_given = fill; 94878910aadSHong Zhang if (ai[mbs] != 0) { 94978910aadSHong Zhang B->info.fill_ratio_needed = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]); 95078910aadSHong Zhang } else { 95178910aadSHong Zhang B->info.fill_ratio_needed = 0.0; 95278910aadSHong Zhang } 95378910aadSHong Zhang if (perm_identity){ 9546ad2eaddSHong Zhang B->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 9556ad2eaddSHong Zhang B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 9566ad2eaddSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering; 95778910aadSHong Zhang } else { 9586ad2eaddSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N; 95978910aadSHong Zhang } 960c05c3958SHong Zhang PetscFunctionReturn(0); 961c05c3958SHong Zhang } 962