15d517e7eSBarry Smith /* 2ec3a8b7bSBarry Smith Factorization code for BAIJ format. 35d517e7eSBarry Smith */ 470f55243SBarry Smith #include "src/mat/impls/baij/seq/baij.h" 54078e994SBarry Smith #include "src/inline/ilu.h" 6ec3a8b7bSBarry Smith 76d84be18SBarry Smith /* ------------------------------------------------------------*/ 86d84be18SBarry Smith /* 94eeb42bcSBarry Smith Version for when blocks are 2 by 2 104eeb42bcSBarry Smith */ 114a2ae208SSatish Balay #undef __FUNCT__ 124a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2" 13af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2(Mat A,MatFactorInfo *info,Mat *B) 144eeb42bcSBarry Smith { 154eeb42bcSBarry Smith Mat C = *B; 164eeb42bcSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; 177cdf2dbbSSatish Balay IS isrow = b->row,isicol = b->icol; 186849ba73SBarry Smith PetscErrorCode ierr; 19690b6cddSBarry Smith PetscInt *r,*ic,i,j,n = a->mbs,*bi = b->i,*bj = b->j; 20690b6cddSBarry Smith PetscInt *ajtmpold,*ajtmp,nz,row; 21690b6cddSBarry Smith PetscInt *diag_offset=b->diag,idx,*ai=a->i,*aj=a->j,*pj; 22329f5518SBarry Smith MatScalar *pv,*v,*rtmp,m1,m2,m3,m4,*pc,*w,*x,x1,x2,x3,x4; 232515f8d2SSatish Balay MatScalar p1,p2,p3,p4; 243f1db9ecSBarry Smith MatScalar *ba = b->a,*aa = a->a; 254eeb42bcSBarry Smith 263a40ed3dSBarry Smith PetscFunctionBegin; 274eeb42bcSBarry Smith ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 284eeb42bcSBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 29b0a32e0cSBarry Smith ierr = PetscMalloc(4*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 304eeb42bcSBarry Smith 314eeb42bcSBarry Smith for (i=0; i<n; i++) { 324078e994SBarry Smith nz = bi[i+1] - bi[i]; 334078e994SBarry Smith ajtmp = bj + bi[i]; 344eeb42bcSBarry Smith for (j=0; j<nz; j++) { 354eeb42bcSBarry Smith x = rtmp+4*ajtmp[j]; x[0] = x[1] = x[2] = x[3] = 0.0; 364eeb42bcSBarry Smith } 374eeb42bcSBarry Smith /* load in initial (unfactored row) */ 384eeb42bcSBarry Smith idx = r[i]; 394078e994SBarry Smith nz = ai[idx+1] - ai[idx]; 404078e994SBarry Smith ajtmpold = aj + ai[idx]; 414078e994SBarry Smith v = aa + 4*ai[idx]; 424eeb42bcSBarry Smith for (j=0; j<nz; j++) { 434eeb42bcSBarry Smith x = rtmp+4*ic[ajtmpold[j]]; 444eeb42bcSBarry Smith x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3]; 454eeb42bcSBarry Smith v += 4; 464eeb42bcSBarry Smith } 474eeb42bcSBarry Smith row = *ajtmp++; 484eeb42bcSBarry Smith while (row < i) { 494eeb42bcSBarry Smith pc = rtmp + 4*row; 504eeb42bcSBarry Smith p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3]; 5188685aaeSLois Curfman McInnes if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0) { 524078e994SBarry Smith pv = ba + 4*diag_offset[row]; 534078e994SBarry Smith pj = bj + diag_offset[row] + 1; 544eeb42bcSBarry Smith x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 554eeb42bcSBarry Smith pc[0] = m1 = p1*x1 + p3*x2; 564eeb42bcSBarry Smith pc[1] = m2 = p2*x1 + p4*x2; 574eeb42bcSBarry Smith pc[2] = m3 = p1*x3 + p3*x4; 584eeb42bcSBarry Smith pc[3] = m4 = p2*x3 + p4*x4; 594078e994SBarry Smith nz = bi[row+1] - diag_offset[row] - 1; 604eeb42bcSBarry Smith pv += 4; 614eeb42bcSBarry Smith for (j=0; j<nz; j++) { 624eeb42bcSBarry Smith x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 634eeb42bcSBarry Smith x = rtmp + 4*pj[j]; 644eeb42bcSBarry Smith x[0] -= m1*x1 + m3*x2; 654eeb42bcSBarry Smith x[1] -= m2*x1 + m4*x2; 664eeb42bcSBarry Smith x[2] -= m1*x3 + m3*x4; 674eeb42bcSBarry Smith x[3] -= m2*x3 + m4*x4; 684eeb42bcSBarry Smith pv += 4; 694eeb42bcSBarry Smith } 70b0a32e0cSBarry Smith PetscLogFlops(16*nz+12); 714eeb42bcSBarry Smith } 724eeb42bcSBarry Smith row = *ajtmp++; 734eeb42bcSBarry Smith } 744eeb42bcSBarry Smith /* finished row so stick it into b->a */ 754078e994SBarry Smith pv = ba + 4*bi[i]; 764078e994SBarry Smith pj = bj + bi[i]; 774078e994SBarry Smith nz = bi[i+1] - bi[i]; 784eeb42bcSBarry Smith for (j=0; j<nz; j++) { 794eeb42bcSBarry Smith x = rtmp+4*pj[j]; 804eeb42bcSBarry Smith pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3]; 814eeb42bcSBarry Smith pv += 4; 824eeb42bcSBarry Smith } 834eeb42bcSBarry Smith /* invert diagonal block */ 844078e994SBarry Smith w = ba + 4*diag_offset[i]; 854cd38560SBarry Smith ierr = Kernel_A_gets_inverse_A_2(w);CHKERRQ(ierr); 864eeb42bcSBarry Smith } 874eeb42bcSBarry Smith 88606d414cSSatish Balay ierr = PetscFree(rtmp);CHKERRQ(ierr); 894eeb42bcSBarry Smith ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 904eeb42bcSBarry Smith ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 914eeb42bcSBarry Smith C->factor = FACTOR_LU; 924eeb42bcSBarry Smith C->assembled = PETSC_TRUE; 93b0a32e0cSBarry Smith PetscLogFlops(1.3333*8*b->mbs); /* from inverting diagonal blocks */ 943a40ed3dSBarry Smith PetscFunctionReturn(0); 954eeb42bcSBarry Smith } 964cd38560SBarry Smith /* 974cd38560SBarry Smith Version for when blocks are 2 by 2 Using natural ordering 984cd38560SBarry Smith */ 994a2ae208SSatish Balay #undef __FUNCT__ 1004a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering" 101af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering(Mat A,MatFactorInfo *info,Mat *B) 1024cd38560SBarry Smith { 1034cd38560SBarry Smith Mat C = *B; 1044cd38560SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; 105dfbe8321SBarry Smith PetscErrorCode ierr; 106690b6cddSBarry Smith PetscInt i,j,n = a->mbs,*bi = b->i,*bj = b->j; 107690b6cddSBarry Smith PetscInt *ajtmpold,*ajtmp,nz,row; 108690b6cddSBarry Smith PetscInt *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj; 109329f5518SBarry Smith MatScalar *pv,*v,*rtmp,*pc,*w,*x; 1102515f8d2SSatish Balay MatScalar p1,p2,p3,p4,m1,m2,m3,m4,x1,x2,x3,x4; 1114cd38560SBarry Smith MatScalar *ba = b->a,*aa = a->a; 1124cd38560SBarry Smith 1134cd38560SBarry Smith PetscFunctionBegin; 114b0a32e0cSBarry Smith ierr = PetscMalloc(4*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 1154cd38560SBarry Smith 1164cd38560SBarry Smith for (i=0; i<n; i++) { 1174cd38560SBarry Smith nz = bi[i+1] - bi[i]; 1184cd38560SBarry Smith ajtmp = bj + bi[i]; 1194cd38560SBarry Smith for (j=0; j<nz; j++) { 1204cd38560SBarry Smith x = rtmp+4*ajtmp[j]; 1214cd38560SBarry Smith x[0] = x[1] = x[2] = x[3] = 0.0; 1224cd38560SBarry Smith } 1234cd38560SBarry Smith /* load in initial (unfactored row) */ 1244cd38560SBarry Smith nz = ai[i+1] - ai[i]; 1254cd38560SBarry Smith ajtmpold = aj + ai[i]; 1264cd38560SBarry Smith v = aa + 4*ai[i]; 1274cd38560SBarry Smith for (j=0; j<nz; j++) { 1284cd38560SBarry Smith x = rtmp+4*ajtmpold[j]; 1294cd38560SBarry Smith x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3]; 1304cd38560SBarry Smith v += 4; 1314cd38560SBarry Smith } 1324cd38560SBarry Smith row = *ajtmp++; 1334cd38560SBarry Smith while (row < i) { 1344cd38560SBarry Smith pc = rtmp + 4*row; 1354cd38560SBarry Smith p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3]; 1364cd38560SBarry Smith if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0) { 1374cd38560SBarry Smith pv = ba + 4*diag_offset[row]; 1384cd38560SBarry Smith pj = bj + diag_offset[row] + 1; 1394cd38560SBarry Smith x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 1404cd38560SBarry Smith pc[0] = m1 = p1*x1 + p3*x2; 1414cd38560SBarry Smith pc[1] = m2 = p2*x1 + p4*x2; 1424cd38560SBarry Smith pc[2] = m3 = p1*x3 + p3*x4; 1434cd38560SBarry Smith pc[3] = m4 = p2*x3 + p4*x4; 1444cd38560SBarry Smith nz = bi[row+1] - diag_offset[row] - 1; 1454cd38560SBarry Smith pv += 4; 1464cd38560SBarry Smith for (j=0; j<nz; j++) { 1474cd38560SBarry Smith x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 1484cd38560SBarry Smith x = rtmp + 4*pj[j]; 1494cd38560SBarry Smith x[0] -= m1*x1 + m3*x2; 1504cd38560SBarry Smith x[1] -= m2*x1 + m4*x2; 1514cd38560SBarry Smith x[2] -= m1*x3 + m3*x4; 1524cd38560SBarry Smith x[3] -= m2*x3 + m4*x4; 1534cd38560SBarry Smith pv += 4; 1544cd38560SBarry Smith } 155b0a32e0cSBarry Smith PetscLogFlops(16*nz+12); 1564cd38560SBarry Smith } 1574cd38560SBarry Smith row = *ajtmp++; 1584cd38560SBarry Smith } 1594cd38560SBarry Smith /* finished row so stick it into b->a */ 1604cd38560SBarry Smith pv = ba + 4*bi[i]; 1614cd38560SBarry Smith pj = bj + bi[i]; 1624cd38560SBarry Smith nz = bi[i+1] - bi[i]; 1634cd38560SBarry Smith for (j=0; j<nz; j++) { 1644cd38560SBarry Smith x = rtmp+4*pj[j]; 1654cd38560SBarry Smith pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3]; 1664cd38560SBarry Smith pv += 4; 1674cd38560SBarry Smith } 1684cd38560SBarry Smith /* invert diagonal block */ 1694cd38560SBarry Smith w = ba + 4*diag_offset[i]; 1704cd38560SBarry Smith ierr = Kernel_A_gets_inverse_A_2(w);CHKERRQ(ierr); 1714cd38560SBarry Smith /*Kernel_A_gets_inverse_A(bs,w,v_pivots,v_work);*/ 1724cd38560SBarry Smith } 1734cd38560SBarry Smith 174606d414cSSatish Balay ierr = PetscFree(rtmp);CHKERRQ(ierr); 1754cd38560SBarry Smith C->factor = FACTOR_LU; 1764cd38560SBarry Smith C->assembled = PETSC_TRUE; 177b0a32e0cSBarry Smith PetscLogFlops(1.3333*8*b->mbs); /* from inverting diagonal blocks */ 1784cd38560SBarry Smith PetscFunctionReturn(0); 1794cd38560SBarry Smith } 1807fc0212eSBarry Smith 1817fc0212eSBarry Smith /* ----------------------------------------------------------- */ 1827fc0212eSBarry Smith /* 1837fc0212eSBarry Smith Version for when blocks are 1 by 1. 1847fc0212eSBarry Smith */ 1854a2ae208SSatish Balay #undef __FUNCT__ 1864a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_1" 187af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqBAIJ_1(Mat A,MatFactorInfo *info,Mat *B) 1887fc0212eSBarry Smith { 1897fc0212eSBarry Smith Mat C = *B; 1907fc0212eSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; 1917cdf2dbbSSatish Balay IS isrow = b->row,isicol = b->icol; 1926849ba73SBarry Smith PetscErrorCode ierr; 193690b6cddSBarry Smith PetscInt *r,*ic,i,j,n = a->mbs,*bi = b->i,*bj = b->j; 194690b6cddSBarry Smith PetscInt *ajtmpold,*ajtmp,nz,row,*ai = a->i,*aj = a->j; 195690b6cddSBarry Smith PetscInt *diag_offset = b->diag,diag,*pj; 196329f5518SBarry Smith MatScalar *pv,*v,*rtmp,multiplier,*pc; 1973f1db9ecSBarry Smith MatScalar *ba = b->a,*aa = a->a; 1987fc0212eSBarry Smith 1993a40ed3dSBarry Smith PetscFunctionBegin; 2007fc0212eSBarry Smith ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 2017fc0212eSBarry Smith ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 202b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 2037fc0212eSBarry Smith 2047fc0212eSBarry Smith for (i=0; i<n; i++) { 2054078e994SBarry Smith nz = bi[i+1] - bi[i]; 2064078e994SBarry Smith ajtmp = bj + bi[i]; 2077fc0212eSBarry Smith for (j=0; j<nz; j++) rtmp[ajtmp[j]] = 0.0; 2087fc0212eSBarry Smith 2097fc0212eSBarry Smith /* load in initial (unfactored row) */ 2104078e994SBarry Smith nz = ai[r[i]+1] - ai[r[i]]; 2114078e994SBarry Smith ajtmpold = aj + ai[r[i]]; 2124078e994SBarry Smith v = aa + ai[r[i]]; 2137fc0212eSBarry Smith for (j=0; j<nz; j++) rtmp[ic[ajtmpold[j]]] = v[j]; 2147fc0212eSBarry Smith 2157fc0212eSBarry Smith row = *ajtmp++; 2167fc0212eSBarry Smith while (row < i) { 2177fc0212eSBarry Smith pc = rtmp + row; 2187fc0212eSBarry Smith if (*pc != 0.0) { 2194078e994SBarry Smith pv = ba + diag_offset[row]; 2204078e994SBarry Smith pj = bj + diag_offset[row] + 1; 2217fc0212eSBarry Smith multiplier = *pc * *pv++; 2227fc0212eSBarry Smith *pc = multiplier; 2234078e994SBarry Smith nz = bi[row+1] - diag_offset[row] - 1; 2247fc0212eSBarry Smith for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j]; 225b0a32e0cSBarry Smith PetscLogFlops(1+2*nz); 2267fc0212eSBarry Smith } 2277fc0212eSBarry Smith row = *ajtmp++; 2287fc0212eSBarry Smith } 2297fc0212eSBarry Smith /* finished row so stick it into b->a */ 2304078e994SBarry Smith pv = ba + bi[i]; 2314078e994SBarry Smith pj = bj + bi[i]; 2324078e994SBarry Smith nz = bi[i+1] - bi[i]; 2337fc0212eSBarry Smith for (j=0; j<nz; j++) {pv[j] = rtmp[pj[j]];} 2344078e994SBarry Smith diag = diag_offset[i] - bi[i]; 2357fc0212eSBarry Smith /* check pivot entry for current row */ 236a73cf429SBarry Smith if (pv[diag] == 0.0) { 23729bbc08cSBarry Smith SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot"); 2387fc0212eSBarry Smith } 2397fc0212eSBarry Smith pv[diag] = 1.0/pv[diag]; 2407fc0212eSBarry Smith } 2417fc0212eSBarry Smith 242606d414cSSatish Balay ierr = PetscFree(rtmp);CHKERRQ(ierr); 2437fc0212eSBarry Smith ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 2447fc0212eSBarry Smith ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 2457fc0212eSBarry Smith C->factor = FACTOR_LU; 2467fc0212eSBarry Smith C->assembled = PETSC_TRUE; 247b0a32e0cSBarry Smith PetscLogFlops(C->n); 2483a40ed3dSBarry Smith PetscFunctionReturn(0); 2497fc0212eSBarry Smith } 2507fc0212eSBarry Smith 251273d9f13SBarry Smith 2525d517e7eSBarry Smith /* ----------------------------------------------------------- */ 2534a2ae208SSatish Balay #undef __FUNCT__ 2544a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqBAIJ" 255dfbe8321SBarry Smith PetscErrorCode MatLUFactor_SeqBAIJ(Mat A,IS row,IS col,MatFactorInfo *info) 2565d517e7eSBarry Smith { 257dfbe8321SBarry Smith PetscErrorCode ierr; 2585d517e7eSBarry Smith Mat C; 2595d517e7eSBarry Smith 2603a40ed3dSBarry Smith PetscFunctionBegin; 261b9b97703SBarry Smith ierr = MatLUFactorSymbolic(A,row,col,info,&C);CHKERRQ(ierr); 262af281ebdSHong Zhang ierr = MatLUFactorNumeric(A,info,&C);CHKERRQ(ierr); 263273d9f13SBarry Smith ierr = MatHeaderCopy(A,C);CHKERRQ(ierr); 264b0a32e0cSBarry Smith PetscLogObjectParent(A,((Mat_SeqBAIJ*)(A->data))->icol); 2653a40ed3dSBarry Smith PetscFunctionReturn(0); 2665d517e7eSBarry Smith } 2674cd38560SBarry Smith 268c05c3958SHong Zhang #include "src/mat/impls/sbaij/seq/sbaij.h" 269c05c3958SHong Zhang #undef __FUNCT__ 270c05c3958SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N" 271af281ebdSHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N(Mat A,MatFactorInfo *info,Mat *B) 272c05c3958SHong Zhang { 273f3086b4bSHong Zhang PetscErrorCode ierr; 27478910aadSHong Zhang Mat C = *B; 27578910aadSHong Zhang Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data; 27678910aadSHong Zhang Mat_SeqSBAIJ *b=(Mat_SeqSBAIJ*)C->data; 27778910aadSHong Zhang IS ip=b->row; 27878910aadSHong Zhang PetscInt *rip,i,j,mbs=a->mbs,bs=A->bs,*bi=b->i,*bj=b->j,*bcol; 27978910aadSHong Zhang PetscInt *ai=a->i,*aj=a->j; 28078910aadSHong Zhang PetscInt k,jmin,jmax,*jl,*il,col,nexti,ili,nz; 28178910aadSHong Zhang MatScalar *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi; 2823cea4cbeSHong Zhang PetscReal zeropivot,rs,shiftnz; 2833cea4cbeSHong Zhang PetscTruth shiftpd; 2843cea4cbeSHong Zhang ChShift_Ctx sctx; 2853cea4cbeSHong Zhang PetscInt newshift; 28678910aadSHong Zhang 287c05c3958SHong Zhang PetscFunctionBegin; 2886ad2eaddSHong Zhang if (bs > 1) { 2896ad2eaddSHong Zhang if (!a->sbaijMat){ 2906ad2eaddSHong Zhang ierr = MatConvert(A,MATSEQSBAIJ,&a->sbaijMat);CHKERRQ(ierr); 2916ad2eaddSHong Zhang } 292f3086b4bSHong Zhang ierr = (a->sbaijMat)->ops->choleskyfactornumeric(a->sbaijMat,info,B);CHKERRQ(ierr); 2936ad2eaddSHong Zhang ierr = MatDestroy(a->sbaijMat);CHKERRQ(ierr); 2946ad2eaddSHong Zhang a->sbaijMat = PETSC_NULL; 2956ad2eaddSHong Zhang PetscFunctionReturn(0); 2966ad2eaddSHong Zhang } 29778910aadSHong Zhang 29878910aadSHong Zhang /* initialization */ 2993cea4cbeSHong Zhang shiftnz = info->shiftnz; 3003cea4cbeSHong Zhang shiftpd = info->shiftpd; 3013cea4cbeSHong Zhang zeropivot = info->zeropivot; 3023cea4cbeSHong Zhang 3036ad2eaddSHong Zhang ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr); 30478910aadSHong Zhang nz = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar); 30578910aadSHong Zhang ierr = PetscMalloc(nz,&il);CHKERRQ(ierr); 30678910aadSHong Zhang jl = il + mbs; 30778910aadSHong Zhang rtmp = (MatScalar*)(jl + mbs); 30878910aadSHong Zhang 3093cea4cbeSHong Zhang sctx.shift_amount = 0; 3103cea4cbeSHong Zhang sctx.nshift = 0; 31178910aadSHong Zhang do { 3123cea4cbeSHong Zhang sctx.chshift = PETSC_FALSE; 31378910aadSHong Zhang for (i=0; i<mbs; i++) { 31478910aadSHong Zhang rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0; 31578910aadSHong Zhang } 31678910aadSHong Zhang 31778910aadSHong Zhang for (k = 0; k<mbs; k++){ 31878910aadSHong Zhang bval = ba + bi[k]; 31978910aadSHong Zhang /* initialize k-th row by the perm[k]-th row of A */ 32078910aadSHong Zhang jmin = ai[rip[k]]; jmax = ai[rip[k]+1]; 32178910aadSHong Zhang for (j = jmin; j < jmax; j++){ 32278910aadSHong Zhang col = rip[aj[j]]; 32378910aadSHong Zhang if (col >= k){ /* only take upper triangular entry */ 32478910aadSHong Zhang rtmp[col] = aa[j]; 32578910aadSHong Zhang *bval++ = 0.0; /* for in-place factorization */ 32678910aadSHong Zhang } 32778910aadSHong Zhang } 3283cea4cbeSHong Zhang 3293cea4cbeSHong Zhang /* shift the diagonal of the matrix */ 3303cea4cbeSHong Zhang if (sctx.nshift) rtmp[k] += sctx.shift_amount; 33178910aadSHong Zhang 33278910aadSHong Zhang /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 33378910aadSHong Zhang dk = rtmp[k]; 33478910aadSHong Zhang i = jl[k]; /* first row to be added to k_th row */ 33578910aadSHong Zhang 33678910aadSHong Zhang while (i < k){ 33778910aadSHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 33878910aadSHong Zhang 33978910aadSHong Zhang /* compute multiplier, update diag(k) and U(i,k) */ 34078910aadSHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 34178910aadSHong Zhang uikdi = - ba[ili]*ba[bi[i]]; /* diagonal(k) */ 34278910aadSHong Zhang dk += uikdi*ba[ili]; 34378910aadSHong Zhang ba[ili] = uikdi; /* -U(i,k) */ 34478910aadSHong Zhang 34578910aadSHong Zhang /* add multiple of row i to k-th row */ 34678910aadSHong Zhang jmin = ili + 1; jmax = bi[i+1]; 34778910aadSHong Zhang if (jmin < jmax){ 34878910aadSHong Zhang for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j]; 34978910aadSHong Zhang /* update il and jl for row i */ 35078910aadSHong Zhang il[i] = jmin; 35178910aadSHong Zhang j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; 35278910aadSHong Zhang } 35378910aadSHong Zhang i = nexti; 35478910aadSHong Zhang } 35578910aadSHong Zhang 3563cea4cbeSHong Zhang /* shift the diagonals when zero pivot is detected */ 3573cea4cbeSHong Zhang /* compute rs=sum of abs(off-diagonal) */ 3583cea4cbeSHong Zhang rs = 0.0; 35978910aadSHong Zhang jmin = bi[k]+1; 36078910aadSHong Zhang nz = bi[k+1] - jmin; 36178910aadSHong Zhang if (nz){ 36278910aadSHong Zhang bcol = bj + jmin; 36378910aadSHong Zhang while (nz--){ 364*89f845c8SHong Zhang rs += PetscAbsScalar(rtmp[*bcol]); 365*89f845c8SHong Zhang bcol++; 36678910aadSHong Zhang } 36778910aadSHong Zhang } 3683cea4cbeSHong Zhang 3693cea4cbeSHong Zhang sctx.rs = rs; 3703cea4cbeSHong Zhang sctx.pv = dk; 3713cea4cbeSHong Zhang ierr = Mat_CholeskyCheckShift(info,&sctx,&newshift);CHKERRQ(ierr); 3723cea4cbeSHong Zhang if (newshift == 1){ 3733cea4cbeSHong Zhang break; /* sctx.shift_amount is updated */ 3743cea4cbeSHong Zhang } else if (newshift == -1){ 3753cea4cbeSHong Zhang SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g * rs %g",k,PetscAbsScalar(dk),zeropivot,rs); 37678910aadSHong Zhang } 37778910aadSHong Zhang 37878910aadSHong Zhang /* copy data into U(k,:) */ 37978910aadSHong Zhang ba[bi[k]] = 1.0/dk; /* U(k,k) */ 38078910aadSHong Zhang jmin = bi[k]+1; jmax = bi[k+1]; 38178910aadSHong Zhang if (jmin < jmax) { 38278910aadSHong Zhang for (j=jmin; j<jmax; j++){ 38378910aadSHong Zhang col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0; 38478910aadSHong Zhang } 38578910aadSHong Zhang /* add the k-th row into il and jl */ 38678910aadSHong Zhang il[k] = jmin; 38778910aadSHong Zhang i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; 38878910aadSHong Zhang } 38978910aadSHong Zhang } 3903cea4cbeSHong Zhang } while (sctx.chshift); 39178910aadSHong Zhang ierr = PetscFree(il);CHKERRQ(ierr); 39278910aadSHong Zhang 39378910aadSHong Zhang ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr); 39478910aadSHong Zhang C->factor = FACTOR_CHOLESKY; 39578910aadSHong Zhang C->assembled = PETSC_TRUE; 39678910aadSHong Zhang C->preallocated = PETSC_TRUE; 39778910aadSHong Zhang PetscLogFlops(C->m); 3983cea4cbeSHong Zhang if (sctx.nshift){ 3993cea4cbeSHong Zhang if (shiftnz) { 4003cea4cbeSHong Zhang PetscLogInfo(0,"MatCholeskyFactorNumeric_SeqBAIJ_1: number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount); 4013cea4cbeSHong Zhang } else if (shiftpd) { 4023cea4cbeSHong Zhang PetscLogInfo(0,"MatCholeskyFactorNumeric_SeqBAIJ_1: number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount); 40378910aadSHong Zhang } 40478910aadSHong Zhang } 405c05c3958SHong Zhang PetscFunctionReturn(0); 406c05c3958SHong Zhang } 4074cd38560SBarry Smith 408c05c3958SHong Zhang #undef __FUNCT__ 409c05c3958SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering" 410af281ebdSHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering(Mat A,MatFactorInfo *info,Mat *fact) 411c05c3958SHong Zhang { 41278910aadSHong Zhang Mat C = *fact; 41378910aadSHong Zhang Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data; 41478910aadSHong Zhang Mat_SeqSBAIJ *b=(Mat_SeqSBAIJ*)C->data; 41578910aadSHong Zhang PetscErrorCode ierr; 4163cea4cbeSHong Zhang PetscInt i,j,am=a->mbs; 41778910aadSHong Zhang PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 4183cea4cbeSHong Zhang PetscInt k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz; 41978910aadSHong Zhang MatScalar *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval; 4203cea4cbeSHong Zhang PetscReal zeropivot,rs,shiftnz; 4213cea4cbeSHong Zhang PetscTruth shiftpd; 4223cea4cbeSHong Zhang ChShift_Ctx sctx; 4233cea4cbeSHong Zhang PetscInt newshift; 42478910aadSHong Zhang 425c05c3958SHong Zhang PetscFunctionBegin; 42678910aadSHong Zhang /* initialization */ 4273cea4cbeSHong Zhang shiftnz = info->shiftnz; 4283cea4cbeSHong Zhang shiftpd = info->shiftpd; 4293cea4cbeSHong Zhang zeropivot = info->zeropivot; 4303cea4cbeSHong Zhang 43178910aadSHong Zhang nz = (2*am+1)*sizeof(PetscInt)+am*sizeof(MatScalar); 43278910aadSHong Zhang ierr = PetscMalloc(nz,&il);CHKERRQ(ierr); 43378910aadSHong Zhang jl = il + am; 43478910aadSHong Zhang rtmp = (MatScalar*)(jl + am); 43578910aadSHong Zhang 4363cea4cbeSHong Zhang sctx.shift_amount = 0; 4373cea4cbeSHong Zhang sctx.nshift = 0; 43878910aadSHong Zhang do { 4393cea4cbeSHong Zhang sctx.chshift = PETSC_FALSE; 44078910aadSHong Zhang for (i=0; i<am; i++) { 44178910aadSHong Zhang rtmp[i] = 0.0; jl[i] = am; il[0] = 0; 44278910aadSHong Zhang } 44378910aadSHong Zhang 44478910aadSHong Zhang for (k = 0; k<am; k++){ 44578910aadSHong Zhang /* initialize k-th row with elements nonzero in row perm(k) of A */ 44678910aadSHong Zhang nz = ai[k+1] - ai[k]; 44778910aadSHong Zhang acol = aj + ai[k]; 44878910aadSHong Zhang aval = aa + ai[k]; 44978910aadSHong Zhang bval = ba + bi[k]; 45078910aadSHong Zhang while (nz -- ){ 45178910aadSHong Zhang if (*acol < k) { /* skip lower triangular entries */ 45278910aadSHong Zhang acol++; aval++; 45378910aadSHong Zhang } else { 45478910aadSHong Zhang rtmp[*acol++] = *aval++; 45578910aadSHong Zhang *bval++ = 0.0; /* for in-place factorization */ 45678910aadSHong Zhang } 45778910aadSHong Zhang } 4583cea4cbeSHong Zhang 4593cea4cbeSHong Zhang /* shift the diagonal of the matrix */ 4603cea4cbeSHong Zhang if (sctx.nshift) rtmp[k] += sctx.shift_amount; 46178910aadSHong Zhang 46278910aadSHong Zhang /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 46378910aadSHong Zhang dk = rtmp[k]; 46478910aadSHong Zhang i = jl[k]; /* first row to be added to k_th row */ 46578910aadSHong Zhang 46678910aadSHong Zhang while (i < k){ 46778910aadSHong Zhang nexti = jl[i]; /* next row to be added to k_th row */ 46878910aadSHong Zhang /* compute multiplier, update D(k) and U(i,k) */ 46978910aadSHong Zhang ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 47078910aadSHong Zhang uikdi = - ba[ili]*ba[bi[i]]; 47178910aadSHong Zhang dk += uikdi*ba[ili]; 47278910aadSHong Zhang ba[ili] = uikdi; /* -U(i,k) */ 47378910aadSHong Zhang 47478910aadSHong Zhang /* add multiple of row i to k-th row ... */ 47578910aadSHong Zhang jmin = ili + 1; 47678910aadSHong Zhang nz = bi[i+1] - jmin; 47778910aadSHong Zhang if (nz > 0){ 47878910aadSHong Zhang bcol = bj + jmin; 47978910aadSHong Zhang bval = ba + jmin; 48078910aadSHong Zhang while (nz --) rtmp[*bcol++] += uikdi*(*bval++); 48178910aadSHong Zhang /* update il and jl for i-th row */ 48278910aadSHong Zhang il[i] = jmin; 48378910aadSHong Zhang j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; 48478910aadSHong Zhang } 48578910aadSHong Zhang i = nexti; 48678910aadSHong Zhang } 48778910aadSHong Zhang 4883cea4cbeSHong Zhang /* shift the diagonals when zero pivot is detected */ 4893cea4cbeSHong Zhang /* compute rs=sum of abs(off-diagonal) */ 4903cea4cbeSHong Zhang rs = 0.0; 49178910aadSHong Zhang jmin = bi[k]+1; 49278910aadSHong Zhang nz = bi[k+1] - jmin; 49378910aadSHong Zhang if (nz){ 49478910aadSHong Zhang bcol = bj + jmin; 49578910aadSHong Zhang while (nz--){ 496*89f845c8SHong Zhang rs += PetscAbsScalar(rtmp[*bcol]); 497*89f845c8SHong Zhang bcol++; 49878910aadSHong Zhang } 49978910aadSHong Zhang } 5003cea4cbeSHong Zhang 5013cea4cbeSHong Zhang sctx.rs = rs; 5023cea4cbeSHong Zhang sctx.pv = dk; 5033cea4cbeSHong Zhang ierr = Mat_CholeskyCheckShift(info,&sctx,&newshift);CHKERRQ(ierr); 5043cea4cbeSHong Zhang if (newshift == 1){ 5053cea4cbeSHong Zhang break; /* sctx.shift_amount is updated */ 5063cea4cbeSHong Zhang } else if (newshift == -1){ 5073cea4cbeSHong Zhang SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g * rs %g",k,PetscAbsScalar(dk),zeropivot,rs); 50878910aadSHong Zhang } 50978910aadSHong Zhang 51078910aadSHong Zhang /* copy data into U(k,:) */ 51178910aadSHong Zhang ba[bi[k]] = 1.0/dk; 51278910aadSHong Zhang jmin = bi[k]+1; 51378910aadSHong Zhang nz = bi[k+1] - jmin; 51478910aadSHong Zhang if (nz){ 51578910aadSHong Zhang bcol = bj + jmin; 51678910aadSHong Zhang bval = ba + jmin; 51778910aadSHong Zhang while (nz--){ 51878910aadSHong Zhang *bval++ = rtmp[*bcol]; 51978910aadSHong Zhang rtmp[*bcol++] = 0.0; 52078910aadSHong Zhang } 52178910aadSHong Zhang /* add k-th row into il and jl */ 52278910aadSHong Zhang il[k] = jmin; 52378910aadSHong Zhang i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; 52478910aadSHong Zhang } 52578910aadSHong Zhang } 5263cea4cbeSHong Zhang } while (sctx.chshift); 52778910aadSHong Zhang ierr = PetscFree(il);CHKERRQ(ierr); 52878910aadSHong Zhang 52978910aadSHong Zhang C->factor = FACTOR_CHOLESKY; 53078910aadSHong Zhang C->assembled = PETSC_TRUE; 53178910aadSHong Zhang C->preallocated = PETSC_TRUE; 53278910aadSHong Zhang PetscLogFlops(C->m); 5333cea4cbeSHong Zhang if (sctx.nshift){ 5343cea4cbeSHong Zhang if (shiftnz) { 5353cea4cbeSHong Zhang PetscLogInfo(0,"MatCholeskyFactorNumeric_SeqBAIJ_1_NaturalOrdering: number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount); 5363cea4cbeSHong Zhang } else if (shiftpd) { 5373cea4cbeSHong Zhang PetscLogInfo(0,"MatCholeskyFactorNumeric_SeqBAIJ_1_NaturalOrdering: number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount); 53878910aadSHong Zhang } 53978910aadSHong Zhang } 540c05c3958SHong Zhang PetscFunctionReturn(0); 541c05c3958SHong Zhang } 542c05c3958SHong Zhang 543c05c3958SHong Zhang #include "petscbt.h" 544c05c3958SHong Zhang #include "src/mat/utils/freespace.h" 545c05c3958SHong Zhang #undef __FUNCT__ 546c05c3958SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqBAIJ" 547c05c3958SHong Zhang PetscErrorCode MatICCFactorSymbolic_SeqBAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact) 548c05c3958SHong Zhang { 54978910aadSHong Zhang Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 55078910aadSHong Zhang Mat_SeqSBAIJ *b; 55178910aadSHong Zhang Mat B; 55278910aadSHong Zhang PetscErrorCode ierr; 55378910aadSHong Zhang PetscTruth perm_identity; 55478910aadSHong Zhang PetscInt reallocs=0,*rip,i,*ai=a->i,*aj=a->j,am=a->mbs,bs=A->bs,*ui; 55578910aadSHong Zhang PetscInt jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow; 55678910aadSHong Zhang PetscInt nlnk,*lnk,*lnk_lvl,ncols,ncols_upper,*cols,*cols_lvl,*uj,**uj_ptr,**uj_lvl_ptr; 55778910aadSHong Zhang PetscReal fill=info->fill,levels=info->levels; 55878910aadSHong Zhang FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 55978910aadSHong Zhang FreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL; 56078910aadSHong Zhang PetscBT lnkbt; 56178910aadSHong Zhang 562c05c3958SHong Zhang PetscFunctionBegin; 5636ad2eaddSHong Zhang if (bs > 1){ 5646ad2eaddSHong Zhang if (!a->sbaijMat){ 5656ad2eaddSHong Zhang ierr = MatConvert(A,MATSEQSBAIJ,&a->sbaijMat);CHKERRQ(ierr); 5666ad2eaddSHong Zhang } 5676ad2eaddSHong Zhang ierr = MatICCFactorSymbolic(a->sbaijMat,perm,info,fact);CHKERRQ(ierr); 5686ad2eaddSHong Zhang B = *fact; 5696ad2eaddSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N; 5706ad2eaddSHong Zhang PetscFunctionReturn(0); 5716ad2eaddSHong Zhang } 5726ad2eaddSHong Zhang 57378910aadSHong Zhang ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 57478910aadSHong Zhang ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 57578910aadSHong Zhang 57678910aadSHong Zhang /* special case that simply copies fill pattern */ 57778910aadSHong Zhang if (!levels && perm_identity) { 57878910aadSHong Zhang ierr = MatMarkDiagonal_SeqBAIJ(A);CHKERRQ(ierr); 57978910aadSHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 58078910aadSHong Zhang for (i=0; i<am; i++) { 58178910aadSHong Zhang ui[i] = ai[i+1] - a->diag[i]; /* ui: rowlengths - changes when !perm_identity */ 58278910aadSHong Zhang } 58378910aadSHong Zhang ierr = MatCreate(PETSC_COMM_SELF,am,am,am,am,fact);CHKERRQ(ierr); 58478910aadSHong Zhang B = *fact; 58578910aadSHong Zhang ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr); 58678910aadSHong Zhang ierr = MatSeqSBAIJSetPreallocation(B,1,0,ui);CHKERRQ(ierr); 58778910aadSHong Zhang 58878910aadSHong Zhang b = (Mat_SeqSBAIJ*)B->data; 58978910aadSHong Zhang uj = b->j; 59078910aadSHong Zhang for (i=0; i<am; i++) { 59178910aadSHong Zhang aj = a->j + a->diag[i]; 59278910aadSHong Zhang for (j=0; j<ui[i]; j++){ 59378910aadSHong Zhang *uj++ = *aj++; 59478910aadSHong Zhang } 59578910aadSHong Zhang b->ilen[i] = ui[i]; 59678910aadSHong Zhang } 59778910aadSHong Zhang ierr = PetscFree(ui);CHKERRQ(ierr); 59878910aadSHong Zhang ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 59978910aadSHong Zhang ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 60078910aadSHong Zhang 60178910aadSHong Zhang B->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 60278910aadSHong Zhang B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 60378910aadSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering; 60478910aadSHong Zhang PetscFunctionReturn(0); 60578910aadSHong Zhang } 60678910aadSHong Zhang 60778910aadSHong Zhang /* initialization */ 60878910aadSHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 60978910aadSHong Zhang ui[0] = 0; 61078910aadSHong Zhang ierr = PetscMalloc((2*am+1)*sizeof(PetscInt),&cols_lvl);CHKERRQ(ierr); 61178910aadSHong Zhang 61278910aadSHong Zhang /* jl: linked list for storing indices of the pivot rows 61378910aadSHong Zhang il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */ 61478910aadSHong Zhang ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt*),&jl);CHKERRQ(ierr); 61578910aadSHong Zhang il = jl + am; 61678910aadSHong Zhang uj_ptr = (PetscInt**)(il + am); 61778910aadSHong Zhang uj_lvl_ptr = (PetscInt**)(uj_ptr + am); 61878910aadSHong Zhang for (i=0; i<am; i++){ 61978910aadSHong Zhang jl[i] = am; il[i] = 0; 62078910aadSHong Zhang } 62178910aadSHong Zhang 62278910aadSHong Zhang /* create and initialize a linked list for storing column indices of the active row k */ 62378910aadSHong Zhang nlnk = am + 1; 62478910aadSHong Zhang ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 62578910aadSHong Zhang 62678910aadSHong Zhang /* initial FreeSpace size is fill*(ai[am]+1) */ 62778910aadSHong Zhang ierr = GetMoreSpace((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr); 62878910aadSHong Zhang current_space = free_space; 62978910aadSHong Zhang ierr = GetMoreSpace((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr); 63078910aadSHong Zhang current_space_lvl = free_space_lvl; 63178910aadSHong Zhang 63278910aadSHong Zhang for (k=0; k<am; k++){ /* for each active row k */ 63378910aadSHong Zhang /* initialize lnk by the column indices of row rip[k] of A */ 63478910aadSHong Zhang nzk = 0; 63578910aadSHong Zhang ncols = ai[rip[k]+1] - ai[rip[k]]; 63678910aadSHong Zhang ncols_upper = 0; 63778910aadSHong Zhang cols = cols_lvl + am; 63878910aadSHong Zhang for (j=0; j<ncols; j++){ 63978910aadSHong Zhang i = rip[*(aj + ai[rip[k]] + j)]; 64078910aadSHong Zhang if (i >= k){ /* only take upper triangular entry */ 64178910aadSHong Zhang cols[ncols_upper] = i; 64278910aadSHong Zhang cols_lvl[ncols_upper] = -1; /* initialize level for nonzero entries */ 64378910aadSHong Zhang ncols_upper++; 64478910aadSHong Zhang } 64578910aadSHong Zhang } 64678910aadSHong Zhang ierr = PetscIncompleteLLAdd(ncols_upper,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 64778910aadSHong Zhang nzk += nlnk; 64878910aadSHong Zhang 64978910aadSHong Zhang /* update lnk by computing fill-in for each pivot row to be merged in */ 65078910aadSHong Zhang prow = jl[k]; /* 1st pivot row */ 65178910aadSHong Zhang 65278910aadSHong Zhang while (prow < k){ 65378910aadSHong Zhang nextprow = jl[prow]; 65478910aadSHong Zhang 65578910aadSHong Zhang /* merge prow into k-th row */ 65678910aadSHong Zhang jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:am-1) */ 65778910aadSHong Zhang jmax = ui[prow+1]; 65878910aadSHong Zhang ncols = jmax-jmin; 65978910aadSHong Zhang i = jmin - ui[prow]; 66078910aadSHong Zhang cols = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */ 66178910aadSHong Zhang for (j=0; j<ncols; j++) cols_lvl[j] = *(uj_lvl_ptr[prow] + i + j); 66278910aadSHong Zhang ierr = PetscIncompleteLLAdd(ncols,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 66378910aadSHong Zhang nzk += nlnk; 66478910aadSHong Zhang 66578910aadSHong Zhang /* update il and jl for prow */ 66678910aadSHong Zhang if (jmin < jmax){ 66778910aadSHong Zhang il[prow] = jmin; 66878910aadSHong Zhang j = *cols; jl[prow] = jl[j]; jl[j] = prow; 66978910aadSHong Zhang } 67078910aadSHong Zhang prow = nextprow; 67178910aadSHong Zhang } 67278910aadSHong Zhang 67378910aadSHong Zhang /* if free space is not available, make more free space */ 67478910aadSHong Zhang if (current_space->local_remaining<nzk) { 67578910aadSHong Zhang i = am - k + 1; /* num of unfactored rows */ 67678910aadSHong Zhang i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 67778910aadSHong Zhang ierr = GetMoreSpace(i,¤t_space);CHKERRQ(ierr); 67878910aadSHong Zhang ierr = GetMoreSpace(i,¤t_space_lvl);CHKERRQ(ierr); 67978910aadSHong Zhang reallocs++; 68078910aadSHong Zhang } 68178910aadSHong Zhang 68278910aadSHong Zhang /* copy data into free_space and free_space_lvl, then initialize lnk */ 68378910aadSHong Zhang ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr); 68478910aadSHong Zhang 68578910aadSHong Zhang /* add the k-th row into il and jl */ 68678910aadSHong Zhang if (nzk-1 > 0){ 68778910aadSHong Zhang i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */ 68878910aadSHong Zhang jl[k] = jl[i]; jl[i] = k; 68978910aadSHong Zhang il[k] = ui[k] + 1; 69078910aadSHong Zhang } 69178910aadSHong Zhang uj_ptr[k] = current_space->array; 69278910aadSHong Zhang uj_lvl_ptr[k] = current_space_lvl->array; 69378910aadSHong Zhang 69478910aadSHong Zhang current_space->array += nzk; 69578910aadSHong Zhang current_space->local_used += nzk; 69678910aadSHong Zhang current_space->local_remaining -= nzk; 69778910aadSHong Zhang 69878910aadSHong Zhang current_space_lvl->array += nzk; 69978910aadSHong Zhang current_space_lvl->local_used += nzk; 70078910aadSHong Zhang current_space_lvl->local_remaining -= nzk; 70178910aadSHong Zhang 70278910aadSHong Zhang ui[k+1] = ui[k] + nzk; 70378910aadSHong Zhang } 70478910aadSHong Zhang 70578910aadSHong Zhang if (ai[am] != 0) { 70678910aadSHong Zhang PetscReal af = ((PetscReal)(2*ui[am]-am))/((PetscReal)ai[am]); 70778910aadSHong Zhang PetscLogInfo(A,"MatICCFactorSymbolic_SeqBAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,fill,af); 70878910aadSHong Zhang PetscLogInfo(A,"MatICCFactorSymbolic_SeqBAIJ:Run with -pc_cholesky_fill %g or use \n",af); 70978910aadSHong Zhang PetscLogInfo(A,"MatICCFactorSymbolic_SeqBAIJ:PCCholeskySetFill(pc,%g) for best performance.\n",af); 71078910aadSHong Zhang } else { 71178910aadSHong Zhang PetscLogInfo(A,"MatICCFactorSymbolic_SeqBAIJ:Empty matrix.\n"); 71278910aadSHong Zhang } 71378910aadSHong Zhang 71478910aadSHong Zhang ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 71578910aadSHong Zhang ierr = PetscFree(jl);CHKERRQ(ierr); 71678910aadSHong Zhang ierr = PetscFree(cols_lvl);CHKERRQ(ierr); 71778910aadSHong Zhang 71878910aadSHong Zhang /* destroy list of free space and other temporary array(s) */ 71978910aadSHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 72078910aadSHong Zhang ierr = MakeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 72178910aadSHong Zhang ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 72278910aadSHong Zhang ierr = DestroySpace(free_space_lvl);CHKERRQ(ierr); 72378910aadSHong Zhang 72478910aadSHong Zhang /* put together the new matrix in MATSEQSBAIJ format */ 72578910aadSHong Zhang ierr = MatCreate(PETSC_COMM_SELF,am,am,am,am,fact);CHKERRQ(ierr); 72678910aadSHong Zhang B = *fact; 72778910aadSHong Zhang ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr); 72878910aadSHong Zhang ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr); 72978910aadSHong Zhang 73078910aadSHong Zhang b = (Mat_SeqSBAIJ*)B->data; 73178910aadSHong Zhang ierr = PetscFree(b->imax);CHKERRQ(ierr); 73278910aadSHong Zhang b->singlemalloc = PETSC_FALSE; 73378910aadSHong Zhang /* the next line frees the default space generated by the Create() */ 73478910aadSHong Zhang ierr = PetscFree(b->a);CHKERRQ(ierr); 73578910aadSHong Zhang ierr = PetscFree(b->ilen);CHKERRQ(ierr); 73678910aadSHong Zhang ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 73778910aadSHong Zhang b->j = uj; 73878910aadSHong Zhang b->i = ui; 73978910aadSHong Zhang b->diag = 0; 74078910aadSHong Zhang b->ilen = 0; 74178910aadSHong Zhang b->imax = 0; 74278910aadSHong Zhang b->row = perm; 74378910aadSHong Zhang b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 74478910aadSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 74578910aadSHong Zhang b->icol = perm; 74678910aadSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 74778910aadSHong Zhang ierr = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 74878910aadSHong Zhang PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar))); 74978910aadSHong Zhang b->maxnz = b->nz = ui[am]; 75078910aadSHong Zhang 75178910aadSHong Zhang B->factor = FACTOR_CHOLESKY; 75278910aadSHong Zhang B->info.factor_mallocs = reallocs; 75378910aadSHong Zhang B->info.fill_ratio_given = fill; 75478910aadSHong Zhang if (ai[am] != 0) { 75578910aadSHong Zhang B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]); 75678910aadSHong Zhang } else { 75778910aadSHong Zhang B->info.fill_ratio_needed = 0.0; 75878910aadSHong Zhang } 75978910aadSHong Zhang if (perm_identity){ 76078910aadSHong Zhang B->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 76178910aadSHong Zhang B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 76278910aadSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering; 76378910aadSHong Zhang } else { 76478910aadSHong Zhang (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N; 76578910aadSHong Zhang } 766c05c3958SHong Zhang PetscFunctionReturn(0); 767c05c3958SHong Zhang } 768c05c3958SHong Zhang 769c05c3958SHong Zhang #undef __FUNCT__ 770c05c3958SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqBAIJ" 771c05c3958SHong Zhang PetscErrorCode MatCholeskyFactorSymbolic_SeqBAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact) 772c05c3958SHong Zhang { 77378910aadSHong Zhang Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 77478910aadSHong Zhang Mat_SeqSBAIJ *b; 77578910aadSHong Zhang Mat B; 77678910aadSHong Zhang PetscErrorCode ierr; 77778910aadSHong Zhang PetscTruth perm_identity; 77878910aadSHong Zhang PetscReal fill = info->fill; 77978910aadSHong Zhang PetscInt *rip,*riip,i,mbs=a->mbs,bs=A->bs,*ai=a->i,*aj=a->j,reallocs=0,prow; 78078910aadSHong Zhang PetscInt *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow; 78178910aadSHong Zhang PetscInt nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr; 78278910aadSHong Zhang FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 78378910aadSHong Zhang PetscBT lnkbt; 78478910aadSHong Zhang IS iperm; 78578910aadSHong Zhang 786c05c3958SHong Zhang PetscFunctionBegin; 7876ad2eaddSHong Zhang if (bs > 1) { /* convert to seqsbaij */ 7886ad2eaddSHong Zhang if (!a->sbaijMat){ 7896ad2eaddSHong Zhang ierr = MatConvert(A,MATSEQSBAIJ,&a->sbaijMat);CHKERRQ(ierr); 7906ad2eaddSHong Zhang } 7916ad2eaddSHong Zhang ierr = MatCholeskyFactorSymbolic(a->sbaijMat,perm,info,fact);CHKERRQ(ierr); 7926ad2eaddSHong Zhang B = *fact; 7936ad2eaddSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N; 7946ad2eaddSHong Zhang PetscFunctionReturn(0); 7956ad2eaddSHong Zhang } 7966ad2eaddSHong Zhang 79778910aadSHong Zhang /* check whether perm is the identity mapping */ 79878910aadSHong Zhang ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 79978910aadSHong Zhang ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 80078910aadSHong Zhang 80178910aadSHong Zhang if (!perm_identity){ 80278910aadSHong Zhang /* check if perm is symmetric! */ 80378910aadSHong Zhang ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr); 80478910aadSHong Zhang ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr); 80578910aadSHong Zhang for (i=0; i<mbs; i++) { 80678910aadSHong Zhang if (rip[i] != riip[i]) SETERRQ(PETSC_ERR_ARG_INCOMP,"Non-symmetric permutation, must use symmetric permutation"); 80778910aadSHong Zhang } 80878910aadSHong Zhang ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr); 80978910aadSHong Zhang ierr = ISDestroy(iperm);CHKERRQ(ierr); 81078910aadSHong Zhang } 81178910aadSHong Zhang 81278910aadSHong Zhang /* initialization */ 81378910aadSHong Zhang ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 81478910aadSHong Zhang ui[0] = 0; 81578910aadSHong Zhang 81678910aadSHong Zhang /* jl: linked list for storing indices of the pivot rows 81778910aadSHong Zhang il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */ 81878910aadSHong Zhang ierr = PetscMalloc((3*mbs+1)*sizeof(PetscInt)+mbs*sizeof(PetscInt*),&jl);CHKERRQ(ierr); 81978910aadSHong Zhang il = jl + mbs; 82078910aadSHong Zhang cols = il + mbs; 82178910aadSHong Zhang ui_ptr = (PetscInt**)(cols + mbs); 82278910aadSHong Zhang for (i=0; i<mbs; i++){ 82378910aadSHong Zhang jl[i] = mbs; il[i] = 0; 82478910aadSHong Zhang } 82578910aadSHong Zhang 82678910aadSHong Zhang /* create and initialize a linked list for storing column indices of the active row k */ 82778910aadSHong Zhang nlnk = mbs + 1; 82878910aadSHong Zhang ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 82978910aadSHong Zhang 83078910aadSHong Zhang /* initial FreeSpace size is fill*(ai[mbs]+1) */ 83178910aadSHong Zhang ierr = GetMoreSpace((PetscInt)(fill*(ai[mbs]+1)),&free_space);CHKERRQ(ierr); 83278910aadSHong Zhang current_space = free_space; 83378910aadSHong Zhang 83478910aadSHong Zhang for (k=0; k<mbs; k++){ /* for each active row k */ 83578910aadSHong Zhang /* initialize lnk by the column indices of row rip[k] of A */ 83678910aadSHong Zhang nzk = 0; 83778910aadSHong Zhang ncols = ai[rip[k]+1] - ai[rip[k]]; 83878910aadSHong Zhang ncols_upper = 0; 83978910aadSHong Zhang for (j=0; j<ncols; j++){ 84078910aadSHong Zhang i = rip[*(aj + ai[rip[k]] + j)]; 84178910aadSHong Zhang if (i >= k){ /* only take upper triangular entry */ 84278910aadSHong Zhang cols[ncols_upper] = i; 84378910aadSHong Zhang ncols_upper++; 84478910aadSHong Zhang } 84578910aadSHong Zhang } 84678910aadSHong Zhang ierr = PetscLLAdd(ncols_upper,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 84778910aadSHong Zhang nzk += nlnk; 84878910aadSHong Zhang 84978910aadSHong Zhang /* update lnk by computing fill-in for each pivot row to be merged in */ 85078910aadSHong Zhang prow = jl[k]; /* 1st pivot row */ 85178910aadSHong Zhang 85278910aadSHong Zhang while (prow < k){ 85378910aadSHong Zhang nextprow = jl[prow]; 85478910aadSHong Zhang /* merge prow into k-th row */ 85578910aadSHong Zhang jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:mbs-1) */ 85678910aadSHong Zhang jmax = ui[prow+1]; 85778910aadSHong Zhang ncols = jmax-jmin; 85878910aadSHong Zhang uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */ 85978910aadSHong Zhang ierr = PetscLLAdd(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 86078910aadSHong Zhang nzk += nlnk; 86178910aadSHong Zhang 86278910aadSHong Zhang /* update il and jl for prow */ 86378910aadSHong Zhang if (jmin < jmax){ 86478910aadSHong Zhang il[prow] = jmin; 86578910aadSHong Zhang j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow; 86678910aadSHong Zhang } 86778910aadSHong Zhang prow = nextprow; 86878910aadSHong Zhang } 86978910aadSHong Zhang 87078910aadSHong Zhang /* if free space is not available, make more free space */ 87178910aadSHong Zhang if (current_space->local_remaining<nzk) { 87278910aadSHong Zhang i = mbs - k + 1; /* num of unfactored rows */ 87378910aadSHong Zhang i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 87478910aadSHong Zhang ierr = GetMoreSpace(i,¤t_space);CHKERRQ(ierr); 87578910aadSHong Zhang reallocs++; 87678910aadSHong Zhang } 87778910aadSHong Zhang 87878910aadSHong Zhang /* copy data into free space, then initialize lnk */ 87978910aadSHong Zhang ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 88078910aadSHong Zhang 88178910aadSHong Zhang /* add the k-th row into il and jl */ 88278910aadSHong Zhang if (nzk-1 > 0){ 88378910aadSHong Zhang i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */ 88478910aadSHong Zhang jl[k] = jl[i]; jl[i] = k; 88578910aadSHong Zhang il[k] = ui[k] + 1; 88678910aadSHong Zhang } 88778910aadSHong Zhang ui_ptr[k] = current_space->array; 88878910aadSHong Zhang current_space->array += nzk; 88978910aadSHong Zhang current_space->local_used += nzk; 89078910aadSHong Zhang current_space->local_remaining -= nzk; 89178910aadSHong Zhang 89278910aadSHong Zhang ui[k+1] = ui[k] + nzk; 89378910aadSHong Zhang } 89478910aadSHong Zhang 89578910aadSHong Zhang if (ai[mbs] != 0) { 89678910aadSHong Zhang PetscReal af = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]); 89778910aadSHong Zhang PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,fill,af); 89878910aadSHong Zhang PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Run with -pc_cholesky_fill %g or use \n",af); 89978910aadSHong Zhang PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:PCCholeskySetFill(pc,%g) for best performance.\n",af); 90078910aadSHong Zhang } else { 90178910aadSHong Zhang PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Empty matrix.\n"); 90278910aadSHong Zhang } 90378910aadSHong Zhang 90478910aadSHong Zhang ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 90578910aadSHong Zhang ierr = PetscFree(jl);CHKERRQ(ierr); 90678910aadSHong Zhang 90778910aadSHong Zhang /* destroy list of free space and other temporary array(s) */ 90878910aadSHong Zhang ierr = PetscMalloc((ui[mbs]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 90978910aadSHong Zhang ierr = MakeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 91078910aadSHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 91178910aadSHong Zhang 91278910aadSHong Zhang /* put together the new matrix in MATSEQSBAIJ format */ 91378910aadSHong Zhang ierr = MatCreate(PETSC_COMM_SELF,mbs,mbs,mbs,mbs,fact);CHKERRQ(ierr); 91478910aadSHong Zhang B = *fact; 91578910aadSHong Zhang ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr); 91678910aadSHong Zhang ierr = MatSeqSBAIJSetPreallocation(B,bs,0,PETSC_NULL);CHKERRQ(ierr); 91778910aadSHong Zhang 91878910aadSHong Zhang b = (Mat_SeqSBAIJ*)B->data; 91978910aadSHong Zhang ierr = PetscFree(b->imax);CHKERRQ(ierr); 92078910aadSHong Zhang b->singlemalloc = PETSC_FALSE; 92178910aadSHong Zhang /* the next line frees the default space generated by the Create() */ 92278910aadSHong Zhang ierr = PetscFree(b->a);CHKERRQ(ierr); 92378910aadSHong Zhang ierr = PetscFree(b->ilen);CHKERRQ(ierr); 92478910aadSHong Zhang ierr = PetscMalloc((ui[mbs]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 92578910aadSHong Zhang b->j = uj; 92678910aadSHong Zhang b->i = ui; 92778910aadSHong Zhang b->diag = 0; 92878910aadSHong Zhang b->ilen = 0; 92978910aadSHong Zhang b->imax = 0; 93078910aadSHong Zhang b->row = perm; 93178910aadSHong Zhang b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 93278910aadSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 93378910aadSHong Zhang b->icol = perm; 93478910aadSHong Zhang ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 93578910aadSHong Zhang ierr = PetscMalloc((mbs+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 93678910aadSHong Zhang PetscLogObjectMemory(B,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar))); 93778910aadSHong Zhang b->maxnz = b->nz = ui[mbs]; 93878910aadSHong Zhang 93978910aadSHong Zhang B->factor = FACTOR_CHOLESKY; 94078910aadSHong Zhang B->info.factor_mallocs = reallocs; 94178910aadSHong Zhang B->info.fill_ratio_given = fill; 94278910aadSHong Zhang if (ai[mbs] != 0) { 94378910aadSHong Zhang B->info.fill_ratio_needed = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]); 94478910aadSHong Zhang } else { 94578910aadSHong Zhang B->info.fill_ratio_needed = 0.0; 94678910aadSHong Zhang } 94778910aadSHong Zhang if (perm_identity){ 9486ad2eaddSHong Zhang B->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 9496ad2eaddSHong Zhang B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 9506ad2eaddSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering; 95178910aadSHong Zhang } else { 9526ad2eaddSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N; 95378910aadSHong Zhang } 954c05c3958SHong Zhang PetscFunctionReturn(0); 955c05c3958SHong Zhang } 956