1be1d678aSKris Buschelman #define PETSCMAT_DLL 2be1d678aSKris 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" 15*719d5645SBarry Smith PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2(Mat B,Mat A,MatFactorInfo *info) 164eeb42bcSBarry Smith { 17*719d5645SBarry 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; 2762bba022SBarry Smith PetscReal shift = info->shiftinblocks; 284eeb42bcSBarry Smith 293a40ed3dSBarry Smith PetscFunctionBegin; 304eeb42bcSBarry Smith ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 314eeb42bcSBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 32b0a32e0cSBarry Smith ierr = PetscMalloc(4*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 334eeb42bcSBarry Smith 344eeb42bcSBarry Smith for (i=0; i<n; i++) { 354078e994SBarry Smith nz = bi[i+1] - bi[i]; 364078e994SBarry Smith ajtmp = bj + bi[i]; 374eeb42bcSBarry Smith for (j=0; j<nz; j++) { 384eeb42bcSBarry Smith x = rtmp+4*ajtmp[j]; x[0] = x[1] = x[2] = x[3] = 0.0; 394eeb42bcSBarry Smith } 404eeb42bcSBarry Smith /* load in initial (unfactored row) */ 414eeb42bcSBarry Smith idx = r[i]; 424078e994SBarry Smith nz = ai[idx+1] - ai[idx]; 434078e994SBarry Smith ajtmpold = aj + ai[idx]; 444078e994SBarry Smith v = aa + 4*ai[idx]; 454eeb42bcSBarry Smith for (j=0; j<nz; j++) { 464eeb42bcSBarry Smith x = rtmp+4*ic[ajtmpold[j]]; 474eeb42bcSBarry Smith x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3]; 484eeb42bcSBarry Smith v += 4; 494eeb42bcSBarry Smith } 504eeb42bcSBarry Smith row = *ajtmp++; 514eeb42bcSBarry Smith while (row < i) { 524eeb42bcSBarry Smith pc = rtmp + 4*row; 534eeb42bcSBarry Smith p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3]; 5488685aaeSLois Curfman McInnes if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0) { 554078e994SBarry Smith pv = ba + 4*diag_offset[row]; 564078e994SBarry Smith pj = bj + diag_offset[row] + 1; 574eeb42bcSBarry Smith x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 584eeb42bcSBarry Smith pc[0] = m1 = p1*x1 + p3*x2; 594eeb42bcSBarry Smith pc[1] = m2 = p2*x1 + p4*x2; 604eeb42bcSBarry Smith pc[2] = m3 = p1*x3 + p3*x4; 614eeb42bcSBarry Smith pc[3] = m4 = p2*x3 + p4*x4; 624078e994SBarry Smith nz = bi[row+1] - diag_offset[row] - 1; 634eeb42bcSBarry Smith pv += 4; 644eeb42bcSBarry Smith for (j=0; j<nz; j++) { 654eeb42bcSBarry Smith x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 664eeb42bcSBarry Smith x = rtmp + 4*pj[j]; 674eeb42bcSBarry Smith x[0] -= m1*x1 + m3*x2; 684eeb42bcSBarry Smith x[1] -= m2*x1 + m4*x2; 694eeb42bcSBarry Smith x[2] -= m1*x3 + m3*x4; 704eeb42bcSBarry Smith x[3] -= m2*x3 + m4*x4; 714eeb42bcSBarry Smith pv += 4; 724eeb42bcSBarry Smith } 73efee365bSSatish Balay ierr = PetscLogFlops(16*nz+12);CHKERRQ(ierr); 744eeb42bcSBarry Smith } 754eeb42bcSBarry Smith row = *ajtmp++; 764eeb42bcSBarry Smith } 774eeb42bcSBarry Smith /* finished row so stick it into b->a */ 784078e994SBarry Smith pv = ba + 4*bi[i]; 794078e994SBarry Smith pj = bj + bi[i]; 804078e994SBarry Smith nz = bi[i+1] - bi[i]; 814eeb42bcSBarry Smith for (j=0; j<nz; j++) { 824eeb42bcSBarry Smith x = rtmp+4*pj[j]; 834eeb42bcSBarry Smith pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3]; 844eeb42bcSBarry Smith pv += 4; 854eeb42bcSBarry Smith } 864eeb42bcSBarry Smith /* invert diagonal block */ 874078e994SBarry Smith w = ba + 4*diag_offset[i]; 8862bba022SBarry Smith ierr = Kernel_A_gets_inverse_A_2(w,shift);CHKERRQ(ierr); 894eeb42bcSBarry Smith } 904eeb42bcSBarry Smith 91606d414cSSatish Balay ierr = PetscFree(rtmp);CHKERRQ(ierr); 924eeb42bcSBarry Smith ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 934eeb42bcSBarry Smith ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 94db4efbfdSBarry Smith C->ops->solve = MatSolve_SeqBAIJ_2; 95db4efbfdSBarry Smith C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2; 964eeb42bcSBarry Smith C->assembled = PETSC_TRUE; 97efee365bSSatish Balay ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 983a40ed3dSBarry Smith PetscFunctionReturn(0); 994eeb42bcSBarry Smith } 1004cd38560SBarry Smith /* 1014cd38560SBarry Smith Version for when blocks are 2 by 2 Using natural ordering 1024cd38560SBarry Smith */ 1034a2ae208SSatish Balay #undef __FUNCT__ 1044a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering" 105*719d5645SBarry Smith PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering(Mat C,Mat A,MatFactorInfo *info) 1064cd38560SBarry Smith { 1074cd38560SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; 108dfbe8321SBarry Smith PetscErrorCode ierr; 109690b6cddSBarry Smith PetscInt i,j,n = a->mbs,*bi = b->i,*bj = b->j; 110690b6cddSBarry Smith PetscInt *ajtmpold,*ajtmp,nz,row; 111690b6cddSBarry Smith PetscInt *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj; 112329f5518SBarry Smith MatScalar *pv,*v,*rtmp,*pc,*w,*x; 1132515f8d2SSatish Balay MatScalar p1,p2,p3,p4,m1,m2,m3,m4,x1,x2,x3,x4; 1144cd38560SBarry Smith MatScalar *ba = b->a,*aa = a->a; 11562bba022SBarry Smith PetscReal shift = info->shiftinblocks; 1164cd38560SBarry Smith 1174cd38560SBarry Smith PetscFunctionBegin; 118b0a32e0cSBarry Smith ierr = PetscMalloc(4*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 1194cd38560SBarry Smith 1204cd38560SBarry Smith for (i=0; i<n; i++) { 1214cd38560SBarry Smith nz = bi[i+1] - bi[i]; 1224cd38560SBarry Smith ajtmp = bj + bi[i]; 1234cd38560SBarry Smith for (j=0; j<nz; j++) { 1244cd38560SBarry Smith x = rtmp+4*ajtmp[j]; 1254cd38560SBarry Smith x[0] = x[1] = x[2] = x[3] = 0.0; 1264cd38560SBarry Smith } 1274cd38560SBarry Smith /* load in initial (unfactored row) */ 1284cd38560SBarry Smith nz = ai[i+1] - ai[i]; 1294cd38560SBarry Smith ajtmpold = aj + ai[i]; 1304cd38560SBarry Smith v = aa + 4*ai[i]; 1314cd38560SBarry Smith for (j=0; j<nz; j++) { 1324cd38560SBarry Smith x = rtmp+4*ajtmpold[j]; 1334cd38560SBarry Smith x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3]; 1344cd38560SBarry Smith v += 4; 1354cd38560SBarry Smith } 1364cd38560SBarry Smith row = *ajtmp++; 1374cd38560SBarry Smith while (row < i) { 1384cd38560SBarry Smith pc = rtmp + 4*row; 1394cd38560SBarry Smith p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3]; 1404cd38560SBarry Smith if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0) { 1414cd38560SBarry Smith pv = ba + 4*diag_offset[row]; 1424cd38560SBarry Smith pj = bj + diag_offset[row] + 1; 1434cd38560SBarry Smith x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 1444cd38560SBarry Smith pc[0] = m1 = p1*x1 + p3*x2; 1454cd38560SBarry Smith pc[1] = m2 = p2*x1 + p4*x2; 1464cd38560SBarry Smith pc[2] = m3 = p1*x3 + p3*x4; 1474cd38560SBarry Smith pc[3] = m4 = p2*x3 + p4*x4; 1484cd38560SBarry Smith nz = bi[row+1] - diag_offset[row] - 1; 1494cd38560SBarry Smith pv += 4; 1504cd38560SBarry Smith for (j=0; j<nz; j++) { 1514cd38560SBarry Smith x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 1524cd38560SBarry Smith x = rtmp + 4*pj[j]; 1534cd38560SBarry Smith x[0] -= m1*x1 + m3*x2; 1544cd38560SBarry Smith x[1] -= m2*x1 + m4*x2; 1554cd38560SBarry Smith x[2] -= m1*x3 + m3*x4; 1564cd38560SBarry Smith x[3] -= m2*x3 + m4*x4; 1574cd38560SBarry Smith pv += 4; 1584cd38560SBarry Smith } 159efee365bSSatish Balay ierr = PetscLogFlops(16*nz+12);CHKERRQ(ierr); 1604cd38560SBarry Smith } 1614cd38560SBarry Smith row = *ajtmp++; 1624cd38560SBarry Smith } 1634cd38560SBarry Smith /* finished row so stick it into b->a */ 1644cd38560SBarry Smith pv = ba + 4*bi[i]; 1654cd38560SBarry Smith pj = bj + bi[i]; 1664cd38560SBarry Smith nz = bi[i+1] - bi[i]; 1674cd38560SBarry Smith for (j=0; j<nz; j++) { 1684cd38560SBarry Smith x = rtmp+4*pj[j]; 1694cd38560SBarry Smith pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3]; 1704cd38560SBarry Smith pv += 4; 1714cd38560SBarry Smith } 1724cd38560SBarry Smith /* invert diagonal block */ 1734cd38560SBarry Smith w = ba + 4*diag_offset[i]; 17462bba022SBarry Smith ierr = Kernel_A_gets_inverse_A_2(w,shift);CHKERRQ(ierr); 1754cd38560SBarry Smith } 1764cd38560SBarry Smith 177606d414cSSatish Balay ierr = PetscFree(rtmp);CHKERRQ(ierr); 178db4efbfdSBarry Smith C->ops->solve = MatSolve_SeqBAIJ_2_NaturalOrdering; 179db4efbfdSBarry Smith C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_NaturalOrdering; 1804cd38560SBarry Smith C->assembled = PETSC_TRUE; 181efee365bSSatish Balay ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 1824cd38560SBarry Smith PetscFunctionReturn(0); 1834cd38560SBarry Smith } 1847fc0212eSBarry Smith 1857fc0212eSBarry Smith /* ----------------------------------------------------------- */ 1867fc0212eSBarry Smith /* 1877fc0212eSBarry Smith Version for when blocks are 1 by 1. 1887fc0212eSBarry Smith */ 1894a2ae208SSatish Balay #undef __FUNCT__ 1904a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_1" 191*719d5645SBarry Smith PetscErrorCode MatLUFactorNumeric_SeqBAIJ_1(Mat C,Mat A,MatFactorInfo *info) 1927fc0212eSBarry Smith { 1937fc0212eSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; 1947cdf2dbbSSatish Balay IS isrow = b->row,isicol = b->icol; 1956849ba73SBarry Smith PetscErrorCode ierr; 196690b6cddSBarry Smith PetscInt *r,*ic,i,j,n = a->mbs,*bi = b->i,*bj = b->j; 197690b6cddSBarry Smith PetscInt *ajtmpold,*ajtmp,nz,row,*ai = a->i,*aj = a->j; 198690b6cddSBarry Smith PetscInt *diag_offset = b->diag,diag,*pj; 199329f5518SBarry Smith MatScalar *pv,*v,*rtmp,multiplier,*pc; 2003f1db9ecSBarry Smith MatScalar *ba = b->a,*aa = a->a; 2017fc0212eSBarry Smith 2023a40ed3dSBarry Smith PetscFunctionBegin; 2037fc0212eSBarry Smith ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 2047fc0212eSBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 205b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 2067fc0212eSBarry Smith 2077fc0212eSBarry Smith for (i=0; i<n; i++) { 2084078e994SBarry Smith nz = bi[i+1] - bi[i]; 2094078e994SBarry Smith ajtmp = bj + bi[i]; 2107fc0212eSBarry Smith for (j=0; j<nz; j++) rtmp[ajtmp[j]] = 0.0; 2117fc0212eSBarry Smith 2127fc0212eSBarry Smith /* load in initial (unfactored row) */ 2134078e994SBarry Smith nz = ai[r[i]+1] - ai[r[i]]; 2144078e994SBarry Smith ajtmpold = aj + ai[r[i]]; 2154078e994SBarry Smith v = aa + ai[r[i]]; 2167fc0212eSBarry Smith for (j=0; j<nz; j++) rtmp[ic[ajtmpold[j]]] = v[j]; 2177fc0212eSBarry Smith 2187fc0212eSBarry Smith row = *ajtmp++; 2197fc0212eSBarry Smith while (row < i) { 2207fc0212eSBarry Smith pc = rtmp + row; 2217fc0212eSBarry Smith if (*pc != 0.0) { 2224078e994SBarry Smith pv = ba + diag_offset[row]; 2234078e994SBarry Smith pj = bj + diag_offset[row] + 1; 2247fc0212eSBarry Smith multiplier = *pc * *pv++; 2257fc0212eSBarry Smith *pc = multiplier; 2264078e994SBarry Smith nz = bi[row+1] - diag_offset[row] - 1; 2277fc0212eSBarry Smith for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j]; 228efee365bSSatish Balay ierr = PetscLogFlops(1+2*nz);CHKERRQ(ierr); 2297fc0212eSBarry Smith } 2307fc0212eSBarry Smith row = *ajtmp++; 2317fc0212eSBarry Smith } 2327fc0212eSBarry Smith /* finished row so stick it into b->a */ 2334078e994SBarry Smith pv = ba + bi[i]; 2344078e994SBarry Smith pj = bj + bi[i]; 2354078e994SBarry Smith nz = bi[i+1] - bi[i]; 2367fc0212eSBarry Smith for (j=0; j<nz; j++) {pv[j] = rtmp[pj[j]];} 2374078e994SBarry Smith diag = diag_offset[i] - bi[i]; 2387fc0212eSBarry Smith /* check pivot entry for current row */ 239a73cf429SBarry Smith if (pv[diag] == 0.0) { 24029bbc08cSBarry Smith SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot"); 2417fc0212eSBarry Smith } 2427fc0212eSBarry Smith pv[diag] = 1.0/pv[diag]; 2437fc0212eSBarry Smith } 2447fc0212eSBarry Smith 245606d414cSSatish Balay ierr = PetscFree(rtmp);CHKERRQ(ierr); 2467fc0212eSBarry Smith ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 2477fc0212eSBarry Smith ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 248db4efbfdSBarry Smith C->ops->solve = MatSolve_SeqBAIJ_1; 249db4efbfdSBarry Smith C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1; 2507fc0212eSBarry Smith C->assembled = PETSC_TRUE; 251d0f46423SBarry Smith ierr = PetscLogFlops(C->cmap->n);CHKERRQ(ierr); 2523a40ed3dSBarry Smith PetscFunctionReturn(0); 2537fc0212eSBarry Smith } 2547fc0212eSBarry Smith 255e631078cSBarry Smith EXTERN_C_BEGIN 256b24902e0SBarry Smith #undef __FUNCT__ 257b24902e0SBarry Smith #define __FUNCT__ "MatGetFactor_seqbaij_petsc" 258b24902e0SBarry Smith PetscErrorCode MatGetFactor_seqbaij_petsc(Mat A,MatFactorType ftype,Mat *B) 259b24902e0SBarry Smith { 260d0f46423SBarry Smith PetscInt n = A->rmap->n; 261b24902e0SBarry Smith PetscErrorCode ierr; 262b24902e0SBarry Smith 263b24902e0SBarry Smith PetscFunctionBegin; 264b24902e0SBarry Smith ierr = MatCreate(((PetscObject)A)->comm,B);CHKERRQ(ierr); 265b24902e0SBarry Smith ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr); 266b24902e0SBarry Smith if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) { 267b24902e0SBarry Smith ierr = MatSetType(*B,MATSEQBAIJ);CHKERRQ(ierr); 268db4efbfdSBarry Smith (*B)->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqBAIJ; 269db4efbfdSBarry Smith (*B)->ops->ilufactorsymbolic = MatILUFactorSymbolic_SeqBAIJ; 270b24902e0SBarry Smith } else if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) { 271b24902e0SBarry Smith ierr = MatSetType(*B,MATSEQSBAIJ);CHKERRQ(ierr); 2725c9eb25fSBarry Smith ierr = MatSeqSBAIJSetPreallocation(*B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 2735c9eb25fSBarry Smith (*B)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqBAIJ; 2745c9eb25fSBarry Smith (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqBAIJ; 275b24902e0SBarry Smith } else SETERRQ(PETSC_ERR_SUP,"Factor type not supported"); 276*719d5645SBarry Smith (*B)->factor = ftype; 277b24902e0SBarry Smith PetscFunctionReturn(0); 278b24902e0SBarry Smith } 279e631078cSBarry Smith EXTERN_C_END 280273d9f13SBarry Smith 281db4efbfdSBarry Smith EXTERN_C_BEGIN 282db4efbfdSBarry Smith #undef __FUNCT__ 283db4efbfdSBarry Smith #define __FUNCT__ "MatGetFactorAvailable_seqaij_petsc" 284db4efbfdSBarry Smith PetscErrorCode MatGetFactorAvailable_seqbaij_petsc(Mat A,MatFactorType ftype,PetscTruth *flg) 285db4efbfdSBarry Smith { 286db4efbfdSBarry Smith PetscFunctionBegin; 287db4efbfdSBarry Smith *flg = PETSC_TRUE; 288db4efbfdSBarry Smith PetscFunctionReturn(0); 289db4efbfdSBarry Smith } 290db4efbfdSBarry Smith EXTERN_C_END 291db4efbfdSBarry Smith 2925d517e7eSBarry Smith /* ----------------------------------------------------------- */ 2934a2ae208SSatish Balay #undef __FUNCT__ 2944a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqBAIJ" 295dfbe8321SBarry Smith PetscErrorCode MatLUFactor_SeqBAIJ(Mat A,IS row,IS col,MatFactorInfo *info) 2965d517e7eSBarry Smith { 297dfbe8321SBarry Smith PetscErrorCode ierr; 2985d517e7eSBarry Smith Mat C; 2995d517e7eSBarry Smith 3003a40ed3dSBarry Smith PetscFunctionBegin; 30143244d56SBarry Smith ierr = MatGetFactor(A,MAT_SOLVER_PETSC,MAT_FACTOR_LU,&C);CHKERRQ(ierr); 302*719d5645SBarry Smith ierr = MatLUFactorSymbolic(C,A,row,col,info);CHKERRQ(ierr); 303*719d5645SBarry Smith ierr = MatLUFactorNumeric(C,A,info);CHKERRQ(ierr); 304db4efbfdSBarry Smith A->ops->solve = C->ops->solve; 305db4efbfdSBarry Smith A->ops->solvetranspose = C->ops->solvetranspose; 306273d9f13SBarry Smith ierr = MatHeaderCopy(A,C);CHKERRQ(ierr); 30752e6d16bSBarry Smith ierr = PetscLogObjectParent(A,((Mat_SeqBAIJ*)(A->data))->icol);CHKERRQ(ierr); 3083a40ed3dSBarry Smith PetscFunctionReturn(0); 3095d517e7eSBarry Smith } 3104cd38560SBarry Smith 311c05c3958SHong Zhang #include "src/mat/impls/sbaij/seq/sbaij.h" 312c05c3958SHong Zhang #undef __FUNCT__ 313c05c3958SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N" 314*719d5645SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N(Mat C,Mat A,MatFactorInfo *info) 315c05c3958SHong Zhang { 316f3086b4bSHong Zhang PetscErrorCode ierr; 31778910aadSHong Zhang Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data; 31878910aadSHong Zhang Mat_SeqSBAIJ *b=(Mat_SeqSBAIJ*)C->data; 31978910aadSHong Zhang IS ip=b->row; 320d0f46423SBarry Smith PetscInt *rip,i,j,mbs=a->mbs,bs=A->rmap->bs,*bi=b->i,*bj=b->j,*bcol; 32178910aadSHong Zhang PetscInt *ai=a->i,*aj=a->j; 32278910aadSHong Zhang PetscInt k,jmin,jmax,*jl,*il,col,nexti,ili,nz; 32378910aadSHong Zhang MatScalar *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi; 3243cea4cbeSHong Zhang PetscReal zeropivot,rs,shiftnz; 325fbf22428SSatish Balay PetscReal shiftpd; 3263cea4cbeSHong Zhang ChShift_Ctx sctx; 3273cea4cbeSHong Zhang PetscInt newshift; 32878910aadSHong Zhang 329c05c3958SHong Zhang PetscFunctionBegin; 3306ad2eaddSHong Zhang if (bs > 1) { 3316ad2eaddSHong Zhang if (!a->sbaijMat){ 332ceb03754SKris Buschelman ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr); 3336ad2eaddSHong Zhang } 334*719d5645SBarry Smith ierr = (a->sbaijMat)->ops->choleskyfactornumeric(C,a->sbaijMat,info);CHKERRQ(ierr); 3356ad2eaddSHong Zhang ierr = MatDestroy(a->sbaijMat);CHKERRQ(ierr); 3366ad2eaddSHong Zhang a->sbaijMat = PETSC_NULL; 3376ad2eaddSHong Zhang PetscFunctionReturn(0); 3386ad2eaddSHong Zhang } 33978910aadSHong Zhang 34078910aadSHong Zhang /* initialization */ 3413cea4cbeSHong Zhang shiftnz = info->shiftnz; 3423cea4cbeSHong Zhang shiftpd = info->shiftpd; 3433cea4cbeSHong Zhang zeropivot = info->zeropivot; 3443cea4cbeSHong Zhang 3456ad2eaddSHong Zhang ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr); 34678910aadSHong Zhang nz = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar); 34778910aadSHong Zhang ierr = PetscMalloc(nz,&il);CHKERRQ(ierr); 34878910aadSHong Zhang jl = il + mbs; 34978910aadSHong Zhang rtmp = (MatScalar*)(jl + mbs); 35078910aadSHong Zhang 3513cea4cbeSHong Zhang sctx.shift_amount = 0; 3523cea4cbeSHong Zhang sctx.nshift = 0; 35378910aadSHong Zhang do { 3543cea4cbeSHong Zhang sctx.chshift = PETSC_FALSE; 35578910aadSHong Zhang for (i=0; i<mbs; i++) { 35678910aadSHong Zhang rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0; 35778910aadSHong Zhang } 35878910aadSHong Zhang 35978910aadSHong Zhang for (k = 0; k<mbs; k++){ 36078910aadSHong Zhang bval = ba + bi[k]; 36178910aadSHong Zhang /* initialize k-th row by the perm[k]-th row of A */ 36278910aadSHong Zhang jmin = ai[rip[k]]; jmax = ai[rip[k]+1]; 36378910aadSHong Zhang for (j = jmin; j < jmax; j++){ 36478910aadSHong Zhang col = rip[aj[j]]; 36578910aadSHong Zhang if (col >= k){ /* only take upper triangular entry */ 36678910aadSHong Zhang rtmp[col] = aa[j]; 36778910aadSHong Zhang *bval++ = 0.0; /* for in-place factorization */ 36878910aadSHong Zhang } 36978910aadSHong Zhang } 3703cea4cbeSHong Zhang 3713cea4cbeSHong Zhang /* shift the diagonal of the matrix */ 3723cea4cbeSHong Zhang if (sctx.nshift) rtmp[k] += sctx.shift_amount; 37378910aadSHong Zhang 37478910aadSHong Zhang /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 37578910aadSHong Zhang dk = rtmp[k]; 37678910aadSHong Zhang i = jl[k]; /* first row to be added to k_th row */ 37778910aadSHong Zhang 37878910aadSHong Zhang while (i < k){ 37978910aadSHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 38078910aadSHong Zhang 38178910aadSHong Zhang /* compute multiplier, update diag(k) and U(i,k) */ 38278910aadSHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 38378910aadSHong Zhang uikdi = - ba[ili]*ba[bi[i]]; /* diagonal(k) */ 38478910aadSHong Zhang dk += uikdi*ba[ili]; 38578910aadSHong Zhang ba[ili] = uikdi; /* -U(i,k) */ 38678910aadSHong Zhang 38778910aadSHong Zhang /* add multiple of row i to k-th row */ 38878910aadSHong Zhang jmin = ili + 1; jmax = bi[i+1]; 38978910aadSHong Zhang if (jmin < jmax){ 39078910aadSHong Zhang for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j]; 39178910aadSHong Zhang /* update il and jl for row i */ 39278910aadSHong Zhang il[i] = jmin; 39378910aadSHong Zhang j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; 39478910aadSHong Zhang } 39578910aadSHong Zhang i = nexti; 39678910aadSHong Zhang } 39778910aadSHong Zhang 3983cea4cbeSHong Zhang /* shift the diagonals when zero pivot is detected */ 3993cea4cbeSHong Zhang /* compute rs=sum of abs(off-diagonal) */ 4003cea4cbeSHong Zhang rs = 0.0; 40178910aadSHong Zhang jmin = bi[k]+1; 40278910aadSHong Zhang nz = bi[k+1] - jmin; 40378910aadSHong Zhang if (nz){ 40478910aadSHong Zhang bcol = bj + jmin; 40578910aadSHong Zhang while (nz--){ 40689f845c8SHong Zhang rs += PetscAbsScalar(rtmp[*bcol]); 40789f845c8SHong Zhang bcol++; 40878910aadSHong Zhang } 40978910aadSHong Zhang } 4103cea4cbeSHong Zhang 4113cea4cbeSHong Zhang sctx.rs = rs; 4123cea4cbeSHong Zhang sctx.pv = dk; 41345938b79SHong Zhang ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr); 41445938b79SHong Zhang if (newshift == 1) break; 41578910aadSHong Zhang 41678910aadSHong Zhang /* copy data into U(k,:) */ 41778910aadSHong Zhang ba[bi[k]] = 1.0/dk; /* U(k,k) */ 41878910aadSHong Zhang jmin = bi[k]+1; jmax = bi[k+1]; 41978910aadSHong Zhang if (jmin < jmax) { 42078910aadSHong Zhang for (j=jmin; j<jmax; j++){ 42178910aadSHong Zhang col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0; 42278910aadSHong Zhang } 42378910aadSHong Zhang /* add the k-th row into il and jl */ 42478910aadSHong Zhang il[k] = jmin; 42578910aadSHong Zhang i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; 42678910aadSHong Zhang } 42778910aadSHong Zhang } 4283cea4cbeSHong Zhang } while (sctx.chshift); 42978910aadSHong Zhang ierr = PetscFree(il);CHKERRQ(ierr); 43078910aadSHong Zhang 43178910aadSHong Zhang ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr); 43278910aadSHong Zhang C->assembled = PETSC_TRUE; 43378910aadSHong Zhang C->preallocated = PETSC_TRUE; 434d0f46423SBarry Smith ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr); 4353cea4cbeSHong Zhang if (sctx.nshift){ 4363cea4cbeSHong Zhang if (shiftnz) { 4371e2582c4SBarry Smith ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 4383cea4cbeSHong Zhang } else if (shiftpd) { 4391e2582c4SBarry Smith ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 44078910aadSHong Zhang } 44178910aadSHong Zhang } 442c05c3958SHong Zhang PetscFunctionReturn(0); 443c05c3958SHong Zhang } 4444cd38560SBarry Smith 445c05c3958SHong Zhang #undef __FUNCT__ 446c05c3958SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering" 447*719d5645SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering(Mat C,Mat A,MatFactorInfo *info) 448c05c3958SHong Zhang { 44978910aadSHong Zhang Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data; 45078910aadSHong Zhang Mat_SeqSBAIJ *b=(Mat_SeqSBAIJ*)C->data; 45178910aadSHong Zhang PetscErrorCode ierr; 4523cea4cbeSHong Zhang PetscInt i,j,am=a->mbs; 45378910aadSHong Zhang PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 4543cea4cbeSHong Zhang PetscInt k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz; 45578910aadSHong Zhang MatScalar *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval; 4563cea4cbeSHong Zhang PetscReal zeropivot,rs,shiftnz; 457fbf22428SSatish Balay PetscReal shiftpd; 4583cea4cbeSHong Zhang ChShift_Ctx sctx; 4593cea4cbeSHong Zhang PetscInt newshift; 46078910aadSHong Zhang 461c05c3958SHong Zhang PetscFunctionBegin; 46278910aadSHong Zhang /* initialization */ 4633cea4cbeSHong Zhang shiftnz = info->shiftnz; 4643cea4cbeSHong Zhang shiftpd = info->shiftpd; 4653cea4cbeSHong Zhang zeropivot = info->zeropivot; 4663cea4cbeSHong Zhang 46778910aadSHong Zhang nz = (2*am+1)*sizeof(PetscInt)+am*sizeof(MatScalar); 46878910aadSHong Zhang ierr = PetscMalloc(nz,&il);CHKERRQ(ierr); 46978910aadSHong Zhang jl = il + am; 47078910aadSHong Zhang rtmp = (MatScalar*)(jl + am); 47178910aadSHong Zhang 4723cea4cbeSHong Zhang sctx.shift_amount = 0; 4733cea4cbeSHong Zhang sctx.nshift = 0; 47478910aadSHong Zhang do { 4753cea4cbeSHong Zhang sctx.chshift = PETSC_FALSE; 47678910aadSHong Zhang for (i=0; i<am; i++) { 47778910aadSHong Zhang rtmp[i] = 0.0; jl[i] = am; il[0] = 0; 47878910aadSHong Zhang } 47978910aadSHong Zhang 48078910aadSHong Zhang for (k = 0; k<am; k++){ 48178910aadSHong Zhang /* initialize k-th row with elements nonzero in row perm(k) of A */ 48278910aadSHong Zhang nz = ai[k+1] - ai[k]; 48378910aadSHong Zhang acol = aj + ai[k]; 48478910aadSHong Zhang aval = aa + ai[k]; 48578910aadSHong Zhang bval = ba + bi[k]; 48678910aadSHong Zhang while (nz -- ){ 48778910aadSHong Zhang if (*acol < k) { /* skip lower triangular entries */ 48878910aadSHong Zhang acol++; aval++; 48978910aadSHong Zhang } else { 49078910aadSHong Zhang rtmp[*acol++] = *aval++; 49178910aadSHong Zhang *bval++ = 0.0; /* for in-place factorization */ 49278910aadSHong Zhang } 49378910aadSHong Zhang } 4943cea4cbeSHong Zhang 4953cea4cbeSHong Zhang /* shift the diagonal of the matrix */ 4963cea4cbeSHong Zhang if (sctx.nshift) rtmp[k] += sctx.shift_amount; 49778910aadSHong Zhang 49878910aadSHong Zhang /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 49978910aadSHong Zhang dk = rtmp[k]; 50078910aadSHong Zhang i = jl[k]; /* first row to be added to k_th row */ 50178910aadSHong Zhang 50278910aadSHong Zhang while (i < k){ 50378910aadSHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 50478910aadSHong Zhang /* compute multiplier, update D(k) and U(i,k) */ 50578910aadSHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 50678910aadSHong Zhang uikdi = - ba[ili]*ba[bi[i]]; 50778910aadSHong Zhang dk += uikdi*ba[ili]; 50878910aadSHong Zhang ba[ili] = uikdi; /* -U(i,k) */ 50978910aadSHong Zhang 51078910aadSHong Zhang /* add multiple of row i to k-th row ... */ 51178910aadSHong Zhang jmin = ili + 1; 51278910aadSHong Zhang nz = bi[i+1] - jmin; 51378910aadSHong Zhang if (nz > 0){ 51478910aadSHong Zhang bcol = bj + jmin; 51578910aadSHong Zhang bval = ba + jmin; 51678910aadSHong Zhang while (nz --) rtmp[*bcol++] += uikdi*(*bval++); 51778910aadSHong Zhang /* update il and jl for i-th row */ 51878910aadSHong Zhang il[i] = jmin; 51978910aadSHong Zhang j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; 52078910aadSHong Zhang } 52178910aadSHong Zhang i = nexti; 52278910aadSHong Zhang } 52378910aadSHong Zhang 5243cea4cbeSHong Zhang /* shift the diagonals when zero pivot is detected */ 5253cea4cbeSHong Zhang /* compute rs=sum of abs(off-diagonal) */ 5263cea4cbeSHong Zhang rs = 0.0; 52778910aadSHong Zhang jmin = bi[k]+1; 52878910aadSHong Zhang nz = bi[k+1] - jmin; 52978910aadSHong Zhang if (nz){ 53078910aadSHong Zhang bcol = bj + jmin; 53178910aadSHong Zhang while (nz--){ 53289f845c8SHong Zhang rs += PetscAbsScalar(rtmp[*bcol]); 53389f845c8SHong Zhang bcol++; 53478910aadSHong Zhang } 53578910aadSHong Zhang } 5363cea4cbeSHong Zhang 5373cea4cbeSHong Zhang sctx.rs = rs; 5383cea4cbeSHong Zhang sctx.pv = dk; 53945938b79SHong Zhang ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr); 54045938b79SHong Zhang if (newshift == 1) break; /* sctx.shift_amount is updated */ 54178910aadSHong Zhang 54278910aadSHong Zhang /* copy data into U(k,:) */ 54378910aadSHong Zhang ba[bi[k]] = 1.0/dk; 54478910aadSHong Zhang jmin = bi[k]+1; 54578910aadSHong Zhang nz = bi[k+1] - jmin; 54678910aadSHong Zhang if (nz){ 54778910aadSHong Zhang bcol = bj + jmin; 54878910aadSHong Zhang bval = ba + jmin; 54978910aadSHong Zhang while (nz--){ 55078910aadSHong Zhang *bval++ = rtmp[*bcol]; 55178910aadSHong Zhang rtmp[*bcol++] = 0.0; 55278910aadSHong Zhang } 55378910aadSHong Zhang /* add k-th row into il and jl */ 55478910aadSHong Zhang il[k] = jmin; 55578910aadSHong Zhang i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; 55678910aadSHong Zhang } 55778910aadSHong Zhang } 5583cea4cbeSHong Zhang } while (sctx.chshift); 55978910aadSHong Zhang ierr = PetscFree(il);CHKERRQ(ierr); 56078910aadSHong Zhang 561*719d5645SBarry Smith C->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 562*719d5645SBarry Smith C->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 56378910aadSHong Zhang C->assembled = PETSC_TRUE; 56478910aadSHong Zhang C->preallocated = PETSC_TRUE; 565d0f46423SBarry Smith ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr); 5663cea4cbeSHong Zhang if (sctx.nshift){ 5673cea4cbeSHong Zhang if (shiftnz) { 5681e2582c4SBarry Smith ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 5693cea4cbeSHong Zhang } else if (shiftpd) { 5701e2582c4SBarry Smith ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 57178910aadSHong Zhang } 57278910aadSHong Zhang } 573c05c3958SHong Zhang PetscFunctionReturn(0); 574c05c3958SHong Zhang } 575c05c3958SHong Zhang 576c05c3958SHong Zhang #include "petscbt.h" 577c05c3958SHong Zhang #include "src/mat/utils/freespace.h" 578c05c3958SHong Zhang #undef __FUNCT__ 579c05c3958SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqBAIJ" 580*719d5645SBarry Smith PetscErrorCode MatICCFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS perm,MatFactorInfo *info) 581c05c3958SHong Zhang { 58278910aadSHong Zhang Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 58378910aadSHong Zhang Mat_SeqSBAIJ *b; 58478910aadSHong Zhang Mat B; 58578910aadSHong Zhang PetscErrorCode ierr; 58678910aadSHong Zhang PetscTruth perm_identity; 587d0f46423SBarry Smith PetscInt reallocs=0,*rip,i,*ai=a->i,*aj=a->j,am=a->mbs,bs=A->rmap->bs,*ui; 58878910aadSHong Zhang PetscInt jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow; 589cfdb8b8aSHong Zhang PetscInt nlnk,*lnk,*lnk_lvl=PETSC_NULL,ncols,ncols_upper,*cols,*cols_lvl,*uj,**uj_ptr,**uj_lvl_ptr; 59078910aadSHong Zhang PetscReal fill=info->fill,levels=info->levels; 591a1a86e44SBarry Smith PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 592a1a86e44SBarry Smith PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL; 59378910aadSHong Zhang PetscBT lnkbt; 59478910aadSHong Zhang 595c05c3958SHong Zhang PetscFunctionBegin; 5966ad2eaddSHong Zhang if (bs > 1){ 5976ad2eaddSHong Zhang if (!a->sbaijMat){ 598ceb03754SKris Buschelman ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr); 5996ad2eaddSHong Zhang } 600*719d5645SBarry Smith (fact)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqSBAIJ; /* undue the change made in MatGetFactor_seqbaij_petsc */ 601*719d5645SBarry Smith ierr = MatICCFactorSymbolic(fact,a->sbaijMat,perm,info);CHKERRQ(ierr); 6026ad2eaddSHong Zhang PetscFunctionReturn(0); 6036ad2eaddSHong Zhang } 6046ad2eaddSHong Zhang 60578910aadSHong Zhang ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 60678910aadSHong Zhang ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 60778910aadSHong Zhang 60878910aadSHong Zhang /* special case that simply copies fill pattern */ 60978910aadSHong Zhang if (!levels && perm_identity) { 61078910aadSHong Zhang ierr = MatMarkDiagonal_SeqBAIJ(A);CHKERRQ(ierr); 61178910aadSHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 61278910aadSHong Zhang for (i=0; i<am; i++) { 61378910aadSHong Zhang ui[i] = ai[i+1] - a->diag[i]; /* ui: rowlengths - changes when !perm_identity */ 61478910aadSHong Zhang } 615*719d5645SBarry Smith B = fact; 61678910aadSHong Zhang ierr = MatSeqSBAIJSetPreallocation(B,1,0,ui);CHKERRQ(ierr); 61778910aadSHong Zhang 618b24902e0SBarry Smith 61978910aadSHong Zhang b = (Mat_SeqSBAIJ*)B->data; 62078910aadSHong Zhang uj = b->j; 62178910aadSHong Zhang for (i=0; i<am; i++) { 62278910aadSHong Zhang aj = a->j + a->diag[i]; 62378910aadSHong Zhang for (j=0; j<ui[i]; j++){ 62478910aadSHong Zhang *uj++ = *aj++; 62578910aadSHong Zhang } 62678910aadSHong Zhang b->ilen[i] = ui[i]; 62778910aadSHong Zhang } 62878910aadSHong Zhang ierr = PetscFree(ui);CHKERRQ(ierr); 629*719d5645SBarry Smith B->factor = MAT_FACTOR_NONE; 63078910aadSHong Zhang ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 63178910aadSHong Zhang ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 632*719d5645SBarry Smith B->factor = MAT_FACTOR_ICC; 63378910aadSHong Zhang 63478910aadSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering; 63578910aadSHong Zhang PetscFunctionReturn(0); 63678910aadSHong Zhang } 63778910aadSHong Zhang 63878910aadSHong Zhang /* initialization */ 63978910aadSHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 64078910aadSHong Zhang ui[0] = 0; 64178910aadSHong Zhang ierr = PetscMalloc((2*am+1)*sizeof(PetscInt),&cols_lvl);CHKERRQ(ierr); 64278910aadSHong Zhang 64378910aadSHong Zhang /* jl: linked list for storing indices of the pivot rows 64478910aadSHong Zhang il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */ 64578910aadSHong Zhang ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt*),&jl);CHKERRQ(ierr); 64678910aadSHong Zhang il = jl + am; 64778910aadSHong Zhang uj_ptr = (PetscInt**)(il + am); 64878910aadSHong Zhang uj_lvl_ptr = (PetscInt**)(uj_ptr + am); 64978910aadSHong Zhang for (i=0; i<am; i++){ 65078910aadSHong Zhang jl[i] = am; il[i] = 0; 65178910aadSHong Zhang } 65278910aadSHong Zhang 65378910aadSHong Zhang /* create and initialize a linked list for storing column indices of the active row k */ 65478910aadSHong Zhang nlnk = am + 1; 65578910aadSHong Zhang ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 65678910aadSHong Zhang 65778910aadSHong Zhang /* initial FreeSpace size is fill*(ai[am]+1) */ 658a1a86e44SBarry Smith ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr); 65978910aadSHong Zhang current_space = free_space; 660a1a86e44SBarry Smith ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr); 66178910aadSHong Zhang current_space_lvl = free_space_lvl; 66278910aadSHong Zhang 66378910aadSHong Zhang for (k=0; k<am; k++){ /* for each active row k */ 66478910aadSHong Zhang /* initialize lnk by the column indices of row rip[k] of A */ 66578910aadSHong Zhang nzk = 0; 66678910aadSHong Zhang ncols = ai[rip[k]+1] - ai[rip[k]]; 66778910aadSHong Zhang ncols_upper = 0; 66878910aadSHong Zhang cols = cols_lvl + am; 66978910aadSHong Zhang for (j=0; j<ncols; j++){ 67078910aadSHong Zhang i = rip[*(aj + ai[rip[k]] + j)]; 67178910aadSHong Zhang if (i >= k){ /* only take upper triangular entry */ 67278910aadSHong Zhang cols[ncols_upper] = i; 67378910aadSHong Zhang cols_lvl[ncols_upper] = -1; /* initialize level for nonzero entries */ 67478910aadSHong Zhang ncols_upper++; 67578910aadSHong Zhang } 67678910aadSHong Zhang } 67778910aadSHong Zhang ierr = PetscIncompleteLLAdd(ncols_upper,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 67878910aadSHong Zhang nzk += nlnk; 67978910aadSHong Zhang 68078910aadSHong Zhang /* update lnk by computing fill-in for each pivot row to be merged in */ 68178910aadSHong Zhang prow = jl[k]; /* 1st pivot row */ 68278910aadSHong Zhang 68378910aadSHong Zhang while (prow < k){ 68478910aadSHong Zhang nextprow = jl[prow]; 68578910aadSHong Zhang 68678910aadSHong Zhang /* merge prow into k-th row */ 68778910aadSHong Zhang jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:am-1) */ 68878910aadSHong Zhang jmax = ui[prow+1]; 68978910aadSHong Zhang ncols = jmax-jmin; 69078910aadSHong Zhang i = jmin - ui[prow]; 69178910aadSHong Zhang cols = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */ 69278910aadSHong Zhang for (j=0; j<ncols; j++) cols_lvl[j] = *(uj_lvl_ptr[prow] + i + j); 6935a8e39fbSHong Zhang ierr = PetscIncompleteLLAddSorted(ncols,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 69478910aadSHong Zhang nzk += nlnk; 69578910aadSHong Zhang 69678910aadSHong Zhang /* update il and jl for prow */ 69778910aadSHong Zhang if (jmin < jmax){ 69878910aadSHong Zhang il[prow] = jmin; 69978910aadSHong Zhang j = *cols; jl[prow] = jl[j]; jl[j] = prow; 70078910aadSHong Zhang } 70178910aadSHong Zhang prow = nextprow; 70278910aadSHong Zhang } 70378910aadSHong Zhang 70478910aadSHong Zhang /* if free space is not available, make more free space */ 70578910aadSHong Zhang if (current_space->local_remaining<nzk) { 70678910aadSHong Zhang i = am - k + 1; /* num of unfactored rows */ 70778910aadSHong Zhang i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 708a1a86e44SBarry Smith ierr = PetscFreeSpaceGet(i,¤t_space);CHKERRQ(ierr); 709a1a86e44SBarry Smith ierr = PetscFreeSpaceGet(i,¤t_space_lvl);CHKERRQ(ierr); 71078910aadSHong Zhang reallocs++; 71178910aadSHong Zhang } 71278910aadSHong Zhang 71378910aadSHong Zhang /* copy data into free_space and free_space_lvl, then initialize lnk */ 71478910aadSHong Zhang ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr); 71578910aadSHong Zhang 71678910aadSHong Zhang /* add the k-th row into il and jl */ 71778910aadSHong Zhang if (nzk-1 > 0){ 71878910aadSHong Zhang i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */ 71978910aadSHong Zhang jl[k] = jl[i]; jl[i] = k; 72078910aadSHong Zhang il[k] = ui[k] + 1; 72178910aadSHong Zhang } 72278910aadSHong Zhang uj_ptr[k] = current_space->array; 72378910aadSHong Zhang uj_lvl_ptr[k] = current_space_lvl->array; 72478910aadSHong Zhang 72578910aadSHong Zhang current_space->array += nzk; 72678910aadSHong Zhang current_space->local_used += nzk; 72778910aadSHong Zhang current_space->local_remaining -= nzk; 72878910aadSHong Zhang 72978910aadSHong Zhang current_space_lvl->array += nzk; 73078910aadSHong Zhang current_space_lvl->local_used += nzk; 73178910aadSHong Zhang current_space_lvl->local_remaining -= nzk; 73278910aadSHong Zhang 73378910aadSHong Zhang ui[k+1] = ui[k] + nzk; 73478910aadSHong Zhang } 73578910aadSHong Zhang 7366cf91177SBarry Smith #if defined(PETSC_USE_INFO) 73778910aadSHong Zhang if (ai[am] != 0) { 73878910aadSHong Zhang PetscReal af = ((PetscReal)(2*ui[am]-am))/((PetscReal)ai[am]); 739ae15b995SBarry Smith ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr); 740ae15b995SBarry Smith ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 741ae15b995SBarry Smith ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr); 74278910aadSHong Zhang } else { 743ae15b995SBarry Smith ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr); 74478910aadSHong Zhang } 74563ba0a88SBarry Smith #endif 74678910aadSHong Zhang 74778910aadSHong Zhang ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 74878910aadSHong Zhang ierr = PetscFree(jl);CHKERRQ(ierr); 74978910aadSHong Zhang ierr = PetscFree(cols_lvl);CHKERRQ(ierr); 75078910aadSHong Zhang 75178910aadSHong Zhang /* destroy list of free space and other temporary array(s) */ 75278910aadSHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 753a1a86e44SBarry Smith ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 75478910aadSHong Zhang ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 755a1a86e44SBarry Smith ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr); 75678910aadSHong Zhang 75778910aadSHong Zhang /* put together the new matrix in MATSEQSBAIJ format */ 758*719d5645SBarry Smith B = fact; 759ab93d7beSBarry Smith ierr = MatSeqSBAIJSetPreallocation(B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 76078910aadSHong Zhang 76178910aadSHong Zhang b = (Mat_SeqSBAIJ*)B->data; 76278910aadSHong Zhang b->singlemalloc = PETSC_FALSE; 763e6b907acSBarry Smith b->free_a = PETSC_TRUE; 764e6b907acSBarry Smith b->free_ij = PETSC_TRUE; 76578910aadSHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 76678910aadSHong Zhang b->j = uj; 76778910aadSHong Zhang b->i = ui; 76878910aadSHong Zhang b->diag = 0; 76978910aadSHong Zhang b->ilen = 0; 77078910aadSHong Zhang b->imax = 0; 77178910aadSHong Zhang b->row = perm; 77278910aadSHong Zhang b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 77378910aadSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 77478910aadSHong Zhang b->icol = perm; 77578910aadSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 77678910aadSHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 77752e6d16bSBarry Smith ierr = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 77878910aadSHong Zhang b->maxnz = b->nz = ui[am]; 77978910aadSHong Zhang 78078910aadSHong Zhang B->info.factor_mallocs = reallocs; 78178910aadSHong Zhang B->info.fill_ratio_given = fill; 78278910aadSHong Zhang if (ai[am] != 0) { 78378910aadSHong Zhang B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]); 78478910aadSHong Zhang } else { 78578910aadSHong Zhang B->info.fill_ratio_needed = 0.0; 78678910aadSHong Zhang } 78778910aadSHong Zhang if (perm_identity){ 78878910aadSHong Zhang B->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 78978910aadSHong Zhang B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 79078910aadSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering; 79178910aadSHong Zhang } else { 792*719d5645SBarry Smith (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N; 79378910aadSHong Zhang } 794c05c3958SHong Zhang PetscFunctionReturn(0); 795c05c3958SHong Zhang } 796c05c3958SHong Zhang 797c05c3958SHong Zhang #undef __FUNCT__ 798c05c3958SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqBAIJ" 799*719d5645SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS perm,MatFactorInfo *info) 800c05c3958SHong Zhang { 80178910aadSHong Zhang Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 80278910aadSHong Zhang Mat_SeqSBAIJ *b; 80378910aadSHong Zhang Mat B; 80478910aadSHong Zhang PetscErrorCode ierr; 80578910aadSHong Zhang PetscTruth perm_identity; 80678910aadSHong Zhang PetscReal fill = info->fill; 807d0f46423SBarry Smith PetscInt *rip,i,mbs=a->mbs,bs=A->rmap->bs,*ai=a->i,*aj=a->j,reallocs=0,prow; 80878910aadSHong Zhang PetscInt *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow; 80978910aadSHong Zhang PetscInt nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr; 810a1a86e44SBarry Smith PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 81178910aadSHong Zhang PetscBT lnkbt; 81278910aadSHong Zhang 813c05c3958SHong Zhang PetscFunctionBegin; 8146ad2eaddSHong Zhang if (bs > 1) { /* convert to seqsbaij */ 8156ad2eaddSHong Zhang if (!a->sbaijMat){ 816ceb03754SKris Buschelman ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr); 8176ad2eaddSHong Zhang } 818*719d5645SBarry Smith (fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqSBAIJ; /* undue the change made in MatGetFactor_seqbaij_petsc */ 819*719d5645SBarry Smith ierr = MatCholeskyFactorSymbolic(fact,a->sbaijMat,perm,info);CHKERRQ(ierr); 8206ad2eaddSHong Zhang PetscFunctionReturn(0); 8216ad2eaddSHong Zhang } 8226ad2eaddSHong Zhang 82378910aadSHong Zhang /* check whether perm is the identity mapping */ 82478910aadSHong Zhang ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 825c84f5b01SHong Zhang if (!perm_identity) SETERRQ(PETSC_ERR_SUP,"Matrix reordering is not supported"); 82678910aadSHong Zhang ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 82778910aadSHong Zhang 82878910aadSHong Zhang /* initialization */ 82978910aadSHong Zhang ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 83078910aadSHong Zhang ui[0] = 0; 83178910aadSHong Zhang 83278910aadSHong Zhang /* jl: linked list for storing indices of the pivot rows 83378910aadSHong Zhang il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */ 83478910aadSHong Zhang ierr = PetscMalloc((3*mbs+1)*sizeof(PetscInt)+mbs*sizeof(PetscInt*),&jl);CHKERRQ(ierr); 83578910aadSHong Zhang il = jl + mbs; 83678910aadSHong Zhang cols = il + mbs; 83778910aadSHong Zhang ui_ptr = (PetscInt**)(cols + mbs); 83878910aadSHong Zhang for (i=0; i<mbs; i++){ 83978910aadSHong Zhang jl[i] = mbs; il[i] = 0; 84078910aadSHong Zhang } 84178910aadSHong Zhang 84278910aadSHong Zhang /* create and initialize a linked list for storing column indices of the active row k */ 84378910aadSHong Zhang nlnk = mbs + 1; 84478910aadSHong Zhang ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 84578910aadSHong Zhang 84678910aadSHong Zhang /* initial FreeSpace size is fill*(ai[mbs]+1) */ 847a1a86e44SBarry Smith ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[mbs]+1)),&free_space);CHKERRQ(ierr); 84878910aadSHong Zhang current_space = free_space; 84978910aadSHong Zhang 85078910aadSHong Zhang for (k=0; k<mbs; k++){ /* for each active row k */ 85178910aadSHong Zhang /* initialize lnk by the column indices of row rip[k] of A */ 85278910aadSHong Zhang nzk = 0; 85378910aadSHong Zhang ncols = ai[rip[k]+1] - ai[rip[k]]; 85478910aadSHong Zhang ncols_upper = 0; 85578910aadSHong Zhang for (j=0; j<ncols; j++){ 85678910aadSHong Zhang i = rip[*(aj + ai[rip[k]] + j)]; 85778910aadSHong Zhang if (i >= k){ /* only take upper triangular entry */ 85878910aadSHong Zhang cols[ncols_upper] = i; 85978910aadSHong Zhang ncols_upper++; 86078910aadSHong Zhang } 86178910aadSHong Zhang } 86278910aadSHong Zhang ierr = PetscLLAdd(ncols_upper,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 86378910aadSHong Zhang nzk += nlnk; 86478910aadSHong Zhang 86578910aadSHong Zhang /* update lnk by computing fill-in for each pivot row to be merged in */ 86678910aadSHong Zhang prow = jl[k]; /* 1st pivot row */ 86778910aadSHong Zhang 86878910aadSHong Zhang while (prow < k){ 86978910aadSHong Zhang nextprow = jl[prow]; 87078910aadSHong Zhang /* merge prow into k-th row */ 87178910aadSHong Zhang jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:mbs-1) */ 87278910aadSHong Zhang jmax = ui[prow+1]; 87378910aadSHong Zhang ncols = jmax-jmin; 87478910aadSHong Zhang uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */ 8755a8e39fbSHong Zhang ierr = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 87678910aadSHong Zhang nzk += nlnk; 87778910aadSHong Zhang 87878910aadSHong Zhang /* update il and jl for prow */ 87978910aadSHong Zhang if (jmin < jmax){ 88078910aadSHong Zhang il[prow] = jmin; 88178910aadSHong Zhang j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow; 88278910aadSHong Zhang } 88378910aadSHong Zhang prow = nextprow; 88478910aadSHong Zhang } 88578910aadSHong Zhang 88678910aadSHong Zhang /* if free space is not available, make more free space */ 88778910aadSHong Zhang if (current_space->local_remaining<nzk) { 88878910aadSHong Zhang i = mbs - k + 1; /* num of unfactored rows */ 88978910aadSHong Zhang i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 890a1a86e44SBarry Smith ierr = PetscFreeSpaceGet(i,¤t_space);CHKERRQ(ierr); 89178910aadSHong Zhang reallocs++; 89278910aadSHong Zhang } 89378910aadSHong Zhang 89478910aadSHong Zhang /* copy data into free space, then initialize lnk */ 89578910aadSHong Zhang ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 89678910aadSHong Zhang 89778910aadSHong Zhang /* add the k-th row into il and jl */ 89878910aadSHong Zhang if (nzk-1 > 0){ 89978910aadSHong Zhang i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */ 90078910aadSHong Zhang jl[k] = jl[i]; jl[i] = k; 90178910aadSHong Zhang il[k] = ui[k] + 1; 90278910aadSHong Zhang } 90378910aadSHong Zhang ui_ptr[k] = current_space->array; 90478910aadSHong Zhang current_space->array += nzk; 90578910aadSHong Zhang current_space->local_used += nzk; 90678910aadSHong Zhang current_space->local_remaining -= nzk; 90778910aadSHong Zhang 90878910aadSHong Zhang ui[k+1] = ui[k] + nzk; 90978910aadSHong Zhang } 91078910aadSHong Zhang 9116cf91177SBarry Smith #if defined(PETSC_USE_INFO) 91278910aadSHong Zhang if (ai[mbs] != 0) { 91378910aadSHong Zhang PetscReal af = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]); 914ae15b995SBarry Smith ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr); 915ae15b995SBarry Smith ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 916ae15b995SBarry Smith ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr); 91778910aadSHong Zhang } else { 918ae15b995SBarry Smith ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr); 91978910aadSHong Zhang } 92063ba0a88SBarry Smith #endif 92178910aadSHong Zhang 92278910aadSHong Zhang ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 92378910aadSHong Zhang ierr = PetscFree(jl);CHKERRQ(ierr); 92478910aadSHong Zhang 92578910aadSHong Zhang /* destroy list of free space and other temporary array(s) */ 92678910aadSHong Zhang ierr = PetscMalloc((ui[mbs]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 927a1a86e44SBarry Smith ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 92878910aadSHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 92978910aadSHong Zhang 93078910aadSHong Zhang /* put together the new matrix in MATSEQSBAIJ format */ 931*719d5645SBarry Smith B = fact; 932ab93d7beSBarry Smith ierr = MatSeqSBAIJSetPreallocation(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 93378910aadSHong Zhang 93478910aadSHong Zhang b = (Mat_SeqSBAIJ*)B->data; 93578910aadSHong Zhang b->singlemalloc = PETSC_FALSE; 936e6b907acSBarry Smith b->free_a = PETSC_TRUE; 937e6b907acSBarry Smith b->free_ij = PETSC_TRUE; 93878910aadSHong Zhang ierr = PetscMalloc((ui[mbs]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 93978910aadSHong Zhang b->j = uj; 94078910aadSHong Zhang b->i = ui; 94178910aadSHong Zhang b->diag = 0; 94278910aadSHong Zhang b->ilen = 0; 94378910aadSHong Zhang b->imax = 0; 94478910aadSHong Zhang b->row = perm; 94578910aadSHong Zhang b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 94678910aadSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 94778910aadSHong Zhang b->icol = perm; 94878910aadSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 94978910aadSHong Zhang ierr = PetscMalloc((mbs+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 95052e6d16bSBarry Smith ierr = PetscLogObjectMemory(B,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 95178910aadSHong Zhang b->maxnz = b->nz = ui[mbs]; 95278910aadSHong Zhang 95378910aadSHong Zhang B->info.factor_mallocs = reallocs; 95478910aadSHong Zhang B->info.fill_ratio_given = fill; 95578910aadSHong Zhang if (ai[mbs] != 0) { 95678910aadSHong Zhang B->info.fill_ratio_needed = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]); 95778910aadSHong Zhang } else { 95878910aadSHong Zhang B->info.fill_ratio_needed = 0.0; 95978910aadSHong Zhang } 96078910aadSHong Zhang if (perm_identity){ 9616ad2eaddSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering; 96278910aadSHong Zhang } else { 9636ad2eaddSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N; 96478910aadSHong Zhang } 965c05c3958SHong Zhang PetscFunctionReturn(0); 966c05c3958SHong Zhang } 967