xref: /petsc/src/mat/impls/baij/seq/baijfact.c (revision 5c9eb25f753630cfd361293e05e7358a00a954ac)
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;
2762bba022SBarry Smith   PetscReal      shift = info->shiftinblocks;
284eeb42bcSBarry Smith 
293a40ed3dSBarry Smith   PetscFunctionBegin;
304eeb42bcSBarry Smith   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
314eeb42bcSBarry Smith   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
32b0a32e0cSBarry Smith   ierr  = PetscMalloc(4*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
334eeb42bcSBarry Smith 
344eeb42bcSBarry Smith   for (i=0; i<n; i++) {
354078e994SBarry Smith     nz    = bi[i+1] - bi[i];
364078e994SBarry Smith     ajtmp = bj + bi[i];
374eeb42bcSBarry Smith     for  (j=0; j<nz; j++) {
384eeb42bcSBarry Smith       x = rtmp+4*ajtmp[j]; x[0] = x[1] = x[2] = x[3] = 0.0;
394eeb42bcSBarry Smith     }
404eeb42bcSBarry Smith     /* load in initial (unfactored row) */
414eeb42bcSBarry Smith     idx      = r[i];
424078e994SBarry Smith     nz       = ai[idx+1] - ai[idx];
434078e994SBarry Smith     ajtmpold = aj + ai[idx];
444078e994SBarry Smith     v        = aa + 4*ai[idx];
454eeb42bcSBarry Smith     for (j=0; j<nz; j++) {
464eeb42bcSBarry Smith       x    = rtmp+4*ic[ajtmpold[j]];
474eeb42bcSBarry Smith       x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3];
484eeb42bcSBarry Smith       v    += 4;
494eeb42bcSBarry Smith     }
504eeb42bcSBarry Smith     row = *ajtmp++;
514eeb42bcSBarry Smith     while (row < i) {
524eeb42bcSBarry Smith       pc = rtmp + 4*row;
534eeb42bcSBarry Smith       p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3];
5488685aaeSLois Curfman McInnes       if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0) {
554078e994SBarry Smith         pv = ba + 4*diag_offset[row];
564078e994SBarry Smith         pj = bj + diag_offset[row] + 1;
574eeb42bcSBarry Smith         x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
584eeb42bcSBarry Smith         pc[0] = m1 = p1*x1 + p3*x2;
594eeb42bcSBarry Smith         pc[1] = m2 = p2*x1 + p4*x2;
604eeb42bcSBarry Smith         pc[2] = m3 = p1*x3 + p3*x4;
614eeb42bcSBarry Smith         pc[3] = m4 = p2*x3 + p4*x4;
624078e994SBarry Smith         nz = bi[row+1] - diag_offset[row] - 1;
634eeb42bcSBarry Smith         pv += 4;
644eeb42bcSBarry Smith         for (j=0; j<nz; j++) {
654eeb42bcSBarry Smith           x1   = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
664eeb42bcSBarry Smith           x    = rtmp + 4*pj[j];
674eeb42bcSBarry Smith           x[0] -= m1*x1 + m3*x2;
684eeb42bcSBarry Smith           x[1] -= m2*x1 + m4*x2;
694eeb42bcSBarry Smith           x[2] -= m1*x3 + m3*x4;
704eeb42bcSBarry Smith           x[3] -= m2*x3 + m4*x4;
714eeb42bcSBarry Smith           pv   += 4;
724eeb42bcSBarry Smith         }
73efee365bSSatish Balay         ierr = PetscLogFlops(16*nz+12);CHKERRQ(ierr);
744eeb42bcSBarry Smith       }
754eeb42bcSBarry Smith       row = *ajtmp++;
764eeb42bcSBarry Smith     }
774eeb42bcSBarry Smith     /* finished row so stick it into b->a */
784078e994SBarry Smith     pv = ba + 4*bi[i];
794078e994SBarry Smith     pj = bj + bi[i];
804078e994SBarry Smith     nz = bi[i+1] - bi[i];
814eeb42bcSBarry Smith     for (j=0; j<nz; j++) {
824eeb42bcSBarry Smith       x     = rtmp+4*pj[j];
834eeb42bcSBarry Smith       pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3];
844eeb42bcSBarry Smith       pv   += 4;
854eeb42bcSBarry Smith     }
864eeb42bcSBarry Smith     /* invert diagonal block */
874078e994SBarry Smith     w = ba + 4*diag_offset[i];
8862bba022SBarry Smith     ierr = Kernel_A_gets_inverse_A_2(w,shift);CHKERRQ(ierr);
894eeb42bcSBarry Smith   }
904eeb42bcSBarry Smith 
91606d414cSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
924eeb42bcSBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
934eeb42bcSBarry Smith   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
94*5c9eb25fSBarry Smith   C->factor = MAT_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;
11562bba022SBarry Smith   PetscReal      shift = info->shiftinblocks;
1164cd38560SBarry Smith 
1174cd38560SBarry Smith   PetscFunctionBegin;
118b0a32e0cSBarry Smith   ierr = PetscMalloc(4*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
1194cd38560SBarry Smith 
1204cd38560SBarry Smith   for (i=0; i<n; i++) {
1214cd38560SBarry Smith     nz    = bi[i+1] - bi[i];
1224cd38560SBarry Smith     ajtmp = bj + bi[i];
1234cd38560SBarry Smith     for  (j=0; j<nz; j++) {
1244cd38560SBarry Smith       x = rtmp+4*ajtmp[j];
1254cd38560SBarry Smith       x[0]  = x[1]  = x[2]  = x[3]  = 0.0;
1264cd38560SBarry Smith     }
1274cd38560SBarry Smith     /* load in initial (unfactored row) */
1284cd38560SBarry Smith     nz       = ai[i+1] - ai[i];
1294cd38560SBarry Smith     ajtmpold = aj + ai[i];
1304cd38560SBarry Smith     v        = aa + 4*ai[i];
1314cd38560SBarry Smith     for (j=0; j<nz; j++) {
1324cd38560SBarry Smith       x    = rtmp+4*ajtmpold[j];
1334cd38560SBarry Smith       x[0]  = v[0];  x[1]  = v[1];  x[2]  = v[2];  x[3]  = v[3];
1344cd38560SBarry Smith       v    += 4;
1354cd38560SBarry Smith     }
1364cd38560SBarry Smith     row = *ajtmp++;
1374cd38560SBarry Smith     while (row < i) {
1384cd38560SBarry Smith       pc  = rtmp + 4*row;
1394cd38560SBarry Smith       p1  = pc[0];  p2  = pc[1];  p3  = pc[2];  p4  = pc[3];
1404cd38560SBarry Smith       if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0) {
1414cd38560SBarry Smith         pv = ba + 4*diag_offset[row];
1424cd38560SBarry Smith         pj = bj + diag_offset[row] + 1;
1434cd38560SBarry Smith         x1  = pv[0];  x2  = pv[1];  x3  = pv[2];  x4  = pv[3];
1444cd38560SBarry Smith         pc[0] = m1 = p1*x1 + p3*x2;
1454cd38560SBarry Smith         pc[1] = m2 = p2*x1 + p4*x2;
1464cd38560SBarry Smith         pc[2] = m3 = p1*x3 + p3*x4;
1474cd38560SBarry Smith         pc[3] = m4 = p2*x3 + p4*x4;
1484cd38560SBarry Smith         nz = bi[row+1] - diag_offset[row] - 1;
1494cd38560SBarry Smith         pv += 4;
1504cd38560SBarry Smith         for (j=0; j<nz; j++) {
1514cd38560SBarry Smith           x1   = pv[0];  x2  = pv[1];   x3 = pv[2];  x4  = pv[3];
1524cd38560SBarry Smith           x    = rtmp + 4*pj[j];
1534cd38560SBarry Smith           x[0] -= m1*x1 + m3*x2;
1544cd38560SBarry Smith           x[1] -= m2*x1 + m4*x2;
1554cd38560SBarry Smith           x[2] -= m1*x3 + m3*x4;
1564cd38560SBarry Smith           x[3] -= m2*x3 + m4*x4;
1574cd38560SBarry Smith           pv   += 4;
1584cd38560SBarry Smith         }
159efee365bSSatish Balay         ierr = PetscLogFlops(16*nz+12);CHKERRQ(ierr);
1604cd38560SBarry Smith       }
1614cd38560SBarry Smith       row = *ajtmp++;
1624cd38560SBarry Smith     }
1634cd38560SBarry Smith     /* finished row so stick it into b->a */
1644cd38560SBarry Smith     pv = ba + 4*bi[i];
1654cd38560SBarry Smith     pj = bj + bi[i];
1664cd38560SBarry Smith     nz = bi[i+1] - bi[i];
1674cd38560SBarry Smith     for (j=0; j<nz; j++) {
1684cd38560SBarry Smith       x      = rtmp+4*pj[j];
1694cd38560SBarry Smith       pv[0]  = x[0];  pv[1]  = x[1];  pv[2]  = x[2];  pv[3]  = x[3];
1704cd38560SBarry Smith       pv   += 4;
1714cd38560SBarry Smith     }
1724cd38560SBarry Smith     /* invert diagonal block */
1734cd38560SBarry Smith     w = ba + 4*diag_offset[i];
17462bba022SBarry Smith     ierr = Kernel_A_gets_inverse_A_2(w,shift);CHKERRQ(ierr);
1754cd38560SBarry Smith   }
1764cd38560SBarry Smith 
177606d414cSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
178*5c9eb25fSBarry Smith   C->factor    = MAT_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);
248*5c9eb25fSBarry Smith   C->factor    = MAT_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 
254b24902e0SBarry Smith #undef __FUNCT__
255b24902e0SBarry Smith #define __FUNCT__ "MatGetFactor_seqbaij_petsc"
256b24902e0SBarry Smith PetscErrorCode MatGetFactor_seqbaij_petsc(Mat A,MatFactorType ftype,Mat *B)
257b24902e0SBarry Smith {
258b24902e0SBarry Smith   PetscInt           n = A->rmap.n;
259b24902e0SBarry Smith   PetscErrorCode     ierr;
260b24902e0SBarry Smith 
261b24902e0SBarry Smith   PetscFunctionBegin;
262b24902e0SBarry Smith   ierr = MatCreate(((PetscObject)A)->comm,B);CHKERRQ(ierr);
263b24902e0SBarry Smith   ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr);
264b24902e0SBarry Smith   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) {
265b24902e0SBarry Smith     ierr = MatSetType(*B,MATSEQBAIJ);CHKERRQ(ierr);
266b24902e0SBarry Smith   } else if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) {
267b24902e0SBarry Smith     ierr = MatSetType(*B,MATSEQSBAIJ);CHKERRQ(ierr);
268*5c9eb25fSBarry Smith     ierr = MatSeqSBAIJSetPreallocation(*B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
269*5c9eb25fSBarry Smith     (*B)->ops->iccfactorsymbolic      = MatICCFactorSymbolic_SeqBAIJ;
270*5c9eb25fSBarry Smith     (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqBAIJ;
271b24902e0SBarry Smith   } else SETERRQ(PETSC_ERR_SUP,"Factor type not supported");
272b24902e0SBarry Smith   PetscFunctionReturn(0);
273b24902e0SBarry Smith }
274b24902e0SBarry Smith 
275273d9f13SBarry Smith 
2765d517e7eSBarry Smith /* ----------------------------------------------------------- */
2774a2ae208SSatish Balay #undef __FUNCT__
2784a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqBAIJ"
279dfbe8321SBarry Smith PetscErrorCode MatLUFactor_SeqBAIJ(Mat A,IS row,IS col,MatFactorInfo *info)
2805d517e7eSBarry Smith {
281dfbe8321SBarry Smith   PetscErrorCode ierr;
2825d517e7eSBarry Smith   Mat            C;
2835d517e7eSBarry Smith 
2843a40ed3dSBarry Smith   PetscFunctionBegin;
285b24902e0SBarry Smith   ierr = MatGetFactor(A,"petsc",MAT_FACTOR_LU,&C);CHKERRQ(ierr);
286b9b97703SBarry Smith   ierr = MatLUFactorSymbolic(A,row,col,info,&C);CHKERRQ(ierr);
287af281ebdSHong Zhang   ierr = MatLUFactorNumeric(A,info,&C);CHKERRQ(ierr);
288273d9f13SBarry Smith   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
28952e6d16bSBarry Smith   ierr = PetscLogObjectParent(A,((Mat_SeqBAIJ*)(A->data))->icol);CHKERRQ(ierr);
2903a40ed3dSBarry Smith   PetscFunctionReturn(0);
2915d517e7eSBarry Smith }
2924cd38560SBarry Smith 
293c05c3958SHong Zhang #include "src/mat/impls/sbaij/seq/sbaij.h"
294c05c3958SHong Zhang #undef __FUNCT__
295c05c3958SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N"
296af281ebdSHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N(Mat A,MatFactorInfo *info,Mat *B)
297c05c3958SHong Zhang {
298f3086b4bSHong Zhang   PetscErrorCode ierr;
29978910aadSHong Zhang   Mat            C = *B;
30078910aadSHong Zhang   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data;
30178910aadSHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
30278910aadSHong Zhang   IS             ip=b->row;
303899cda47SBarry Smith   PetscInt       *rip,i,j,mbs=a->mbs,bs=A->rmap.bs,*bi=b->i,*bj=b->j,*bcol;
30478910aadSHong Zhang   PetscInt       *ai=a->i,*aj=a->j;
30578910aadSHong Zhang   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
30678910aadSHong Zhang   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
3073cea4cbeSHong Zhang   PetscReal      zeropivot,rs,shiftnz;
308fbf22428SSatish Balay   PetscReal      shiftpd;
3093cea4cbeSHong Zhang   ChShift_Ctx    sctx;
3103cea4cbeSHong Zhang   PetscInt       newshift;
31178910aadSHong Zhang 
312c05c3958SHong Zhang   PetscFunctionBegin;
3136ad2eaddSHong Zhang   if (bs > 1) {
3146ad2eaddSHong Zhang     if (!a->sbaijMat){
315ceb03754SKris Buschelman       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
3166ad2eaddSHong Zhang     }
317f3086b4bSHong Zhang     ierr = (a->sbaijMat)->ops->choleskyfactornumeric(a->sbaijMat,info,B);CHKERRQ(ierr);
3186ad2eaddSHong Zhang     ierr = MatDestroy(a->sbaijMat);CHKERRQ(ierr);
3196ad2eaddSHong Zhang     a->sbaijMat = PETSC_NULL;
3206ad2eaddSHong Zhang     PetscFunctionReturn(0);
3216ad2eaddSHong Zhang   }
32278910aadSHong Zhang 
32378910aadSHong Zhang   /* initialization */
3243cea4cbeSHong Zhang   shiftnz   = info->shiftnz;
3253cea4cbeSHong Zhang   shiftpd   = info->shiftpd;
3263cea4cbeSHong Zhang   zeropivot = info->zeropivot;
3273cea4cbeSHong Zhang 
3286ad2eaddSHong Zhang   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
32978910aadSHong Zhang   nz   = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar);
33078910aadSHong Zhang   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
33178910aadSHong Zhang   jl   = il + mbs;
33278910aadSHong Zhang   rtmp = (MatScalar*)(jl + mbs);
33378910aadSHong Zhang 
3343cea4cbeSHong Zhang   sctx.shift_amount = 0;
3353cea4cbeSHong Zhang   sctx.nshift       = 0;
33678910aadSHong Zhang   do {
3373cea4cbeSHong Zhang     sctx.chshift = PETSC_FALSE;
33878910aadSHong Zhang     for (i=0; i<mbs; i++) {
33978910aadSHong Zhang       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
34078910aadSHong Zhang     }
34178910aadSHong Zhang 
34278910aadSHong Zhang     for (k = 0; k<mbs; k++){
34378910aadSHong Zhang       bval = ba + bi[k];
34478910aadSHong Zhang       /* initialize k-th row by the perm[k]-th row of A */
34578910aadSHong Zhang       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
34678910aadSHong Zhang       for (j = jmin; j < jmax; j++){
34778910aadSHong Zhang         col = rip[aj[j]];
34878910aadSHong Zhang         if (col >= k){ /* only take upper triangular entry */
34978910aadSHong Zhang           rtmp[col] = aa[j];
35078910aadSHong Zhang           *bval++  = 0.0; /* for in-place factorization */
35178910aadSHong Zhang         }
35278910aadSHong Zhang       }
3533cea4cbeSHong Zhang 
3543cea4cbeSHong Zhang       /* shift the diagonal of the matrix */
3553cea4cbeSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
35678910aadSHong Zhang 
35778910aadSHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
35878910aadSHong Zhang       dk = rtmp[k];
35978910aadSHong Zhang       i = jl[k]; /* first row to be added to k_th row  */
36078910aadSHong Zhang 
36178910aadSHong Zhang       while (i < k){
36278910aadSHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
36378910aadSHong Zhang 
36478910aadSHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
36578910aadSHong Zhang         ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
36678910aadSHong Zhang         uikdi = - ba[ili]*ba[bi[i]];  /* diagonal(k) */
36778910aadSHong Zhang         dk += uikdi*ba[ili];
36878910aadSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
36978910aadSHong Zhang 
37078910aadSHong Zhang         /* add multiple of row i to k-th row */
37178910aadSHong Zhang         jmin = ili + 1; jmax = bi[i+1];
37278910aadSHong Zhang         if (jmin < jmax){
37378910aadSHong Zhang           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
37478910aadSHong Zhang           /* update il and jl for row i */
37578910aadSHong Zhang           il[i] = jmin;
37678910aadSHong Zhang           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
37778910aadSHong Zhang         }
37878910aadSHong Zhang         i = nexti;
37978910aadSHong Zhang       }
38078910aadSHong Zhang 
3813cea4cbeSHong Zhang       /* shift the diagonals when zero pivot is detected */
3823cea4cbeSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
3833cea4cbeSHong Zhang       rs   = 0.0;
38478910aadSHong Zhang       jmin = bi[k]+1;
38578910aadSHong Zhang       nz   = bi[k+1] - jmin;
38678910aadSHong Zhang       if (nz){
38778910aadSHong Zhang         bcol = bj + jmin;
38878910aadSHong Zhang         while (nz--){
38989f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
39089f845c8SHong Zhang           bcol++;
39178910aadSHong Zhang         }
39278910aadSHong Zhang       }
3933cea4cbeSHong Zhang 
3943cea4cbeSHong Zhang       sctx.rs = rs;
3953cea4cbeSHong Zhang       sctx.pv = dk;
39645938b79SHong Zhang       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
39745938b79SHong Zhang       if (newshift == 1) break;
39878910aadSHong Zhang 
39978910aadSHong Zhang       /* copy data into U(k,:) */
40078910aadSHong Zhang       ba[bi[k]] = 1.0/dk; /* U(k,k) */
40178910aadSHong Zhang       jmin = bi[k]+1; jmax = bi[k+1];
40278910aadSHong Zhang       if (jmin < jmax) {
40378910aadSHong Zhang         for (j=jmin; j<jmax; j++){
40478910aadSHong Zhang           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
40578910aadSHong Zhang         }
40678910aadSHong Zhang         /* add the k-th row into il and jl */
40778910aadSHong Zhang         il[k] = jmin;
40878910aadSHong Zhang         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
40978910aadSHong Zhang       }
41078910aadSHong Zhang     }
4113cea4cbeSHong Zhang   } while (sctx.chshift);
41278910aadSHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
41378910aadSHong Zhang 
41478910aadSHong Zhang   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
415*5c9eb25fSBarry Smith   C->factor       = MAT_FACTOR_CHOLESKY;
41678910aadSHong Zhang   C->assembled    = PETSC_TRUE;
41778910aadSHong Zhang   C->preallocated = PETSC_TRUE;
418899cda47SBarry Smith   ierr = PetscLogFlops(C->rmap.N);CHKERRQ(ierr);
4193cea4cbeSHong Zhang   if (sctx.nshift){
4203cea4cbeSHong Zhang     if (shiftnz) {
4211e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
4223cea4cbeSHong Zhang     } else if (shiftpd) {
4231e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
42478910aadSHong Zhang     }
42578910aadSHong Zhang   }
426c05c3958SHong Zhang   PetscFunctionReturn(0);
427c05c3958SHong Zhang }
4284cd38560SBarry Smith 
429c05c3958SHong Zhang #undef __FUNCT__
430c05c3958SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering"
431af281ebdSHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering(Mat A,MatFactorInfo *info,Mat *fact)
432c05c3958SHong Zhang {
43378910aadSHong Zhang   Mat            C = *fact;
43478910aadSHong Zhang   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data;
43578910aadSHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
43678910aadSHong Zhang   PetscErrorCode ierr;
4373cea4cbeSHong Zhang   PetscInt       i,j,am=a->mbs;
43878910aadSHong Zhang   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
4393cea4cbeSHong Zhang   PetscInt       k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz;
44078910aadSHong Zhang   MatScalar      *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval;
4413cea4cbeSHong Zhang   PetscReal      zeropivot,rs,shiftnz;
442fbf22428SSatish Balay   PetscReal      shiftpd;
4433cea4cbeSHong Zhang   ChShift_Ctx    sctx;
4443cea4cbeSHong Zhang   PetscInt       newshift;
44578910aadSHong Zhang 
446c05c3958SHong Zhang   PetscFunctionBegin;
44778910aadSHong Zhang   /* initialization */
4483cea4cbeSHong Zhang   shiftnz   = info->shiftnz;
4493cea4cbeSHong Zhang   shiftpd   = info->shiftpd;
4503cea4cbeSHong Zhang   zeropivot = info->zeropivot;
4513cea4cbeSHong Zhang 
45278910aadSHong Zhang   nz   = (2*am+1)*sizeof(PetscInt)+am*sizeof(MatScalar);
45378910aadSHong Zhang   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
45478910aadSHong Zhang   jl   = il + am;
45578910aadSHong Zhang   rtmp = (MatScalar*)(jl + am);
45678910aadSHong Zhang 
4573cea4cbeSHong Zhang   sctx.shift_amount = 0;
4583cea4cbeSHong Zhang   sctx.nshift       = 0;
45978910aadSHong Zhang   do {
4603cea4cbeSHong Zhang     sctx.chshift = PETSC_FALSE;
46178910aadSHong Zhang     for (i=0; i<am; i++) {
46278910aadSHong Zhang       rtmp[i] = 0.0; jl[i] = am; il[0] = 0;
46378910aadSHong Zhang     }
46478910aadSHong Zhang 
46578910aadSHong Zhang     for (k = 0; k<am; k++){
46678910aadSHong Zhang     /* initialize k-th row with elements nonzero in row perm(k) of A */
46778910aadSHong Zhang       nz   = ai[k+1] - ai[k];
46878910aadSHong Zhang       acol = aj + ai[k];
46978910aadSHong Zhang       aval = aa + ai[k];
47078910aadSHong Zhang       bval = ba + bi[k];
47178910aadSHong Zhang       while (nz -- ){
47278910aadSHong Zhang         if (*acol < k) { /* skip lower triangular entries */
47378910aadSHong Zhang           acol++; aval++;
47478910aadSHong Zhang         } else {
47578910aadSHong Zhang           rtmp[*acol++] = *aval++;
47678910aadSHong Zhang           *bval++       = 0.0; /* for in-place factorization */
47778910aadSHong Zhang         }
47878910aadSHong Zhang       }
4793cea4cbeSHong Zhang 
4803cea4cbeSHong Zhang       /* shift the diagonal of the matrix */
4813cea4cbeSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
48278910aadSHong Zhang 
48378910aadSHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
48478910aadSHong Zhang       dk = rtmp[k];
48578910aadSHong Zhang       i  = jl[k]; /* first row to be added to k_th row  */
48678910aadSHong Zhang 
48778910aadSHong Zhang       while (i < k){
48878910aadSHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
48978910aadSHong Zhang         /* compute multiplier, update D(k) and U(i,k) */
49078910aadSHong Zhang         ili   = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
49178910aadSHong Zhang         uikdi = - ba[ili]*ba[bi[i]];
49278910aadSHong Zhang         dk   += uikdi*ba[ili];
49378910aadSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
49478910aadSHong Zhang 
49578910aadSHong Zhang         /* add multiple of row i to k-th row ... */
49678910aadSHong Zhang         jmin = ili + 1;
49778910aadSHong Zhang         nz   = bi[i+1] - jmin;
49878910aadSHong Zhang         if (nz > 0){
49978910aadSHong Zhang           bcol = bj + jmin;
50078910aadSHong Zhang           bval = ba + jmin;
50178910aadSHong Zhang           while (nz --) rtmp[*bcol++] += uikdi*(*bval++);
50278910aadSHong Zhang           /* update il and jl for i-th row */
50378910aadSHong Zhang           il[i] = jmin;
50478910aadSHong Zhang           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
50578910aadSHong Zhang         }
50678910aadSHong Zhang         i = nexti;
50778910aadSHong Zhang       }
50878910aadSHong Zhang 
5093cea4cbeSHong Zhang       /* shift the diagonals when zero pivot is detected */
5103cea4cbeSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
5113cea4cbeSHong Zhang       rs   = 0.0;
51278910aadSHong Zhang       jmin = bi[k]+1;
51378910aadSHong Zhang       nz   = bi[k+1] - jmin;
51478910aadSHong Zhang       if (nz){
51578910aadSHong Zhang         bcol = bj + jmin;
51678910aadSHong Zhang         while (nz--){
51789f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
51889f845c8SHong Zhang           bcol++;
51978910aadSHong Zhang         }
52078910aadSHong Zhang       }
5213cea4cbeSHong Zhang 
5223cea4cbeSHong Zhang       sctx.rs = rs;
5233cea4cbeSHong Zhang       sctx.pv = dk;
52445938b79SHong Zhang       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
52545938b79SHong Zhang       if (newshift == 1) break;    /* sctx.shift_amount is updated */
52678910aadSHong Zhang 
52778910aadSHong Zhang       /* copy data into U(k,:) */
52878910aadSHong Zhang       ba[bi[k]] = 1.0/dk;
52978910aadSHong Zhang       jmin      = bi[k]+1;
53078910aadSHong Zhang       nz        = bi[k+1] - jmin;
53178910aadSHong Zhang       if (nz){
53278910aadSHong Zhang         bcol = bj + jmin;
53378910aadSHong Zhang         bval = ba + jmin;
53478910aadSHong Zhang         while (nz--){
53578910aadSHong Zhang           *bval++       = rtmp[*bcol];
53678910aadSHong Zhang           rtmp[*bcol++] = 0.0;
53778910aadSHong Zhang         }
53878910aadSHong Zhang         /* add k-th row into il and jl */
53978910aadSHong Zhang         il[k] = jmin;
54078910aadSHong Zhang         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
54178910aadSHong Zhang       }
54278910aadSHong Zhang     }
5433cea4cbeSHong Zhang   } while (sctx.chshift);
54478910aadSHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
54578910aadSHong Zhang 
546*5c9eb25fSBarry Smith   C->factor       = MAT_FACTOR_CHOLESKY;
54778910aadSHong Zhang   C->assembled    = PETSC_TRUE;
54878910aadSHong Zhang   C->preallocated = PETSC_TRUE;
549899cda47SBarry Smith   ierr = PetscLogFlops(C->rmap.N);CHKERRQ(ierr);
5503cea4cbeSHong Zhang     if (sctx.nshift){
5513cea4cbeSHong Zhang     if (shiftnz) {
5521e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
5533cea4cbeSHong Zhang     } else if (shiftpd) {
5541e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
55578910aadSHong Zhang     }
55678910aadSHong Zhang   }
557c05c3958SHong Zhang   PetscFunctionReturn(0);
558c05c3958SHong Zhang }
559c05c3958SHong Zhang 
560c05c3958SHong Zhang #include "petscbt.h"
561c05c3958SHong Zhang #include "src/mat/utils/freespace.h"
562c05c3958SHong Zhang #undef __FUNCT__
563c05c3958SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqBAIJ"
564c05c3958SHong Zhang PetscErrorCode MatICCFactorSymbolic_SeqBAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
565c05c3958SHong Zhang {
56678910aadSHong Zhang   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data;
56778910aadSHong Zhang   Mat_SeqSBAIJ       *b;
56878910aadSHong Zhang   Mat                B;
56978910aadSHong Zhang   PetscErrorCode     ierr;
57078910aadSHong Zhang   PetscTruth         perm_identity;
571899cda47SBarry Smith   PetscInt           reallocs=0,*rip,i,*ai=a->i,*aj=a->j,am=a->mbs,bs=A->rmap.bs,*ui;
57278910aadSHong Zhang   PetscInt           jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow;
573cfdb8b8aSHong Zhang   PetscInt           nlnk,*lnk,*lnk_lvl=PETSC_NULL,ncols,ncols_upper,*cols,*cols_lvl,*uj,**uj_ptr,**uj_lvl_ptr;
57478910aadSHong Zhang   PetscReal          fill=info->fill,levels=info->levels;
575a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
576a1a86e44SBarry Smith   PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
57778910aadSHong Zhang   PetscBT            lnkbt;
57878910aadSHong Zhang 
579c05c3958SHong Zhang   PetscFunctionBegin;
5806ad2eaddSHong Zhang   if (bs > 1){
5816ad2eaddSHong Zhang     if (!a->sbaijMat){
582ceb03754SKris Buschelman       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
5836ad2eaddSHong Zhang     }
584*5c9eb25fSBarry Smith     (*fact)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqSBAIJ;  /* undue the change made in MatGetFactor_seqbaij_petsc */
5856ad2eaddSHong Zhang     ierr = MatICCFactorSymbolic(a->sbaijMat,perm,info,fact);CHKERRQ(ierr);
5866ad2eaddSHong Zhang     PetscFunctionReturn(0);
5876ad2eaddSHong Zhang   }
5886ad2eaddSHong Zhang 
58978910aadSHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
59078910aadSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
59178910aadSHong Zhang 
59278910aadSHong Zhang   /* special case that simply copies fill pattern */
59378910aadSHong Zhang   if (!levels && perm_identity) {
59478910aadSHong Zhang     ierr = MatMarkDiagonal_SeqBAIJ(A);CHKERRQ(ierr);
59578910aadSHong Zhang     ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
59678910aadSHong Zhang     for (i=0; i<am; i++) {
59778910aadSHong Zhang       ui[i] = ai[i+1] - a->diag[i]; /* ui: rowlengths - changes when !perm_identity */
59878910aadSHong Zhang     }
59978910aadSHong Zhang     B = *fact;
60078910aadSHong Zhang     ierr = MatSeqSBAIJSetPreallocation(B,1,0,ui);CHKERRQ(ierr);
60178910aadSHong Zhang 
602b24902e0SBarry Smith 
60378910aadSHong Zhang     b  = (Mat_SeqSBAIJ*)B->data;
60478910aadSHong Zhang     uj = b->j;
60578910aadSHong Zhang     for (i=0; i<am; i++) {
60678910aadSHong Zhang       aj = a->j + a->diag[i];
60778910aadSHong Zhang       for (j=0; j<ui[i]; j++){
60878910aadSHong Zhang         *uj++ = *aj++;
60978910aadSHong Zhang       }
61078910aadSHong Zhang       b->ilen[i] = ui[i];
61178910aadSHong Zhang     }
61278910aadSHong Zhang     ierr = PetscFree(ui);CHKERRQ(ierr);
61378910aadSHong Zhang     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
61478910aadSHong Zhang     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
61578910aadSHong Zhang 
61678910aadSHong Zhang     B->ops->solve                 = MatSolve_SeqSBAIJ_1_NaturalOrdering;
61778910aadSHong Zhang     B->ops->solvetranspose        = MatSolve_SeqSBAIJ_1_NaturalOrdering;
61878910aadSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
61978910aadSHong Zhang     PetscFunctionReturn(0);
62078910aadSHong Zhang   }
62178910aadSHong Zhang 
62278910aadSHong Zhang   /* initialization */
62378910aadSHong Zhang   ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
62478910aadSHong Zhang   ui[0] = 0;
62578910aadSHong Zhang   ierr  = PetscMalloc((2*am+1)*sizeof(PetscInt),&cols_lvl);CHKERRQ(ierr);
62678910aadSHong Zhang 
62778910aadSHong Zhang   /* jl: linked list for storing indices of the pivot rows
62878910aadSHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
62978910aadSHong Zhang   ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt*),&jl);CHKERRQ(ierr);
63078910aadSHong Zhang   il         = jl + am;
63178910aadSHong Zhang   uj_ptr     = (PetscInt**)(il + am);
63278910aadSHong Zhang   uj_lvl_ptr = (PetscInt**)(uj_ptr + am);
63378910aadSHong Zhang   for (i=0; i<am; i++){
63478910aadSHong Zhang     jl[i] = am; il[i] = 0;
63578910aadSHong Zhang   }
63678910aadSHong Zhang 
63778910aadSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
63878910aadSHong Zhang   nlnk = am + 1;
63978910aadSHong Zhang   ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
64078910aadSHong Zhang 
64178910aadSHong Zhang   /* initial FreeSpace size is fill*(ai[am]+1) */
642a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
64378910aadSHong Zhang   current_space = free_space;
644a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr);
64578910aadSHong Zhang   current_space_lvl = free_space_lvl;
64678910aadSHong Zhang 
64778910aadSHong Zhang   for (k=0; k<am; k++){  /* for each active row k */
64878910aadSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
64978910aadSHong Zhang     nzk   = 0;
65078910aadSHong Zhang     ncols = ai[rip[k]+1] - ai[rip[k]];
65178910aadSHong Zhang     ncols_upper = 0;
65278910aadSHong Zhang     cols        = cols_lvl + am;
65378910aadSHong Zhang     for (j=0; j<ncols; j++){
65478910aadSHong Zhang       i = rip[*(aj + ai[rip[k]] + j)];
65578910aadSHong Zhang       if (i >= k){ /* only take upper triangular entry */
65678910aadSHong Zhang         cols[ncols_upper] = i;
65778910aadSHong Zhang         cols_lvl[ncols_upper] = -1;  /* initialize level for nonzero entries */
65878910aadSHong Zhang         ncols_upper++;
65978910aadSHong Zhang       }
66078910aadSHong Zhang     }
66178910aadSHong Zhang     ierr = PetscIncompleteLLAdd(ncols_upper,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
66278910aadSHong Zhang     nzk += nlnk;
66378910aadSHong Zhang 
66478910aadSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
66578910aadSHong Zhang     prow = jl[k]; /* 1st pivot row */
66678910aadSHong Zhang 
66778910aadSHong Zhang     while (prow < k){
66878910aadSHong Zhang       nextprow = jl[prow];
66978910aadSHong Zhang 
67078910aadSHong Zhang       /* merge prow into k-th row */
67178910aadSHong Zhang       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
67278910aadSHong Zhang       jmax = ui[prow+1];
67378910aadSHong Zhang       ncols = jmax-jmin;
67478910aadSHong Zhang       i     = jmin - ui[prow];
67578910aadSHong Zhang       cols = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
67678910aadSHong Zhang       for (j=0; j<ncols; j++) cols_lvl[j] = *(uj_lvl_ptr[prow] + i + j);
6775a8e39fbSHong Zhang       ierr = PetscIncompleteLLAddSorted(ncols,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
67878910aadSHong Zhang       nzk += nlnk;
67978910aadSHong Zhang 
68078910aadSHong Zhang       /* update il and jl for prow */
68178910aadSHong Zhang       if (jmin < jmax){
68278910aadSHong Zhang         il[prow] = jmin;
68378910aadSHong Zhang         j = *cols; jl[prow] = jl[j]; jl[j] = prow;
68478910aadSHong Zhang       }
68578910aadSHong Zhang       prow = nextprow;
68678910aadSHong Zhang     }
68778910aadSHong Zhang 
68878910aadSHong Zhang     /* if free space is not available, make more free space */
68978910aadSHong Zhang     if (current_space->local_remaining<nzk) {
69078910aadSHong Zhang       i = am - k + 1; /* num of unfactored rows */
69178910aadSHong Zhang       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
692a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
693a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(i,&current_space_lvl);CHKERRQ(ierr);
69478910aadSHong Zhang       reallocs++;
69578910aadSHong Zhang     }
69678910aadSHong Zhang 
69778910aadSHong Zhang     /* copy data into free_space and free_space_lvl, then initialize lnk */
69878910aadSHong Zhang     ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
69978910aadSHong Zhang 
70078910aadSHong Zhang     /* add the k-th row into il and jl */
70178910aadSHong Zhang     if (nzk-1 > 0){
70278910aadSHong Zhang       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
70378910aadSHong Zhang       jl[k] = jl[i]; jl[i] = k;
70478910aadSHong Zhang       il[k] = ui[k] + 1;
70578910aadSHong Zhang     }
70678910aadSHong Zhang     uj_ptr[k]     = current_space->array;
70778910aadSHong Zhang     uj_lvl_ptr[k] = current_space_lvl->array;
70878910aadSHong Zhang 
70978910aadSHong Zhang     current_space->array           += nzk;
71078910aadSHong Zhang     current_space->local_used      += nzk;
71178910aadSHong Zhang     current_space->local_remaining -= nzk;
71278910aadSHong Zhang 
71378910aadSHong Zhang     current_space_lvl->array           += nzk;
71478910aadSHong Zhang     current_space_lvl->local_used      += nzk;
71578910aadSHong Zhang     current_space_lvl->local_remaining -= nzk;
71678910aadSHong Zhang 
71778910aadSHong Zhang     ui[k+1] = ui[k] + nzk;
71878910aadSHong Zhang   }
71978910aadSHong Zhang 
7206cf91177SBarry Smith #if defined(PETSC_USE_INFO)
72178910aadSHong Zhang   if (ai[am] != 0) {
72278910aadSHong Zhang     PetscReal af = ((PetscReal)(2*ui[am]-am))/((PetscReal)ai[am]);
723ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
724ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
725ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
72678910aadSHong Zhang   } else {
727ae15b995SBarry Smith     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
72878910aadSHong Zhang   }
72963ba0a88SBarry Smith #endif
73078910aadSHong Zhang 
73178910aadSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
73278910aadSHong Zhang   ierr = PetscFree(jl);CHKERRQ(ierr);
73378910aadSHong Zhang   ierr = PetscFree(cols_lvl);CHKERRQ(ierr);
73478910aadSHong Zhang 
73578910aadSHong Zhang   /* destroy list of free space and other temporary array(s) */
73678910aadSHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
737a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
73878910aadSHong Zhang   ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
739a1a86e44SBarry Smith   ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
74078910aadSHong Zhang 
74178910aadSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
74278910aadSHong Zhang   B = *fact;
743ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation(B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
74478910aadSHong Zhang 
74578910aadSHong Zhang   b = (Mat_SeqSBAIJ*)B->data;
74678910aadSHong Zhang   b->singlemalloc = PETSC_FALSE;
747e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
748e6b907acSBarry Smith   b->free_ij       = PETSC_TRUE;
74978910aadSHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
75078910aadSHong Zhang   b->j    = uj;
75178910aadSHong Zhang   b->i    = ui;
75278910aadSHong Zhang   b->diag = 0;
75378910aadSHong Zhang   b->ilen = 0;
75478910aadSHong Zhang   b->imax = 0;
75578910aadSHong Zhang   b->row  = perm;
75678910aadSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
75778910aadSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
75878910aadSHong Zhang   b->icol = perm;
75978910aadSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
76078910aadSHong Zhang   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
76152e6d16bSBarry Smith   ierr    = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
76278910aadSHong Zhang   b->maxnz = b->nz = ui[am];
76378910aadSHong Zhang 
764*5c9eb25fSBarry Smith   B->factor                 = MAT_FACTOR_CHOLESKY;
76578910aadSHong Zhang   B->info.factor_mallocs    = reallocs;
76678910aadSHong Zhang   B->info.fill_ratio_given  = fill;
76778910aadSHong Zhang   if (ai[am] != 0) {
76878910aadSHong Zhang     B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
76978910aadSHong Zhang   } else {
77078910aadSHong Zhang     B->info.fill_ratio_needed = 0.0;
77178910aadSHong Zhang   }
77278910aadSHong Zhang   if (perm_identity){
77378910aadSHong Zhang     B->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
77478910aadSHong Zhang     B->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
77578910aadSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
77678910aadSHong Zhang   } else {
77778910aadSHong Zhang     (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
77878910aadSHong Zhang   }
779c05c3958SHong Zhang   PetscFunctionReturn(0);
780c05c3958SHong Zhang }
781c05c3958SHong Zhang 
782c05c3958SHong Zhang #undef __FUNCT__
783c05c3958SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqBAIJ"
784c05c3958SHong Zhang PetscErrorCode MatCholeskyFactorSymbolic_SeqBAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
785c05c3958SHong Zhang {
78678910aadSHong Zhang   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data;
78778910aadSHong Zhang   Mat_SeqSBAIJ       *b;
78878910aadSHong Zhang   Mat                B;
78978910aadSHong Zhang   PetscErrorCode     ierr;
79078910aadSHong Zhang   PetscTruth         perm_identity;
79178910aadSHong Zhang   PetscReal          fill = info->fill;
7926f60792eSHong Zhang   PetscInt           *rip,i,mbs=a->mbs,bs=A->rmap.bs,*ai=a->i,*aj=a->j,reallocs=0,prow;
79378910aadSHong Zhang   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
79478910aadSHong Zhang   PetscInt           nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr;
795a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
79678910aadSHong Zhang   PetscBT            lnkbt;
79778910aadSHong Zhang 
798c05c3958SHong Zhang   PetscFunctionBegin;
7996ad2eaddSHong Zhang   if (bs > 1) { /* convert to seqsbaij */
8006ad2eaddSHong Zhang     if (!a->sbaijMat){
801ceb03754SKris Buschelman       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
8026ad2eaddSHong Zhang     }
803*5c9eb25fSBarry Smith     (*fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqSBAIJ; /* undue the change made in MatGetFactor_seqbaij_petsc */
8046ad2eaddSHong Zhang     ierr = MatCholeskyFactorSymbolic(a->sbaijMat,perm,info,fact);CHKERRQ(ierr);
8056ad2eaddSHong Zhang     PetscFunctionReturn(0);
8066ad2eaddSHong Zhang   }
8076ad2eaddSHong Zhang 
80878910aadSHong Zhang   /* check whether perm is the identity mapping */
80978910aadSHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
810c84f5b01SHong Zhang   if (!perm_identity) SETERRQ(PETSC_ERR_SUP,"Matrix reordering is not supported");
81178910aadSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
81278910aadSHong Zhang 
81378910aadSHong Zhang   /* initialization */
81478910aadSHong Zhang   ierr  = PetscMalloc((mbs+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
81578910aadSHong Zhang   ui[0] = 0;
81678910aadSHong Zhang 
81778910aadSHong Zhang   /* jl: linked list for storing indices of the pivot rows
81878910aadSHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
81978910aadSHong Zhang   ierr = PetscMalloc((3*mbs+1)*sizeof(PetscInt)+mbs*sizeof(PetscInt*),&jl);CHKERRQ(ierr);
82078910aadSHong Zhang   il     = jl + mbs;
82178910aadSHong Zhang   cols   = il + mbs;
82278910aadSHong Zhang   ui_ptr = (PetscInt**)(cols + mbs);
82378910aadSHong Zhang   for (i=0; i<mbs; i++){
82478910aadSHong Zhang     jl[i] = mbs; il[i] = 0;
82578910aadSHong Zhang   }
82678910aadSHong Zhang 
82778910aadSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
82878910aadSHong Zhang   nlnk = mbs + 1;
82978910aadSHong Zhang   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
83078910aadSHong Zhang 
83178910aadSHong Zhang   /* initial FreeSpace size is fill*(ai[mbs]+1) */
832a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[mbs]+1)),&free_space);CHKERRQ(ierr);
83378910aadSHong Zhang   current_space = free_space;
83478910aadSHong Zhang 
83578910aadSHong Zhang   for (k=0; k<mbs; k++){  /* for each active row k */
83678910aadSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
83778910aadSHong Zhang     nzk   = 0;
83878910aadSHong Zhang     ncols = ai[rip[k]+1] - ai[rip[k]];
83978910aadSHong Zhang     ncols_upper = 0;
84078910aadSHong Zhang     for (j=0; j<ncols; j++){
84178910aadSHong Zhang       i = rip[*(aj + ai[rip[k]] + j)];
84278910aadSHong Zhang       if (i >= k){ /* only take upper triangular entry */
84378910aadSHong Zhang         cols[ncols_upper] = i;
84478910aadSHong Zhang         ncols_upper++;
84578910aadSHong Zhang       }
84678910aadSHong Zhang     }
84778910aadSHong Zhang     ierr = PetscLLAdd(ncols_upper,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
84878910aadSHong Zhang     nzk += nlnk;
84978910aadSHong Zhang 
85078910aadSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
85178910aadSHong Zhang     prow = jl[k]; /* 1st pivot row */
85278910aadSHong Zhang 
85378910aadSHong Zhang     while (prow < k){
85478910aadSHong Zhang       nextprow = jl[prow];
85578910aadSHong Zhang       /* merge prow into k-th row */
85678910aadSHong Zhang       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
85778910aadSHong Zhang       jmax = ui[prow+1];
85878910aadSHong Zhang       ncols = jmax-jmin;
85978910aadSHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
8605a8e39fbSHong Zhang       ierr = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
86178910aadSHong Zhang       nzk += nlnk;
86278910aadSHong Zhang 
86378910aadSHong Zhang       /* update il and jl for prow */
86478910aadSHong Zhang       if (jmin < jmax){
86578910aadSHong Zhang         il[prow] = jmin;
86678910aadSHong Zhang         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
86778910aadSHong Zhang       }
86878910aadSHong Zhang       prow = nextprow;
86978910aadSHong Zhang     }
87078910aadSHong Zhang 
87178910aadSHong Zhang     /* if free space is not available, make more free space */
87278910aadSHong Zhang     if (current_space->local_remaining<nzk) {
87378910aadSHong Zhang       i = mbs - k + 1; /* num of unfactored rows */
87478910aadSHong Zhang       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
875a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
87678910aadSHong Zhang       reallocs++;
87778910aadSHong Zhang     }
87878910aadSHong Zhang 
87978910aadSHong Zhang     /* copy data into free space, then initialize lnk */
88078910aadSHong Zhang     ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
88178910aadSHong Zhang 
88278910aadSHong Zhang     /* add the k-th row into il and jl */
88378910aadSHong Zhang     if (nzk-1 > 0){
88478910aadSHong Zhang       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
88578910aadSHong Zhang       jl[k] = jl[i]; jl[i] = k;
88678910aadSHong Zhang       il[k] = ui[k] + 1;
88778910aadSHong Zhang     }
88878910aadSHong Zhang     ui_ptr[k] = current_space->array;
88978910aadSHong Zhang     current_space->array           += nzk;
89078910aadSHong Zhang     current_space->local_used      += nzk;
89178910aadSHong Zhang     current_space->local_remaining -= nzk;
89278910aadSHong Zhang 
89378910aadSHong Zhang     ui[k+1] = ui[k] + nzk;
89478910aadSHong Zhang   }
89578910aadSHong Zhang 
8966cf91177SBarry Smith #if defined(PETSC_USE_INFO)
89778910aadSHong Zhang   if (ai[mbs] != 0) {
89878910aadSHong Zhang     PetscReal af = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
899ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
900ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
901ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
90278910aadSHong Zhang   } else {
903ae15b995SBarry Smith     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
90478910aadSHong Zhang   }
90563ba0a88SBarry Smith #endif
90678910aadSHong Zhang 
90778910aadSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
90878910aadSHong Zhang   ierr = PetscFree(jl);CHKERRQ(ierr);
90978910aadSHong Zhang 
91078910aadSHong Zhang   /* destroy list of free space and other temporary array(s) */
91178910aadSHong Zhang   ierr = PetscMalloc((ui[mbs]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
912a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
91378910aadSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
91478910aadSHong Zhang 
91578910aadSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
91678910aadSHong Zhang   B    = *fact;
917ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
91878910aadSHong Zhang 
91978910aadSHong Zhang   b = (Mat_SeqSBAIJ*)B->data;
92078910aadSHong Zhang   b->singlemalloc = PETSC_FALSE;
921e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
922e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
92378910aadSHong Zhang   ierr = PetscMalloc((ui[mbs]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
92478910aadSHong Zhang   b->j    = uj;
92578910aadSHong Zhang   b->i    = ui;
92678910aadSHong Zhang   b->diag = 0;
92778910aadSHong Zhang   b->ilen = 0;
92878910aadSHong Zhang   b->imax = 0;
92978910aadSHong Zhang   b->row  = perm;
93078910aadSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
93178910aadSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
93278910aadSHong Zhang   b->icol = perm;
93378910aadSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
93478910aadSHong Zhang   ierr    = PetscMalloc((mbs+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
93552e6d16bSBarry Smith   ierr    = PetscLogObjectMemory(B,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
93678910aadSHong Zhang   b->maxnz = b->nz = ui[mbs];
93778910aadSHong Zhang 
938*5c9eb25fSBarry Smith   B->factor                 = MAT_FACTOR_CHOLESKY;
93978910aadSHong Zhang   B->info.factor_mallocs    = reallocs;
94078910aadSHong Zhang   B->info.fill_ratio_given  = fill;
94178910aadSHong Zhang   if (ai[mbs] != 0) {
94278910aadSHong Zhang     B->info.fill_ratio_needed = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
94378910aadSHong Zhang   } else {
94478910aadSHong Zhang     B->info.fill_ratio_needed = 0.0;
94578910aadSHong Zhang   }
94678910aadSHong Zhang   if (perm_identity){
9476ad2eaddSHong Zhang     B->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
9486ad2eaddSHong Zhang     B->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
9496ad2eaddSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
95078910aadSHong Zhang   } else {
9516ad2eaddSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
95278910aadSHong Zhang   }
953c05c3958SHong Zhang   PetscFunctionReturn(0);
954c05c3958SHong Zhang }
955