xref: /petsc/src/mat/impls/baij/seq/baijfact.c (revision 62bba022857b4b9c504bd48e984f55e96e1c91e3)
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"
15af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2(Mat A,MatFactorInfo *info,Mat *B)
164eeb42bcSBarry Smith {
174eeb42bcSBarry Smith   Mat            C = *B;
184eeb42bcSBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data;
197cdf2dbbSSatish Balay   IS             isrow = b->row,isicol = b->icol;
206849ba73SBarry Smith   PetscErrorCode ierr;
21690b6cddSBarry Smith   PetscInt       *r,*ic,i,j,n = a->mbs,*bi = b->i,*bj = b->j;
22690b6cddSBarry Smith   PetscInt       *ajtmpold,*ajtmp,nz,row;
23690b6cddSBarry Smith   PetscInt       *diag_offset=b->diag,idx,*ai=a->i,*aj=a->j,*pj;
24329f5518SBarry Smith   MatScalar      *pv,*v,*rtmp,m1,m2,m3,m4,*pc,*w,*x,x1,x2,x3,x4;
252515f8d2SSatish Balay   MatScalar      p1,p2,p3,p4;
263f1db9ecSBarry Smith   MatScalar      *ba = b->a,*aa = a->a;
27*62bba022SBarry 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];
88*62bba022SBarry 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);
944eeb42bcSBarry Smith   C->factor = FACTOR_LU;
954eeb42bcSBarry Smith   C->assembled = PETSC_TRUE;
96efee365bSSatish Balay   ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
973a40ed3dSBarry Smith   PetscFunctionReturn(0);
984eeb42bcSBarry Smith }
994cd38560SBarry Smith /*
1004cd38560SBarry Smith       Version for when blocks are 2 by 2 Using natural ordering
1014cd38560SBarry Smith */
1024a2ae208SSatish Balay #undef __FUNCT__
1034a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering"
104af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering(Mat A,MatFactorInfo *info,Mat *B)
1054cd38560SBarry Smith {
1064cd38560SBarry Smith   Mat            C = *B;
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;
115*62bba022SBarry 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];
174*62bba022SBarry Smith     ierr = Kernel_A_gets_inverse_A_2(w,shift);CHKERRQ(ierr);
1754cd38560SBarry Smith   }
1764cd38560SBarry Smith 
177606d414cSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
1784cd38560SBarry Smith   C->factor    = FACTOR_LU;
1794cd38560SBarry Smith   C->assembled = PETSC_TRUE;
180efee365bSSatish Balay   ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
1814cd38560SBarry Smith   PetscFunctionReturn(0);
1824cd38560SBarry Smith }
1837fc0212eSBarry Smith 
1847fc0212eSBarry Smith /* ----------------------------------------------------------- */
1857fc0212eSBarry Smith /*
1867fc0212eSBarry Smith      Version for when blocks are 1 by 1.
1877fc0212eSBarry Smith */
1884a2ae208SSatish Balay #undef __FUNCT__
1894a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_1"
190af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqBAIJ_1(Mat A,MatFactorInfo *info,Mat *B)
1917fc0212eSBarry Smith {
1927fc0212eSBarry Smith   Mat            C = *B;
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);
2487fc0212eSBarry Smith   C->factor    = FACTOR_LU;
2497fc0212eSBarry Smith   C->assembled = PETSC_TRUE;
250899cda47SBarry Smith   ierr = PetscLogFlops(C->cmap.n);CHKERRQ(ierr);
2513a40ed3dSBarry Smith   PetscFunctionReturn(0);
2527fc0212eSBarry Smith }
2537fc0212eSBarry Smith 
254273d9f13SBarry Smith 
2555d517e7eSBarry Smith /* ----------------------------------------------------------- */
2564a2ae208SSatish Balay #undef __FUNCT__
2574a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqBAIJ"
258dfbe8321SBarry Smith PetscErrorCode MatLUFactor_SeqBAIJ(Mat A,IS row,IS col,MatFactorInfo *info)
2595d517e7eSBarry Smith {
260dfbe8321SBarry Smith   PetscErrorCode ierr;
2615d517e7eSBarry Smith   Mat            C;
2625d517e7eSBarry Smith 
2633a40ed3dSBarry Smith   PetscFunctionBegin;
264b9b97703SBarry Smith   ierr = MatLUFactorSymbolic(A,row,col,info,&C);CHKERRQ(ierr);
265af281ebdSHong Zhang   ierr = MatLUFactorNumeric(A,info,&C);CHKERRQ(ierr);
266273d9f13SBarry Smith   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
26752e6d16bSBarry Smith   ierr = PetscLogObjectParent(A,((Mat_SeqBAIJ*)(A->data))->icol);CHKERRQ(ierr);
2683a40ed3dSBarry Smith   PetscFunctionReturn(0);
2695d517e7eSBarry Smith }
2704cd38560SBarry Smith 
271c05c3958SHong Zhang #include "src/mat/impls/sbaij/seq/sbaij.h"
272c05c3958SHong Zhang #undef __FUNCT__
273c05c3958SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N"
274af281ebdSHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N(Mat A,MatFactorInfo *info,Mat *B)
275c05c3958SHong Zhang {
276f3086b4bSHong Zhang   PetscErrorCode ierr;
27778910aadSHong Zhang   Mat            C = *B;
27878910aadSHong Zhang   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data;
27978910aadSHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
28078910aadSHong Zhang   IS             ip=b->row;
281899cda47SBarry Smith   PetscInt       *rip,i,j,mbs=a->mbs,bs=A->rmap.bs,*bi=b->i,*bj=b->j,*bcol;
28278910aadSHong Zhang   PetscInt       *ai=a->i,*aj=a->j;
28378910aadSHong Zhang   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
28478910aadSHong Zhang   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
2853cea4cbeSHong Zhang   PetscReal      zeropivot,rs,shiftnz;
286fbf22428SSatish Balay   PetscReal      shiftpd;
2873cea4cbeSHong Zhang   ChShift_Ctx    sctx;
2883cea4cbeSHong Zhang   PetscInt       newshift;
28978910aadSHong Zhang 
290c05c3958SHong Zhang   PetscFunctionBegin;
2916ad2eaddSHong Zhang   if (bs > 1) {
2926ad2eaddSHong Zhang     if (!a->sbaijMat){
293ceb03754SKris Buschelman       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
2946ad2eaddSHong Zhang     }
295f3086b4bSHong Zhang     ierr = (a->sbaijMat)->ops->choleskyfactornumeric(a->sbaijMat,info,B);CHKERRQ(ierr);
2966ad2eaddSHong Zhang     ierr = MatDestroy(a->sbaijMat);CHKERRQ(ierr);
2976ad2eaddSHong Zhang     a->sbaijMat = PETSC_NULL;
2986ad2eaddSHong Zhang     PetscFunctionReturn(0);
2996ad2eaddSHong Zhang   }
30078910aadSHong Zhang 
30178910aadSHong Zhang   /* initialization */
3023cea4cbeSHong Zhang   shiftnz   = info->shiftnz;
3033cea4cbeSHong Zhang   shiftpd   = info->shiftpd;
3043cea4cbeSHong Zhang   zeropivot = info->zeropivot;
3053cea4cbeSHong Zhang 
3066ad2eaddSHong Zhang   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
30778910aadSHong Zhang   nz   = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar);
30878910aadSHong Zhang   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
30978910aadSHong Zhang   jl   = il + mbs;
31078910aadSHong Zhang   rtmp = (MatScalar*)(jl + mbs);
31178910aadSHong Zhang 
3123cea4cbeSHong Zhang   sctx.shift_amount = 0;
3133cea4cbeSHong Zhang   sctx.nshift       = 0;
31478910aadSHong Zhang   do {
3153cea4cbeSHong Zhang     sctx.chshift = PETSC_FALSE;
31678910aadSHong Zhang     for (i=0; i<mbs; i++) {
31778910aadSHong Zhang       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
31878910aadSHong Zhang     }
31978910aadSHong Zhang 
32078910aadSHong Zhang     for (k = 0; k<mbs; k++){
32178910aadSHong Zhang       bval = ba + bi[k];
32278910aadSHong Zhang       /* initialize k-th row by the perm[k]-th row of A */
32378910aadSHong Zhang       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
32478910aadSHong Zhang       for (j = jmin; j < jmax; j++){
32578910aadSHong Zhang         col = rip[aj[j]];
32678910aadSHong Zhang         if (col >= k){ /* only take upper triangular entry */
32778910aadSHong Zhang           rtmp[col] = aa[j];
32878910aadSHong Zhang           *bval++  = 0.0; /* for in-place factorization */
32978910aadSHong Zhang         }
33078910aadSHong Zhang       }
3313cea4cbeSHong Zhang 
3323cea4cbeSHong Zhang       /* shift the diagonal of the matrix */
3333cea4cbeSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
33478910aadSHong Zhang 
33578910aadSHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
33678910aadSHong Zhang       dk = rtmp[k];
33778910aadSHong Zhang       i = jl[k]; /* first row to be added to k_th row  */
33878910aadSHong Zhang 
33978910aadSHong Zhang       while (i < k){
34078910aadSHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
34178910aadSHong Zhang 
34278910aadSHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
34378910aadSHong Zhang         ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
34478910aadSHong Zhang         uikdi = - ba[ili]*ba[bi[i]];  /* diagonal(k) */
34578910aadSHong Zhang         dk += uikdi*ba[ili];
34678910aadSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
34778910aadSHong Zhang 
34878910aadSHong Zhang         /* add multiple of row i to k-th row */
34978910aadSHong Zhang         jmin = ili + 1; jmax = bi[i+1];
35078910aadSHong Zhang         if (jmin < jmax){
35178910aadSHong Zhang           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
35278910aadSHong Zhang           /* update il and jl for row i */
35378910aadSHong Zhang           il[i] = jmin;
35478910aadSHong Zhang           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
35578910aadSHong Zhang         }
35678910aadSHong Zhang         i = nexti;
35778910aadSHong Zhang       }
35878910aadSHong Zhang 
3593cea4cbeSHong Zhang       /* shift the diagonals when zero pivot is detected */
3603cea4cbeSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
3613cea4cbeSHong Zhang       rs   = 0.0;
36278910aadSHong Zhang       jmin = bi[k]+1;
36378910aadSHong Zhang       nz   = bi[k+1] - jmin;
36478910aadSHong Zhang       if (nz){
36578910aadSHong Zhang         bcol = bj + jmin;
36678910aadSHong Zhang         while (nz--){
36789f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
36889f845c8SHong Zhang           bcol++;
36978910aadSHong Zhang         }
37078910aadSHong Zhang       }
3713cea4cbeSHong Zhang 
3723cea4cbeSHong Zhang       sctx.rs = rs;
3733cea4cbeSHong Zhang       sctx.pv = dk;
37445938b79SHong Zhang       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
37545938b79SHong Zhang       if (newshift == 1) break;
37678910aadSHong Zhang 
37778910aadSHong Zhang       /* copy data into U(k,:) */
37878910aadSHong Zhang       ba[bi[k]] = 1.0/dk; /* U(k,k) */
37978910aadSHong Zhang       jmin = bi[k]+1; jmax = bi[k+1];
38078910aadSHong Zhang       if (jmin < jmax) {
38178910aadSHong Zhang         for (j=jmin; j<jmax; j++){
38278910aadSHong Zhang           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
38378910aadSHong Zhang         }
38478910aadSHong Zhang         /* add the k-th row into il and jl */
38578910aadSHong Zhang         il[k] = jmin;
38678910aadSHong Zhang         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
38778910aadSHong Zhang       }
38878910aadSHong Zhang     }
3893cea4cbeSHong Zhang   } while (sctx.chshift);
39078910aadSHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
39178910aadSHong Zhang 
39278910aadSHong Zhang   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
39378910aadSHong Zhang   C->factor       = FACTOR_CHOLESKY;
39478910aadSHong Zhang   C->assembled    = PETSC_TRUE;
39578910aadSHong Zhang   C->preallocated = PETSC_TRUE;
396899cda47SBarry Smith   ierr = PetscLogFlops(C->rmap.N);CHKERRQ(ierr);
3973cea4cbeSHong Zhang   if (sctx.nshift){
3983cea4cbeSHong Zhang     if (shiftnz) {
3991e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
4003cea4cbeSHong Zhang     } else if (shiftpd) {
4011e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
40278910aadSHong Zhang     }
40378910aadSHong Zhang   }
404c05c3958SHong Zhang   PetscFunctionReturn(0);
405c05c3958SHong Zhang }
4064cd38560SBarry Smith 
407c05c3958SHong Zhang #undef __FUNCT__
408c05c3958SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering"
409af281ebdSHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering(Mat A,MatFactorInfo *info,Mat *fact)
410c05c3958SHong Zhang {
41178910aadSHong Zhang   Mat            C = *fact;
41278910aadSHong Zhang   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data;
41378910aadSHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
41478910aadSHong Zhang   PetscErrorCode ierr;
4153cea4cbeSHong Zhang   PetscInt       i,j,am=a->mbs;
41678910aadSHong Zhang   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
4173cea4cbeSHong Zhang   PetscInt       k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz;
41878910aadSHong Zhang   MatScalar      *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval;
4193cea4cbeSHong Zhang   PetscReal      zeropivot,rs,shiftnz;
420fbf22428SSatish Balay   PetscReal      shiftpd;
4213cea4cbeSHong Zhang   ChShift_Ctx    sctx;
4223cea4cbeSHong Zhang   PetscInt       newshift;
42378910aadSHong Zhang 
424c05c3958SHong Zhang   PetscFunctionBegin;
42578910aadSHong Zhang   /* initialization */
4263cea4cbeSHong Zhang   shiftnz   = info->shiftnz;
4273cea4cbeSHong Zhang   shiftpd   = info->shiftpd;
4283cea4cbeSHong Zhang   zeropivot = info->zeropivot;
4293cea4cbeSHong Zhang 
43078910aadSHong Zhang   nz   = (2*am+1)*sizeof(PetscInt)+am*sizeof(MatScalar);
43178910aadSHong Zhang   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
43278910aadSHong Zhang   jl   = il + am;
43378910aadSHong Zhang   rtmp = (MatScalar*)(jl + am);
43478910aadSHong Zhang 
4353cea4cbeSHong Zhang   sctx.shift_amount = 0;
4363cea4cbeSHong Zhang   sctx.nshift       = 0;
43778910aadSHong Zhang   do {
4383cea4cbeSHong Zhang     sctx.chshift = PETSC_FALSE;
43978910aadSHong Zhang     for (i=0; i<am; i++) {
44078910aadSHong Zhang       rtmp[i] = 0.0; jl[i] = am; il[0] = 0;
44178910aadSHong Zhang     }
44278910aadSHong Zhang 
44378910aadSHong Zhang     for (k = 0; k<am; k++){
44478910aadSHong Zhang     /* initialize k-th row with elements nonzero in row perm(k) of A */
44578910aadSHong Zhang       nz   = ai[k+1] - ai[k];
44678910aadSHong Zhang       acol = aj + ai[k];
44778910aadSHong Zhang       aval = aa + ai[k];
44878910aadSHong Zhang       bval = ba + bi[k];
44978910aadSHong Zhang       while (nz -- ){
45078910aadSHong Zhang         if (*acol < k) { /* skip lower triangular entries */
45178910aadSHong Zhang           acol++; aval++;
45278910aadSHong Zhang         } else {
45378910aadSHong Zhang           rtmp[*acol++] = *aval++;
45478910aadSHong Zhang           *bval++       = 0.0; /* for in-place factorization */
45578910aadSHong Zhang         }
45678910aadSHong Zhang       }
4573cea4cbeSHong Zhang 
4583cea4cbeSHong Zhang       /* shift the diagonal of the matrix */
4593cea4cbeSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
46078910aadSHong Zhang 
46178910aadSHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
46278910aadSHong Zhang       dk = rtmp[k];
46378910aadSHong Zhang       i  = jl[k]; /* first row to be added to k_th row  */
46478910aadSHong Zhang 
46578910aadSHong Zhang       while (i < k){
46678910aadSHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
46778910aadSHong Zhang         /* compute multiplier, update D(k) and U(i,k) */
46878910aadSHong Zhang         ili   = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
46978910aadSHong Zhang         uikdi = - ba[ili]*ba[bi[i]];
47078910aadSHong Zhang         dk   += uikdi*ba[ili];
47178910aadSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
47278910aadSHong Zhang 
47378910aadSHong Zhang         /* add multiple of row i to k-th row ... */
47478910aadSHong Zhang         jmin = ili + 1;
47578910aadSHong Zhang         nz   = bi[i+1] - jmin;
47678910aadSHong Zhang         if (nz > 0){
47778910aadSHong Zhang           bcol = bj + jmin;
47878910aadSHong Zhang           bval = ba + jmin;
47978910aadSHong Zhang           while (nz --) rtmp[*bcol++] += uikdi*(*bval++);
48078910aadSHong Zhang           /* update il and jl for i-th row */
48178910aadSHong Zhang           il[i] = jmin;
48278910aadSHong Zhang           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
48378910aadSHong Zhang         }
48478910aadSHong Zhang         i = nexti;
48578910aadSHong Zhang       }
48678910aadSHong Zhang 
4873cea4cbeSHong Zhang       /* shift the diagonals when zero pivot is detected */
4883cea4cbeSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
4893cea4cbeSHong Zhang       rs   = 0.0;
49078910aadSHong Zhang       jmin = bi[k]+1;
49178910aadSHong Zhang       nz   = bi[k+1] - jmin;
49278910aadSHong Zhang       if (nz){
49378910aadSHong Zhang         bcol = bj + jmin;
49478910aadSHong Zhang         while (nz--){
49589f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
49689f845c8SHong Zhang           bcol++;
49778910aadSHong Zhang         }
49878910aadSHong Zhang       }
4993cea4cbeSHong Zhang 
5003cea4cbeSHong Zhang       sctx.rs = rs;
5013cea4cbeSHong Zhang       sctx.pv = dk;
50245938b79SHong Zhang       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
50345938b79SHong Zhang       if (newshift == 1) break;    /* sctx.shift_amount is updated */
50478910aadSHong Zhang 
50578910aadSHong Zhang       /* copy data into U(k,:) */
50678910aadSHong Zhang       ba[bi[k]] = 1.0/dk;
50778910aadSHong Zhang       jmin      = bi[k]+1;
50878910aadSHong Zhang       nz        = bi[k+1] - jmin;
50978910aadSHong Zhang       if (nz){
51078910aadSHong Zhang         bcol = bj + jmin;
51178910aadSHong Zhang         bval = ba + jmin;
51278910aadSHong Zhang         while (nz--){
51378910aadSHong Zhang           *bval++       = rtmp[*bcol];
51478910aadSHong Zhang           rtmp[*bcol++] = 0.0;
51578910aadSHong Zhang         }
51678910aadSHong Zhang         /* add k-th row into il and jl */
51778910aadSHong Zhang         il[k] = jmin;
51878910aadSHong Zhang         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
51978910aadSHong Zhang       }
52078910aadSHong Zhang     }
5213cea4cbeSHong Zhang   } while (sctx.chshift);
52278910aadSHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
52378910aadSHong Zhang 
52478910aadSHong Zhang   C->factor       = FACTOR_CHOLESKY;
52578910aadSHong Zhang   C->assembled    = PETSC_TRUE;
52678910aadSHong Zhang   C->preallocated = PETSC_TRUE;
527899cda47SBarry Smith   ierr = PetscLogFlops(C->rmap.N);CHKERRQ(ierr);
5283cea4cbeSHong Zhang     if (sctx.nshift){
5293cea4cbeSHong Zhang     if (shiftnz) {
5301e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
5313cea4cbeSHong Zhang     } else if (shiftpd) {
5321e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
53378910aadSHong Zhang     }
53478910aadSHong Zhang   }
535c05c3958SHong Zhang   PetscFunctionReturn(0);
536c05c3958SHong Zhang }
537c05c3958SHong Zhang 
538c05c3958SHong Zhang #include "petscbt.h"
539c05c3958SHong Zhang #include "src/mat/utils/freespace.h"
540c05c3958SHong Zhang #undef __FUNCT__
541c05c3958SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqBAIJ"
542c05c3958SHong Zhang PetscErrorCode MatICCFactorSymbolic_SeqBAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
543c05c3958SHong Zhang {
54478910aadSHong Zhang   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data;
54578910aadSHong Zhang   Mat_SeqSBAIJ       *b;
54678910aadSHong Zhang   Mat                B;
54778910aadSHong Zhang   PetscErrorCode     ierr;
54878910aadSHong Zhang   PetscTruth         perm_identity;
549899cda47SBarry Smith   PetscInt           reallocs=0,*rip,i,*ai=a->i,*aj=a->j,am=a->mbs,bs=A->rmap.bs,*ui;
55078910aadSHong Zhang   PetscInt           jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow;
551cfdb8b8aSHong Zhang   PetscInt           nlnk,*lnk,*lnk_lvl=PETSC_NULL,ncols,ncols_upper,*cols,*cols_lvl,*uj,**uj_ptr,**uj_lvl_ptr;
55278910aadSHong Zhang   PetscReal          fill=info->fill,levels=info->levels;
553a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
554a1a86e44SBarry Smith   PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
55578910aadSHong Zhang   PetscBT            lnkbt;
55678910aadSHong Zhang 
557c05c3958SHong Zhang   PetscFunctionBegin;
5586ad2eaddSHong Zhang   if (bs > 1){
5596ad2eaddSHong Zhang     if (!a->sbaijMat){
560ceb03754SKris Buschelman       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
5616ad2eaddSHong Zhang     }
5626ad2eaddSHong Zhang     ierr = MatICCFactorSymbolic(a->sbaijMat,perm,info,fact);CHKERRQ(ierr);
5636ad2eaddSHong Zhang     B = *fact;
5646ad2eaddSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
5656ad2eaddSHong Zhang     PetscFunctionReturn(0);
5666ad2eaddSHong Zhang   }
5676ad2eaddSHong Zhang 
56878910aadSHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
56978910aadSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
57078910aadSHong Zhang 
57178910aadSHong Zhang   /* special case that simply copies fill pattern */
57278910aadSHong Zhang   if (!levels && perm_identity) {
57378910aadSHong Zhang     ierr = MatMarkDiagonal_SeqBAIJ(A);CHKERRQ(ierr);
57478910aadSHong Zhang     ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
57578910aadSHong Zhang     for (i=0; i<am; i++) {
57678910aadSHong Zhang       ui[i] = ai[i+1] - a->diag[i]; /* ui: rowlengths - changes when !perm_identity */
57778910aadSHong Zhang     }
578f69a0ea3SMatthew Knepley     ierr = MatCreate(PETSC_COMM_SELF,fact);CHKERRQ(ierr);
579f69a0ea3SMatthew Knepley     ierr = MatSetSizes(*fact,am,am,am,am);CHKERRQ(ierr);
58078910aadSHong Zhang     B = *fact;
58178910aadSHong Zhang     ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
58278910aadSHong Zhang     ierr = MatSeqSBAIJSetPreallocation(B,1,0,ui);CHKERRQ(ierr);
58378910aadSHong Zhang 
58478910aadSHong Zhang     b  = (Mat_SeqSBAIJ*)B->data;
58578910aadSHong Zhang     uj = b->j;
58678910aadSHong Zhang     for (i=0; i<am; i++) {
58778910aadSHong Zhang       aj = a->j + a->diag[i];
58878910aadSHong Zhang       for (j=0; j<ui[i]; j++){
58978910aadSHong Zhang         *uj++ = *aj++;
59078910aadSHong Zhang       }
59178910aadSHong Zhang       b->ilen[i] = ui[i];
59278910aadSHong Zhang     }
59378910aadSHong Zhang     ierr = PetscFree(ui);CHKERRQ(ierr);
59478910aadSHong Zhang     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
59578910aadSHong Zhang     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
59678910aadSHong Zhang 
59778910aadSHong Zhang     B->ops->solve                 = MatSolve_SeqSBAIJ_1_NaturalOrdering;
59878910aadSHong Zhang     B->ops->solvetranspose        = MatSolve_SeqSBAIJ_1_NaturalOrdering;
59978910aadSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
60078910aadSHong Zhang     PetscFunctionReturn(0);
60178910aadSHong Zhang   }
60278910aadSHong Zhang 
60378910aadSHong Zhang   /* initialization */
60478910aadSHong Zhang   ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
60578910aadSHong Zhang   ui[0] = 0;
60678910aadSHong Zhang   ierr  = PetscMalloc((2*am+1)*sizeof(PetscInt),&cols_lvl);CHKERRQ(ierr);
60778910aadSHong Zhang 
60878910aadSHong Zhang   /* jl: linked list for storing indices of the pivot rows
60978910aadSHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
61078910aadSHong Zhang   ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt*),&jl);CHKERRQ(ierr);
61178910aadSHong Zhang   il         = jl + am;
61278910aadSHong Zhang   uj_ptr     = (PetscInt**)(il + am);
61378910aadSHong Zhang   uj_lvl_ptr = (PetscInt**)(uj_ptr + am);
61478910aadSHong Zhang   for (i=0; i<am; i++){
61578910aadSHong Zhang     jl[i] = am; il[i] = 0;
61678910aadSHong Zhang   }
61778910aadSHong Zhang 
61878910aadSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
61978910aadSHong Zhang   nlnk = am + 1;
62078910aadSHong Zhang   ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
62178910aadSHong Zhang 
62278910aadSHong Zhang   /* initial FreeSpace size is fill*(ai[am]+1) */
623a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
62478910aadSHong Zhang   current_space = free_space;
625a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr);
62678910aadSHong Zhang   current_space_lvl = free_space_lvl;
62778910aadSHong Zhang 
62878910aadSHong Zhang   for (k=0; k<am; k++){  /* for each active row k */
62978910aadSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
63078910aadSHong Zhang     nzk   = 0;
63178910aadSHong Zhang     ncols = ai[rip[k]+1] - ai[rip[k]];
63278910aadSHong Zhang     ncols_upper = 0;
63378910aadSHong Zhang     cols        = cols_lvl + am;
63478910aadSHong Zhang     for (j=0; j<ncols; j++){
63578910aadSHong Zhang       i = rip[*(aj + ai[rip[k]] + j)];
63678910aadSHong Zhang       if (i >= k){ /* only take upper triangular entry */
63778910aadSHong Zhang         cols[ncols_upper] = i;
63878910aadSHong Zhang         cols_lvl[ncols_upper] = -1;  /* initialize level for nonzero entries */
63978910aadSHong Zhang         ncols_upper++;
64078910aadSHong Zhang       }
64178910aadSHong Zhang     }
64278910aadSHong Zhang     ierr = PetscIncompleteLLAdd(ncols_upper,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
64378910aadSHong Zhang     nzk += nlnk;
64478910aadSHong Zhang 
64578910aadSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
64678910aadSHong Zhang     prow = jl[k]; /* 1st pivot row */
64778910aadSHong Zhang 
64878910aadSHong Zhang     while (prow < k){
64978910aadSHong Zhang       nextprow = jl[prow];
65078910aadSHong Zhang 
65178910aadSHong Zhang       /* merge prow into k-th row */
65278910aadSHong Zhang       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
65378910aadSHong Zhang       jmax = ui[prow+1];
65478910aadSHong Zhang       ncols = jmax-jmin;
65578910aadSHong Zhang       i     = jmin - ui[prow];
65678910aadSHong Zhang       cols = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
65778910aadSHong Zhang       for (j=0; j<ncols; j++) cols_lvl[j] = *(uj_lvl_ptr[prow] + i + j);
6585a8e39fbSHong Zhang       ierr = PetscIncompleteLLAddSorted(ncols,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
65978910aadSHong Zhang       nzk += nlnk;
66078910aadSHong Zhang 
66178910aadSHong Zhang       /* update il and jl for prow */
66278910aadSHong Zhang       if (jmin < jmax){
66378910aadSHong Zhang         il[prow] = jmin;
66478910aadSHong Zhang         j = *cols; jl[prow] = jl[j]; jl[j] = prow;
66578910aadSHong Zhang       }
66678910aadSHong Zhang       prow = nextprow;
66778910aadSHong Zhang     }
66878910aadSHong Zhang 
66978910aadSHong Zhang     /* if free space is not available, make more free space */
67078910aadSHong Zhang     if (current_space->local_remaining<nzk) {
67178910aadSHong Zhang       i = am - k + 1; /* num of unfactored rows */
67278910aadSHong Zhang       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
673a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
674a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(i,&current_space_lvl);CHKERRQ(ierr);
67578910aadSHong Zhang       reallocs++;
67678910aadSHong Zhang     }
67778910aadSHong Zhang 
67878910aadSHong Zhang     /* copy data into free_space and free_space_lvl, then initialize lnk */
67978910aadSHong Zhang     ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
68078910aadSHong Zhang 
68178910aadSHong Zhang     /* add the k-th row into il and jl */
68278910aadSHong Zhang     if (nzk-1 > 0){
68378910aadSHong Zhang       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
68478910aadSHong Zhang       jl[k] = jl[i]; jl[i] = k;
68578910aadSHong Zhang       il[k] = ui[k] + 1;
68678910aadSHong Zhang     }
68778910aadSHong Zhang     uj_ptr[k]     = current_space->array;
68878910aadSHong Zhang     uj_lvl_ptr[k] = current_space_lvl->array;
68978910aadSHong Zhang 
69078910aadSHong Zhang     current_space->array           += nzk;
69178910aadSHong Zhang     current_space->local_used      += nzk;
69278910aadSHong Zhang     current_space->local_remaining -= nzk;
69378910aadSHong Zhang 
69478910aadSHong Zhang     current_space_lvl->array           += nzk;
69578910aadSHong Zhang     current_space_lvl->local_used      += nzk;
69678910aadSHong Zhang     current_space_lvl->local_remaining -= nzk;
69778910aadSHong Zhang 
69878910aadSHong Zhang     ui[k+1] = ui[k] + nzk;
69978910aadSHong Zhang   }
70078910aadSHong Zhang 
7016cf91177SBarry Smith #if defined(PETSC_USE_INFO)
70278910aadSHong Zhang   if (ai[am] != 0) {
70378910aadSHong Zhang     PetscReal af = ((PetscReal)(2*ui[am]-am))/((PetscReal)ai[am]);
704ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
705ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
706ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
70778910aadSHong Zhang   } else {
708ae15b995SBarry Smith     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
70978910aadSHong Zhang   }
71063ba0a88SBarry Smith #endif
71178910aadSHong Zhang 
71278910aadSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
71378910aadSHong Zhang   ierr = PetscFree(jl);CHKERRQ(ierr);
71478910aadSHong Zhang   ierr = PetscFree(cols_lvl);CHKERRQ(ierr);
71578910aadSHong Zhang 
71678910aadSHong Zhang   /* destroy list of free space and other temporary array(s) */
71778910aadSHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
718a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
71978910aadSHong Zhang   ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
720a1a86e44SBarry Smith   ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
72178910aadSHong Zhang 
72278910aadSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
723f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,fact);CHKERRQ(ierr);
724f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*fact,am,am,am,am);CHKERRQ(ierr);
72578910aadSHong Zhang   B = *fact;
72678910aadSHong Zhang   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
727ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation(B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
72878910aadSHong Zhang 
72978910aadSHong Zhang   b = (Mat_SeqSBAIJ*)B->data;
73078910aadSHong Zhang   b->singlemalloc = PETSC_FALSE;
731e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
732e6b907acSBarry Smith   b->free_ij       = PETSC_TRUE;
73378910aadSHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
73478910aadSHong Zhang   b->j    = uj;
73578910aadSHong Zhang   b->i    = ui;
73678910aadSHong Zhang   b->diag = 0;
73778910aadSHong Zhang   b->ilen = 0;
73878910aadSHong Zhang   b->imax = 0;
73978910aadSHong Zhang   b->row  = perm;
74078910aadSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
74178910aadSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
74278910aadSHong Zhang   b->icol = perm;
74378910aadSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
74478910aadSHong Zhang   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
74552e6d16bSBarry Smith   ierr    = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
74678910aadSHong Zhang   b->maxnz = b->nz = ui[am];
74778910aadSHong Zhang 
74878910aadSHong Zhang   B->factor                 = FACTOR_CHOLESKY;
74978910aadSHong Zhang   B->info.factor_mallocs    = reallocs;
75078910aadSHong Zhang   B->info.fill_ratio_given  = fill;
75178910aadSHong Zhang   if (ai[am] != 0) {
75278910aadSHong Zhang     B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
75378910aadSHong Zhang   } else {
75478910aadSHong Zhang     B->info.fill_ratio_needed = 0.0;
75578910aadSHong Zhang   }
75678910aadSHong Zhang   if (perm_identity){
75778910aadSHong Zhang     B->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
75878910aadSHong Zhang     B->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
75978910aadSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
76078910aadSHong Zhang   } else {
76178910aadSHong Zhang     (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
76278910aadSHong Zhang   }
763c05c3958SHong Zhang   PetscFunctionReturn(0);
764c05c3958SHong Zhang }
765c05c3958SHong Zhang 
766c05c3958SHong Zhang #undef __FUNCT__
767c05c3958SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqBAIJ"
768c05c3958SHong Zhang PetscErrorCode MatCholeskyFactorSymbolic_SeqBAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
769c05c3958SHong Zhang {
77078910aadSHong Zhang   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data;
77178910aadSHong Zhang   Mat_SeqSBAIJ       *b;
77278910aadSHong Zhang   Mat                B;
77378910aadSHong Zhang   PetscErrorCode     ierr;
77478910aadSHong Zhang   PetscTruth         perm_identity;
77578910aadSHong Zhang   PetscReal          fill = info->fill;
7766f60792eSHong Zhang   PetscInt           *rip,i,mbs=a->mbs,bs=A->rmap.bs,*ai=a->i,*aj=a->j,reallocs=0,prow;
77778910aadSHong Zhang   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
77878910aadSHong Zhang   PetscInt           nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr;
779a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
78078910aadSHong Zhang   PetscBT            lnkbt;
78178910aadSHong Zhang 
782c05c3958SHong Zhang   PetscFunctionBegin;
7836ad2eaddSHong Zhang   if (bs > 1) { /* convert to seqsbaij */
7846ad2eaddSHong Zhang     if (!a->sbaijMat){
785ceb03754SKris Buschelman       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
7866ad2eaddSHong Zhang     }
7876ad2eaddSHong Zhang     ierr = MatCholeskyFactorSymbolic(a->sbaijMat,perm,info,fact);CHKERRQ(ierr);
7886ad2eaddSHong Zhang     B    = *fact;
7896ad2eaddSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
7906ad2eaddSHong Zhang     PetscFunctionReturn(0);
7916ad2eaddSHong Zhang   }
7926ad2eaddSHong Zhang 
79378910aadSHong Zhang   /* check whether perm is the identity mapping */
79478910aadSHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
795c84f5b01SHong Zhang   if (!perm_identity) SETERRQ(PETSC_ERR_SUP,"Matrix reordering is not supported");
79678910aadSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
79778910aadSHong Zhang 
79878910aadSHong Zhang   /* initialization */
79978910aadSHong Zhang   ierr  = PetscMalloc((mbs+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
80078910aadSHong Zhang   ui[0] = 0;
80178910aadSHong Zhang 
80278910aadSHong Zhang   /* jl: linked list for storing indices of the pivot rows
80378910aadSHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
80478910aadSHong Zhang   ierr = PetscMalloc((3*mbs+1)*sizeof(PetscInt)+mbs*sizeof(PetscInt*),&jl);CHKERRQ(ierr);
80578910aadSHong Zhang   il     = jl + mbs;
80678910aadSHong Zhang   cols   = il + mbs;
80778910aadSHong Zhang   ui_ptr = (PetscInt**)(cols + mbs);
80878910aadSHong Zhang   for (i=0; i<mbs; i++){
80978910aadSHong Zhang     jl[i] = mbs; il[i] = 0;
81078910aadSHong Zhang   }
81178910aadSHong Zhang 
81278910aadSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
81378910aadSHong Zhang   nlnk = mbs + 1;
81478910aadSHong Zhang   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
81578910aadSHong Zhang 
81678910aadSHong Zhang   /* initial FreeSpace size is fill*(ai[mbs]+1) */
817a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[mbs]+1)),&free_space);CHKERRQ(ierr);
81878910aadSHong Zhang   current_space = free_space;
81978910aadSHong Zhang 
82078910aadSHong Zhang   for (k=0; k<mbs; k++){  /* for each active row k */
82178910aadSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
82278910aadSHong Zhang     nzk   = 0;
82378910aadSHong Zhang     ncols = ai[rip[k]+1] - ai[rip[k]];
82478910aadSHong Zhang     ncols_upper = 0;
82578910aadSHong Zhang     for (j=0; j<ncols; j++){
82678910aadSHong Zhang       i = rip[*(aj + ai[rip[k]] + j)];
82778910aadSHong Zhang       if (i >= k){ /* only take upper triangular entry */
82878910aadSHong Zhang         cols[ncols_upper] = i;
82978910aadSHong Zhang         ncols_upper++;
83078910aadSHong Zhang       }
83178910aadSHong Zhang     }
83278910aadSHong Zhang     ierr = PetscLLAdd(ncols_upper,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
83378910aadSHong Zhang     nzk += nlnk;
83478910aadSHong Zhang 
83578910aadSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
83678910aadSHong Zhang     prow = jl[k]; /* 1st pivot row */
83778910aadSHong Zhang 
83878910aadSHong Zhang     while (prow < k){
83978910aadSHong Zhang       nextprow = jl[prow];
84078910aadSHong Zhang       /* merge prow into k-th row */
84178910aadSHong Zhang       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
84278910aadSHong Zhang       jmax = ui[prow+1];
84378910aadSHong Zhang       ncols = jmax-jmin;
84478910aadSHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
8455a8e39fbSHong Zhang       ierr = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
84678910aadSHong Zhang       nzk += nlnk;
84778910aadSHong Zhang 
84878910aadSHong Zhang       /* update il and jl for prow */
84978910aadSHong Zhang       if (jmin < jmax){
85078910aadSHong Zhang         il[prow] = jmin;
85178910aadSHong Zhang         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
85278910aadSHong Zhang       }
85378910aadSHong Zhang       prow = nextprow;
85478910aadSHong Zhang     }
85578910aadSHong Zhang 
85678910aadSHong Zhang     /* if free space is not available, make more free space */
85778910aadSHong Zhang     if (current_space->local_remaining<nzk) {
85878910aadSHong Zhang       i = mbs - k + 1; /* num of unfactored rows */
85978910aadSHong Zhang       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
860a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
86178910aadSHong Zhang       reallocs++;
86278910aadSHong Zhang     }
86378910aadSHong Zhang 
86478910aadSHong Zhang     /* copy data into free space, then initialize lnk */
86578910aadSHong Zhang     ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
86678910aadSHong Zhang 
86778910aadSHong Zhang     /* add the k-th row into il and jl */
86878910aadSHong Zhang     if (nzk-1 > 0){
86978910aadSHong Zhang       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
87078910aadSHong Zhang       jl[k] = jl[i]; jl[i] = k;
87178910aadSHong Zhang       il[k] = ui[k] + 1;
87278910aadSHong Zhang     }
87378910aadSHong Zhang     ui_ptr[k] = current_space->array;
87478910aadSHong Zhang     current_space->array           += nzk;
87578910aadSHong Zhang     current_space->local_used      += nzk;
87678910aadSHong Zhang     current_space->local_remaining -= nzk;
87778910aadSHong Zhang 
87878910aadSHong Zhang     ui[k+1] = ui[k] + nzk;
87978910aadSHong Zhang   }
88078910aadSHong Zhang 
8816cf91177SBarry Smith #if defined(PETSC_USE_INFO)
88278910aadSHong Zhang   if (ai[mbs] != 0) {
88378910aadSHong Zhang     PetscReal af = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
884ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
885ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
886ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
88778910aadSHong Zhang   } else {
888ae15b995SBarry Smith     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
88978910aadSHong Zhang   }
89063ba0a88SBarry Smith #endif
89178910aadSHong Zhang 
89278910aadSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
89378910aadSHong Zhang   ierr = PetscFree(jl);CHKERRQ(ierr);
89478910aadSHong Zhang 
89578910aadSHong Zhang   /* destroy list of free space and other temporary array(s) */
89678910aadSHong Zhang   ierr = PetscMalloc((ui[mbs]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
897a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
89878910aadSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
89978910aadSHong Zhang 
90078910aadSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
901f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,fact);CHKERRQ(ierr);
902f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*fact,mbs,mbs,mbs,mbs);CHKERRQ(ierr);
90378910aadSHong Zhang   B    = *fact;
90478910aadSHong Zhang   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
905ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
90678910aadSHong Zhang 
90778910aadSHong Zhang   b = (Mat_SeqSBAIJ*)B->data;
90878910aadSHong Zhang   b->singlemalloc = PETSC_FALSE;
909e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
910e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
91178910aadSHong Zhang   ierr = PetscMalloc((ui[mbs]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
91278910aadSHong Zhang   b->j    = uj;
91378910aadSHong Zhang   b->i    = ui;
91478910aadSHong Zhang   b->diag = 0;
91578910aadSHong Zhang   b->ilen = 0;
91678910aadSHong Zhang   b->imax = 0;
91778910aadSHong Zhang   b->row  = perm;
91878910aadSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
91978910aadSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
92078910aadSHong Zhang   b->icol = perm;
92178910aadSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
92278910aadSHong Zhang   ierr    = PetscMalloc((mbs+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
92352e6d16bSBarry Smith   ierr    = PetscLogObjectMemory(B,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
92478910aadSHong Zhang   b->maxnz = b->nz = ui[mbs];
92578910aadSHong Zhang 
92678910aadSHong Zhang   B->factor                 = FACTOR_CHOLESKY;
92778910aadSHong Zhang   B->info.factor_mallocs    = reallocs;
92878910aadSHong Zhang   B->info.fill_ratio_given  = fill;
92978910aadSHong Zhang   if (ai[mbs] != 0) {
93078910aadSHong Zhang     B->info.fill_ratio_needed = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
93178910aadSHong Zhang   } else {
93278910aadSHong Zhang     B->info.fill_ratio_needed = 0.0;
93378910aadSHong Zhang   }
93478910aadSHong Zhang   if (perm_identity){
9356ad2eaddSHong Zhang     B->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
9366ad2eaddSHong Zhang     B->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
9376ad2eaddSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
93878910aadSHong Zhang   } else {
9396ad2eaddSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
94078910aadSHong Zhang   }
941c05c3958SHong Zhang   PetscFunctionReturn(0);
942c05c3958SHong Zhang }
943