xref: /petsc/src/mat/impls/baij/seq/baijfact.c (revision 5d0c19d75c660d4fec594a5399ec8d8ba29c54a8)
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"
15719d5645SBarry Smith PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2(Mat B,Mat A,MatFactorInfo *info)
164eeb42bcSBarry Smith {
17719d5645SBarry 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;
21*5d0c19d7SBarry Smith   const PetscInt *r,*ic;
22*5d0c19d7SBarry Smith   PetscInt       i,j,n = a->mbs,*bi = b->i,*bj = b->j;
23690b6cddSBarry Smith   PetscInt       *ajtmpold,*ajtmp,nz,row;
24690b6cddSBarry Smith   PetscInt       *diag_offset=b->diag,idx,*ai=a->i,*aj=a->j,*pj;
25329f5518SBarry Smith   MatScalar      *pv,*v,*rtmp,m1,m2,m3,m4,*pc,*w,*x,x1,x2,x3,x4;
262515f8d2SSatish Balay   MatScalar      p1,p2,p3,p4;
273f1db9ecSBarry Smith   MatScalar      *ba = b->a,*aa = a->a;
2862bba022SBarry Smith   PetscReal      shift = info->shiftinblocks;
294eeb42bcSBarry Smith 
303a40ed3dSBarry Smith   PetscFunctionBegin;
314eeb42bcSBarry Smith   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
324eeb42bcSBarry Smith   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
33b0a32e0cSBarry Smith   ierr  = PetscMalloc(4*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
344eeb42bcSBarry Smith 
354eeb42bcSBarry Smith   for (i=0; i<n; i++) {
364078e994SBarry Smith     nz    = bi[i+1] - bi[i];
374078e994SBarry Smith     ajtmp = bj + bi[i];
384eeb42bcSBarry Smith     for  (j=0; j<nz; j++) {
394eeb42bcSBarry Smith       x = rtmp+4*ajtmp[j]; x[0] = x[1] = x[2] = x[3] = 0.0;
404eeb42bcSBarry Smith     }
414eeb42bcSBarry Smith     /* load in initial (unfactored row) */
424eeb42bcSBarry Smith     idx      = r[i];
434078e994SBarry Smith     nz       = ai[idx+1] - ai[idx];
444078e994SBarry Smith     ajtmpold = aj + ai[idx];
454078e994SBarry Smith     v        = aa + 4*ai[idx];
464eeb42bcSBarry Smith     for (j=0; j<nz; j++) {
474eeb42bcSBarry Smith       x    = rtmp+4*ic[ajtmpold[j]];
484eeb42bcSBarry Smith       x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3];
494eeb42bcSBarry Smith       v    += 4;
504eeb42bcSBarry Smith     }
514eeb42bcSBarry Smith     row = *ajtmp++;
524eeb42bcSBarry Smith     while (row < i) {
534eeb42bcSBarry Smith       pc = rtmp + 4*row;
544eeb42bcSBarry Smith       p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3];
5588685aaeSLois Curfman McInnes       if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0) {
564078e994SBarry Smith         pv = ba + 4*diag_offset[row];
574078e994SBarry Smith         pj = bj + diag_offset[row] + 1;
584eeb42bcSBarry Smith         x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
594eeb42bcSBarry Smith         pc[0] = m1 = p1*x1 + p3*x2;
604eeb42bcSBarry Smith         pc[1] = m2 = p2*x1 + p4*x2;
614eeb42bcSBarry Smith         pc[2] = m3 = p1*x3 + p3*x4;
624eeb42bcSBarry Smith         pc[3] = m4 = p2*x3 + p4*x4;
634078e994SBarry Smith         nz = bi[row+1] - diag_offset[row] - 1;
644eeb42bcSBarry Smith         pv += 4;
654eeb42bcSBarry Smith         for (j=0; j<nz; j++) {
664eeb42bcSBarry Smith           x1   = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
674eeb42bcSBarry Smith           x    = rtmp + 4*pj[j];
684eeb42bcSBarry Smith           x[0] -= m1*x1 + m3*x2;
694eeb42bcSBarry Smith           x[1] -= m2*x1 + m4*x2;
704eeb42bcSBarry Smith           x[2] -= m1*x3 + m3*x4;
714eeb42bcSBarry Smith           x[3] -= m2*x3 + m4*x4;
724eeb42bcSBarry Smith           pv   += 4;
734eeb42bcSBarry Smith         }
74efee365bSSatish Balay         ierr = PetscLogFlops(16*nz+12);CHKERRQ(ierr);
754eeb42bcSBarry Smith       }
764eeb42bcSBarry Smith       row = *ajtmp++;
774eeb42bcSBarry Smith     }
784eeb42bcSBarry Smith     /* finished row so stick it into b->a */
794078e994SBarry Smith     pv = ba + 4*bi[i];
804078e994SBarry Smith     pj = bj + bi[i];
814078e994SBarry Smith     nz = bi[i+1] - bi[i];
824eeb42bcSBarry Smith     for (j=0; j<nz; j++) {
834eeb42bcSBarry Smith       x     = rtmp+4*pj[j];
844eeb42bcSBarry Smith       pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3];
854eeb42bcSBarry Smith       pv   += 4;
864eeb42bcSBarry Smith     }
874eeb42bcSBarry Smith     /* invert diagonal block */
884078e994SBarry Smith     w = ba + 4*diag_offset[i];
8962bba022SBarry Smith     ierr = Kernel_A_gets_inverse_A_2(w,shift);CHKERRQ(ierr);
904eeb42bcSBarry Smith   }
914eeb42bcSBarry Smith 
92606d414cSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
934eeb42bcSBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
944eeb42bcSBarry Smith   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
95db4efbfdSBarry Smith   C->ops->solve          = MatSolve_SeqBAIJ_2;
96db4efbfdSBarry Smith   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2;
974eeb42bcSBarry Smith   C->assembled = PETSC_TRUE;
98efee365bSSatish Balay   ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
993a40ed3dSBarry Smith   PetscFunctionReturn(0);
1004eeb42bcSBarry Smith }
1014cd38560SBarry Smith /*
1024cd38560SBarry Smith       Version for when blocks are 2 by 2 Using natural ordering
1034cd38560SBarry Smith */
1044a2ae208SSatish Balay #undef __FUNCT__
1054a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering"
106719d5645SBarry Smith PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering(Mat C,Mat A,MatFactorInfo *info)
1074cd38560SBarry Smith {
1084cd38560SBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data;
109dfbe8321SBarry Smith   PetscErrorCode ierr;
110690b6cddSBarry Smith   PetscInt       i,j,n = a->mbs,*bi = b->i,*bj = b->j;
111690b6cddSBarry Smith   PetscInt       *ajtmpold,*ajtmp,nz,row;
112690b6cddSBarry Smith   PetscInt       *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj;
113329f5518SBarry Smith   MatScalar      *pv,*v,*rtmp,*pc,*w,*x;
1142515f8d2SSatish Balay   MatScalar      p1,p2,p3,p4,m1,m2,m3,m4,x1,x2,x3,x4;
1154cd38560SBarry Smith   MatScalar      *ba = b->a,*aa = a->a;
11662bba022SBarry Smith   PetscReal      shift = info->shiftinblocks;
1174cd38560SBarry Smith 
1184cd38560SBarry Smith   PetscFunctionBegin;
119b0a32e0cSBarry Smith   ierr = PetscMalloc(4*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
1204cd38560SBarry Smith 
1214cd38560SBarry Smith   for (i=0; i<n; i++) {
1224cd38560SBarry Smith     nz    = bi[i+1] - bi[i];
1234cd38560SBarry Smith     ajtmp = bj + bi[i];
1244cd38560SBarry Smith     for  (j=0; j<nz; j++) {
1254cd38560SBarry Smith       x = rtmp+4*ajtmp[j];
1264cd38560SBarry Smith       x[0]  = x[1]  = x[2]  = x[3]  = 0.0;
1274cd38560SBarry Smith     }
1284cd38560SBarry Smith     /* load in initial (unfactored row) */
1294cd38560SBarry Smith     nz       = ai[i+1] - ai[i];
1304cd38560SBarry Smith     ajtmpold = aj + ai[i];
1314cd38560SBarry Smith     v        = aa + 4*ai[i];
1324cd38560SBarry Smith     for (j=0; j<nz; j++) {
1334cd38560SBarry Smith       x    = rtmp+4*ajtmpold[j];
1344cd38560SBarry Smith       x[0]  = v[0];  x[1]  = v[1];  x[2]  = v[2];  x[3]  = v[3];
1354cd38560SBarry Smith       v    += 4;
1364cd38560SBarry Smith     }
1374cd38560SBarry Smith     row = *ajtmp++;
1384cd38560SBarry Smith     while (row < i) {
1394cd38560SBarry Smith       pc  = rtmp + 4*row;
1404cd38560SBarry Smith       p1  = pc[0];  p2  = pc[1];  p3  = pc[2];  p4  = pc[3];
1414cd38560SBarry Smith       if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0) {
1424cd38560SBarry Smith         pv = ba + 4*diag_offset[row];
1434cd38560SBarry Smith         pj = bj + diag_offset[row] + 1;
1444cd38560SBarry Smith         x1  = pv[0];  x2  = pv[1];  x3  = pv[2];  x4  = pv[3];
1454cd38560SBarry Smith         pc[0] = m1 = p1*x1 + p3*x2;
1464cd38560SBarry Smith         pc[1] = m2 = p2*x1 + p4*x2;
1474cd38560SBarry Smith         pc[2] = m3 = p1*x3 + p3*x4;
1484cd38560SBarry Smith         pc[3] = m4 = p2*x3 + p4*x4;
1494cd38560SBarry Smith         nz = bi[row+1] - diag_offset[row] - 1;
1504cd38560SBarry Smith         pv += 4;
1514cd38560SBarry Smith         for (j=0; j<nz; j++) {
1524cd38560SBarry Smith           x1   = pv[0];  x2  = pv[1];   x3 = pv[2];  x4  = pv[3];
1534cd38560SBarry Smith           x    = rtmp + 4*pj[j];
1544cd38560SBarry Smith           x[0] -= m1*x1 + m3*x2;
1554cd38560SBarry Smith           x[1] -= m2*x1 + m4*x2;
1564cd38560SBarry Smith           x[2] -= m1*x3 + m3*x4;
1574cd38560SBarry Smith           x[3] -= m2*x3 + m4*x4;
1584cd38560SBarry Smith           pv   += 4;
1594cd38560SBarry Smith         }
160efee365bSSatish Balay         ierr = PetscLogFlops(16*nz+12);CHKERRQ(ierr);
1614cd38560SBarry Smith       }
1624cd38560SBarry Smith       row = *ajtmp++;
1634cd38560SBarry Smith     }
1644cd38560SBarry Smith     /* finished row so stick it into b->a */
1654cd38560SBarry Smith     pv = ba + 4*bi[i];
1664cd38560SBarry Smith     pj = bj + bi[i];
1674cd38560SBarry Smith     nz = bi[i+1] - bi[i];
1684cd38560SBarry Smith     for (j=0; j<nz; j++) {
1694cd38560SBarry Smith       x      = rtmp+4*pj[j];
1704cd38560SBarry Smith       pv[0]  = x[0];  pv[1]  = x[1];  pv[2]  = x[2];  pv[3]  = x[3];
1714cd38560SBarry Smith       pv   += 4;
1724cd38560SBarry Smith     }
1734cd38560SBarry Smith     /* invert diagonal block */
1744cd38560SBarry Smith     w = ba + 4*diag_offset[i];
17562bba022SBarry Smith     ierr = Kernel_A_gets_inverse_A_2(w,shift);CHKERRQ(ierr);
1764cd38560SBarry Smith   }
1774cd38560SBarry Smith 
178606d414cSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
179db4efbfdSBarry Smith   C->ops->solve          = MatSolve_SeqBAIJ_2_NaturalOrdering;
180db4efbfdSBarry Smith   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_NaturalOrdering;
1814cd38560SBarry Smith   C->assembled = PETSC_TRUE;
182efee365bSSatish Balay   ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
1834cd38560SBarry Smith   PetscFunctionReturn(0);
1844cd38560SBarry Smith }
1857fc0212eSBarry Smith 
1867fc0212eSBarry Smith /* ----------------------------------------------------------- */
1877fc0212eSBarry Smith /*
1887fc0212eSBarry Smith      Version for when blocks are 1 by 1.
1897fc0212eSBarry Smith */
1904a2ae208SSatish Balay #undef __FUNCT__
1914a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_1"
192719d5645SBarry Smith PetscErrorCode MatLUFactorNumeric_SeqBAIJ_1(Mat C,Mat A,MatFactorInfo *info)
1937fc0212eSBarry Smith {
1947fc0212eSBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data;
1957cdf2dbbSSatish Balay   IS             isrow = b->row,isicol = b->icol;
1966849ba73SBarry Smith   PetscErrorCode ierr;
197*5d0c19d7SBarry Smith   const PetscInt *r,*ic;
198*5d0c19d7SBarry Smith   PetscInt       i,j,n = a->mbs,*bi = b->i,*bj = b->j;
199690b6cddSBarry Smith   PetscInt       *ajtmpold,*ajtmp,nz,row,*ai = a->i,*aj = a->j;
200690b6cddSBarry Smith   PetscInt       *diag_offset = b->diag,diag,*pj;
201329f5518SBarry Smith   MatScalar      *pv,*v,*rtmp,multiplier,*pc;
2023f1db9ecSBarry Smith   MatScalar      *ba = b->a,*aa = a->a;
2037fc0212eSBarry Smith 
2043a40ed3dSBarry Smith   PetscFunctionBegin;
2057fc0212eSBarry Smith   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
2067fc0212eSBarry Smith   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
207b0a32e0cSBarry Smith   ierr  = PetscMalloc((n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
2087fc0212eSBarry Smith 
2097fc0212eSBarry Smith   for (i=0; i<n; i++) {
2104078e994SBarry Smith     nz    = bi[i+1] - bi[i];
2114078e994SBarry Smith     ajtmp = bj + bi[i];
2127fc0212eSBarry Smith     for  (j=0; j<nz; j++) rtmp[ajtmp[j]] = 0.0;
2137fc0212eSBarry Smith 
2147fc0212eSBarry Smith     /* load in initial (unfactored row) */
2154078e994SBarry Smith     nz       = ai[r[i]+1] - ai[r[i]];
2164078e994SBarry Smith     ajtmpold = aj + ai[r[i]];
2174078e994SBarry Smith     v        = aa + ai[r[i]];
2187fc0212eSBarry Smith     for (j=0; j<nz; j++) rtmp[ic[ajtmpold[j]]] =  v[j];
2197fc0212eSBarry Smith 
2207fc0212eSBarry Smith     row = *ajtmp++;
2217fc0212eSBarry Smith     while (row < i) {
2227fc0212eSBarry Smith       pc = rtmp + row;
2237fc0212eSBarry Smith       if (*pc != 0.0) {
2244078e994SBarry Smith         pv         = ba + diag_offset[row];
2254078e994SBarry Smith         pj         = bj + diag_offset[row] + 1;
2267fc0212eSBarry Smith         multiplier = *pc * *pv++;
2277fc0212eSBarry Smith         *pc        = multiplier;
2284078e994SBarry Smith         nz         = bi[row+1] - diag_offset[row] - 1;
2297fc0212eSBarry Smith         for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j];
230efee365bSSatish Balay         ierr = PetscLogFlops(1+2*nz);CHKERRQ(ierr);
2317fc0212eSBarry Smith       }
2327fc0212eSBarry Smith       row = *ajtmp++;
2337fc0212eSBarry Smith     }
2347fc0212eSBarry Smith     /* finished row so stick it into b->a */
2354078e994SBarry Smith     pv = ba + bi[i];
2364078e994SBarry Smith     pj = bj + bi[i];
2374078e994SBarry Smith     nz = bi[i+1] - bi[i];
2387fc0212eSBarry Smith     for (j=0; j<nz; j++) {pv[j] = rtmp[pj[j]];}
2394078e994SBarry Smith     diag = diag_offset[i] - bi[i];
2407fc0212eSBarry Smith     /* check pivot entry for current row */
241a73cf429SBarry Smith     if (pv[diag] == 0.0) {
24229bbc08cSBarry Smith       SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot");
2437fc0212eSBarry Smith     }
2447fc0212eSBarry Smith     pv[diag] = 1.0/pv[diag];
2457fc0212eSBarry Smith   }
2467fc0212eSBarry Smith 
247606d414cSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
2487fc0212eSBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
2497fc0212eSBarry Smith   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
250db4efbfdSBarry Smith   C->ops->solve          = MatSolve_SeqBAIJ_1;
251db4efbfdSBarry Smith   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1;
2527fc0212eSBarry Smith   C->assembled = PETSC_TRUE;
253d0f46423SBarry Smith   ierr = PetscLogFlops(C->cmap->n);CHKERRQ(ierr);
2543a40ed3dSBarry Smith   PetscFunctionReturn(0);
2557fc0212eSBarry Smith }
2567fc0212eSBarry Smith 
257e631078cSBarry Smith EXTERN_C_BEGIN
258b24902e0SBarry Smith #undef __FUNCT__
259b24902e0SBarry Smith #define __FUNCT__ "MatGetFactor_seqbaij_petsc"
260b24902e0SBarry Smith PetscErrorCode MatGetFactor_seqbaij_petsc(Mat A,MatFactorType ftype,Mat *B)
261b24902e0SBarry Smith {
262d0f46423SBarry Smith   PetscInt           n = A->rmap->n;
263b24902e0SBarry Smith   PetscErrorCode     ierr;
264b24902e0SBarry Smith 
265b24902e0SBarry Smith   PetscFunctionBegin;
266b24902e0SBarry Smith   ierr = MatCreate(((PetscObject)A)->comm,B);CHKERRQ(ierr);
267b24902e0SBarry Smith   ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr);
268b24902e0SBarry Smith   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) {
269b24902e0SBarry Smith     ierr = MatSetType(*B,MATSEQBAIJ);CHKERRQ(ierr);
270db4efbfdSBarry Smith     (*B)->ops->lufactorsymbolic  = MatLUFactorSymbolic_SeqBAIJ;
271db4efbfdSBarry Smith     (*B)->ops->ilufactorsymbolic = MatILUFactorSymbolic_SeqBAIJ;
272b24902e0SBarry Smith   } else if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) {
273b24902e0SBarry Smith     ierr = MatSetType(*B,MATSEQSBAIJ);CHKERRQ(ierr);
2745c9eb25fSBarry Smith     ierr = MatSeqSBAIJSetPreallocation(*B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
2755c9eb25fSBarry Smith     (*B)->ops->iccfactorsymbolic      = MatICCFactorSymbolic_SeqBAIJ;
2765c9eb25fSBarry Smith     (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqBAIJ;
277b24902e0SBarry Smith   } else SETERRQ(PETSC_ERR_SUP,"Factor type not supported");
278719d5645SBarry Smith   (*B)->factor = ftype;
279b24902e0SBarry Smith   PetscFunctionReturn(0);
280b24902e0SBarry Smith }
281e631078cSBarry Smith EXTERN_C_END
282273d9f13SBarry Smith 
283db4efbfdSBarry Smith EXTERN_C_BEGIN
284db4efbfdSBarry Smith #undef __FUNCT__
285db4efbfdSBarry Smith #define __FUNCT__ "MatGetFactorAvailable_seqaij_petsc"
286db4efbfdSBarry Smith PetscErrorCode MatGetFactorAvailable_seqbaij_petsc(Mat A,MatFactorType ftype,PetscTruth *flg)
287db4efbfdSBarry Smith {
288db4efbfdSBarry Smith   PetscFunctionBegin;
289db4efbfdSBarry Smith   *flg = PETSC_TRUE;
290db4efbfdSBarry Smith   PetscFunctionReturn(0);
291db4efbfdSBarry Smith }
292db4efbfdSBarry Smith EXTERN_C_END
293db4efbfdSBarry Smith 
2945d517e7eSBarry Smith /* ----------------------------------------------------------- */
2954a2ae208SSatish Balay #undef __FUNCT__
2964a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqBAIJ"
297dfbe8321SBarry Smith PetscErrorCode MatLUFactor_SeqBAIJ(Mat A,IS row,IS col,MatFactorInfo *info)
2985d517e7eSBarry Smith {
299dfbe8321SBarry Smith   PetscErrorCode ierr;
3005d517e7eSBarry Smith   Mat            C;
3015d517e7eSBarry Smith 
3023a40ed3dSBarry Smith   PetscFunctionBegin;
30343244d56SBarry Smith   ierr = MatGetFactor(A,MAT_SOLVER_PETSC,MAT_FACTOR_LU,&C);CHKERRQ(ierr);
304719d5645SBarry Smith   ierr = MatLUFactorSymbolic(C,A,row,col,info);CHKERRQ(ierr);
305719d5645SBarry Smith   ierr = MatLUFactorNumeric(C,A,info);CHKERRQ(ierr);
306db4efbfdSBarry Smith   A->ops->solve            = C->ops->solve;
307db4efbfdSBarry Smith   A->ops->solvetranspose   = C->ops->solvetranspose;
308273d9f13SBarry Smith   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
30952e6d16bSBarry Smith   ierr = PetscLogObjectParent(A,((Mat_SeqBAIJ*)(A->data))->icol);CHKERRQ(ierr);
3103a40ed3dSBarry Smith   PetscFunctionReturn(0);
3115d517e7eSBarry Smith }
3124cd38560SBarry Smith 
313c05c3958SHong Zhang #include "src/mat/impls/sbaij/seq/sbaij.h"
314c05c3958SHong Zhang #undef __FUNCT__
315c05c3958SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N"
316719d5645SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N(Mat C,Mat A,MatFactorInfo *info)
317c05c3958SHong Zhang {
318f3086b4bSHong Zhang   PetscErrorCode ierr;
31978910aadSHong Zhang   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data;
32078910aadSHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
32178910aadSHong Zhang   IS             ip=b->row;
322*5d0c19d7SBarry Smith   const PetscInt *rip;
323*5d0c19d7SBarry Smith   PetscInt       i,j,mbs=a->mbs,bs=A->rmap->bs,*bi=b->i,*bj=b->j,*bcol;
32478910aadSHong Zhang   PetscInt       *ai=a->i,*aj=a->j;
32578910aadSHong Zhang   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
32678910aadSHong Zhang   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
3273cea4cbeSHong Zhang   PetscReal      zeropivot,rs,shiftnz;
328fbf22428SSatish Balay   PetscReal      shiftpd;
3293cea4cbeSHong Zhang   ChShift_Ctx    sctx;
3303cea4cbeSHong Zhang   PetscInt       newshift;
33178910aadSHong Zhang 
332c05c3958SHong Zhang   PetscFunctionBegin;
3336ad2eaddSHong Zhang   if (bs > 1) {
3346ad2eaddSHong Zhang     if (!a->sbaijMat){
335ceb03754SKris Buschelman       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
3366ad2eaddSHong Zhang     }
337719d5645SBarry Smith     ierr = (a->sbaijMat)->ops->choleskyfactornumeric(C,a->sbaijMat,info);CHKERRQ(ierr);
3386ad2eaddSHong Zhang     ierr = MatDestroy(a->sbaijMat);CHKERRQ(ierr);
3396ad2eaddSHong Zhang     a->sbaijMat = PETSC_NULL;
3406ad2eaddSHong Zhang     PetscFunctionReturn(0);
3416ad2eaddSHong Zhang   }
34278910aadSHong Zhang 
34378910aadSHong Zhang   /* initialization */
3443cea4cbeSHong Zhang   shiftnz   = info->shiftnz;
3453cea4cbeSHong Zhang   shiftpd   = info->shiftpd;
3463cea4cbeSHong Zhang   zeropivot = info->zeropivot;
3473cea4cbeSHong Zhang 
3486ad2eaddSHong Zhang   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
34978910aadSHong Zhang   nz   = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar);
35078910aadSHong Zhang   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
35178910aadSHong Zhang   jl   = il + mbs;
35278910aadSHong Zhang   rtmp = (MatScalar*)(jl + mbs);
35378910aadSHong Zhang 
3543cea4cbeSHong Zhang   sctx.shift_amount = 0;
3553cea4cbeSHong Zhang   sctx.nshift       = 0;
35678910aadSHong Zhang   do {
3573cea4cbeSHong Zhang     sctx.chshift = PETSC_FALSE;
35878910aadSHong Zhang     for (i=0; i<mbs; i++) {
35978910aadSHong Zhang       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
36078910aadSHong Zhang     }
36178910aadSHong Zhang 
36278910aadSHong Zhang     for (k = 0; k<mbs; k++){
36378910aadSHong Zhang       bval = ba + bi[k];
36478910aadSHong Zhang       /* initialize k-th row by the perm[k]-th row of A */
36578910aadSHong Zhang       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
36678910aadSHong Zhang       for (j = jmin; j < jmax; j++){
36778910aadSHong Zhang         col = rip[aj[j]];
36878910aadSHong Zhang         if (col >= k){ /* only take upper triangular entry */
36978910aadSHong Zhang           rtmp[col] = aa[j];
37078910aadSHong Zhang           *bval++  = 0.0; /* for in-place factorization */
37178910aadSHong Zhang         }
37278910aadSHong Zhang       }
3733cea4cbeSHong Zhang 
3743cea4cbeSHong Zhang       /* shift the diagonal of the matrix */
3753cea4cbeSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
37678910aadSHong Zhang 
37778910aadSHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
37878910aadSHong Zhang       dk = rtmp[k];
37978910aadSHong Zhang       i = jl[k]; /* first row to be added to k_th row  */
38078910aadSHong Zhang 
38178910aadSHong Zhang       while (i < k){
38278910aadSHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
38378910aadSHong Zhang 
38478910aadSHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
38578910aadSHong Zhang         ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
38678910aadSHong Zhang         uikdi = - ba[ili]*ba[bi[i]];  /* diagonal(k) */
38778910aadSHong Zhang         dk += uikdi*ba[ili];
38878910aadSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
38978910aadSHong Zhang 
39078910aadSHong Zhang         /* add multiple of row i to k-th row */
39178910aadSHong Zhang         jmin = ili + 1; jmax = bi[i+1];
39278910aadSHong Zhang         if (jmin < jmax){
39378910aadSHong Zhang           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
39478910aadSHong Zhang           /* update il and jl for row i */
39578910aadSHong Zhang           il[i] = jmin;
39678910aadSHong Zhang           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
39778910aadSHong Zhang         }
39878910aadSHong Zhang         i = nexti;
39978910aadSHong Zhang       }
40078910aadSHong Zhang 
4013cea4cbeSHong Zhang       /* shift the diagonals when zero pivot is detected */
4023cea4cbeSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
4033cea4cbeSHong Zhang       rs   = 0.0;
40478910aadSHong Zhang       jmin = bi[k]+1;
40578910aadSHong Zhang       nz   = bi[k+1] - jmin;
40678910aadSHong Zhang       if (nz){
40778910aadSHong Zhang         bcol = bj + jmin;
40878910aadSHong Zhang         while (nz--){
40989f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
41089f845c8SHong Zhang           bcol++;
41178910aadSHong Zhang         }
41278910aadSHong Zhang       }
4133cea4cbeSHong Zhang 
4143cea4cbeSHong Zhang       sctx.rs = rs;
4153cea4cbeSHong Zhang       sctx.pv = dk;
41645938b79SHong Zhang       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
41745938b79SHong Zhang       if (newshift == 1) break;
41878910aadSHong Zhang 
41978910aadSHong Zhang       /* copy data into U(k,:) */
42078910aadSHong Zhang       ba[bi[k]] = 1.0/dk; /* U(k,k) */
42178910aadSHong Zhang       jmin = bi[k]+1; jmax = bi[k+1];
42278910aadSHong Zhang       if (jmin < jmax) {
42378910aadSHong Zhang         for (j=jmin; j<jmax; j++){
42478910aadSHong Zhang           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
42578910aadSHong Zhang         }
42678910aadSHong Zhang         /* add the k-th row into il and jl */
42778910aadSHong Zhang         il[k] = jmin;
42878910aadSHong Zhang         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
42978910aadSHong Zhang       }
43078910aadSHong Zhang     }
4313cea4cbeSHong Zhang   } while (sctx.chshift);
43278910aadSHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
43378910aadSHong Zhang 
43478910aadSHong Zhang   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
43578910aadSHong Zhang   C->assembled    = PETSC_TRUE;
43678910aadSHong Zhang   C->preallocated = PETSC_TRUE;
437d0f46423SBarry Smith   ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr);
4383cea4cbeSHong Zhang   if (sctx.nshift){
4393cea4cbeSHong Zhang     if (shiftnz) {
4401e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
4413cea4cbeSHong Zhang     } else if (shiftpd) {
4421e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
44378910aadSHong Zhang     }
44478910aadSHong Zhang   }
445c05c3958SHong Zhang   PetscFunctionReturn(0);
446c05c3958SHong Zhang }
4474cd38560SBarry Smith 
448c05c3958SHong Zhang #undef __FUNCT__
449c05c3958SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering"
450719d5645SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering(Mat C,Mat A,MatFactorInfo *info)
451c05c3958SHong Zhang {
45278910aadSHong Zhang   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data;
45378910aadSHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
45478910aadSHong Zhang   PetscErrorCode ierr;
4553cea4cbeSHong Zhang   PetscInt       i,j,am=a->mbs;
45678910aadSHong Zhang   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
4573cea4cbeSHong Zhang   PetscInt       k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz;
45878910aadSHong Zhang   MatScalar      *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval;
4593cea4cbeSHong Zhang   PetscReal      zeropivot,rs,shiftnz;
460fbf22428SSatish Balay   PetscReal      shiftpd;
4613cea4cbeSHong Zhang   ChShift_Ctx    sctx;
4623cea4cbeSHong Zhang   PetscInt       newshift;
46378910aadSHong Zhang 
464c05c3958SHong Zhang   PetscFunctionBegin;
46578910aadSHong Zhang   /* initialization */
4663cea4cbeSHong Zhang   shiftnz   = info->shiftnz;
4673cea4cbeSHong Zhang   shiftpd   = info->shiftpd;
4683cea4cbeSHong Zhang   zeropivot = info->zeropivot;
4693cea4cbeSHong Zhang 
47078910aadSHong Zhang   nz   = (2*am+1)*sizeof(PetscInt)+am*sizeof(MatScalar);
47178910aadSHong Zhang   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
47278910aadSHong Zhang   jl   = il + am;
47378910aadSHong Zhang   rtmp = (MatScalar*)(jl + am);
47478910aadSHong Zhang 
4753cea4cbeSHong Zhang   sctx.shift_amount = 0;
4763cea4cbeSHong Zhang   sctx.nshift       = 0;
47778910aadSHong Zhang   do {
4783cea4cbeSHong Zhang     sctx.chshift = PETSC_FALSE;
47978910aadSHong Zhang     for (i=0; i<am; i++) {
48078910aadSHong Zhang       rtmp[i] = 0.0; jl[i] = am; il[0] = 0;
48178910aadSHong Zhang     }
48278910aadSHong Zhang 
48378910aadSHong Zhang     for (k = 0; k<am; k++){
48478910aadSHong Zhang     /* initialize k-th row with elements nonzero in row perm(k) of A */
48578910aadSHong Zhang       nz   = ai[k+1] - ai[k];
48678910aadSHong Zhang       acol = aj + ai[k];
48778910aadSHong Zhang       aval = aa + ai[k];
48878910aadSHong Zhang       bval = ba + bi[k];
48978910aadSHong Zhang       while (nz -- ){
49078910aadSHong Zhang         if (*acol < k) { /* skip lower triangular entries */
49178910aadSHong Zhang           acol++; aval++;
49278910aadSHong Zhang         } else {
49378910aadSHong Zhang           rtmp[*acol++] = *aval++;
49478910aadSHong Zhang           *bval++       = 0.0; /* for in-place factorization */
49578910aadSHong Zhang         }
49678910aadSHong Zhang       }
4973cea4cbeSHong Zhang 
4983cea4cbeSHong Zhang       /* shift the diagonal of the matrix */
4993cea4cbeSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
50078910aadSHong Zhang 
50178910aadSHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
50278910aadSHong Zhang       dk = rtmp[k];
50378910aadSHong Zhang       i  = jl[k]; /* first row to be added to k_th row  */
50478910aadSHong Zhang 
50578910aadSHong Zhang       while (i < k){
50678910aadSHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
50778910aadSHong Zhang         /* compute multiplier, update D(k) and U(i,k) */
50878910aadSHong Zhang         ili   = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
50978910aadSHong Zhang         uikdi = - ba[ili]*ba[bi[i]];
51078910aadSHong Zhang         dk   += uikdi*ba[ili];
51178910aadSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
51278910aadSHong Zhang 
51378910aadSHong Zhang         /* add multiple of row i to k-th row ... */
51478910aadSHong Zhang         jmin = ili + 1;
51578910aadSHong Zhang         nz   = bi[i+1] - jmin;
51678910aadSHong Zhang         if (nz > 0){
51778910aadSHong Zhang           bcol = bj + jmin;
51878910aadSHong Zhang           bval = ba + jmin;
51978910aadSHong Zhang           while (nz --) rtmp[*bcol++] += uikdi*(*bval++);
52078910aadSHong Zhang           /* update il and jl for i-th row */
52178910aadSHong Zhang           il[i] = jmin;
52278910aadSHong Zhang           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
52378910aadSHong Zhang         }
52478910aadSHong Zhang         i = nexti;
52578910aadSHong Zhang       }
52678910aadSHong Zhang 
5273cea4cbeSHong Zhang       /* shift the diagonals when zero pivot is detected */
5283cea4cbeSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
5293cea4cbeSHong Zhang       rs   = 0.0;
53078910aadSHong Zhang       jmin = bi[k]+1;
53178910aadSHong Zhang       nz   = bi[k+1] - jmin;
53278910aadSHong Zhang       if (nz){
53378910aadSHong Zhang         bcol = bj + jmin;
53478910aadSHong Zhang         while (nz--){
53589f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
53689f845c8SHong Zhang           bcol++;
53778910aadSHong Zhang         }
53878910aadSHong Zhang       }
5393cea4cbeSHong Zhang 
5403cea4cbeSHong Zhang       sctx.rs = rs;
5413cea4cbeSHong Zhang       sctx.pv = dk;
54245938b79SHong Zhang       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
54345938b79SHong Zhang       if (newshift == 1) break;    /* sctx.shift_amount is updated */
54478910aadSHong Zhang 
54578910aadSHong Zhang       /* copy data into U(k,:) */
54678910aadSHong Zhang       ba[bi[k]] = 1.0/dk;
54778910aadSHong Zhang       jmin      = bi[k]+1;
54878910aadSHong Zhang       nz        = bi[k+1] - jmin;
54978910aadSHong Zhang       if (nz){
55078910aadSHong Zhang         bcol = bj + jmin;
55178910aadSHong Zhang         bval = ba + jmin;
55278910aadSHong Zhang         while (nz--){
55378910aadSHong Zhang           *bval++       = rtmp[*bcol];
55478910aadSHong Zhang           rtmp[*bcol++] = 0.0;
55578910aadSHong Zhang         }
55678910aadSHong Zhang         /* add k-th row into il and jl */
55778910aadSHong Zhang         il[k] = jmin;
55878910aadSHong Zhang         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
55978910aadSHong Zhang       }
56078910aadSHong Zhang     }
5613cea4cbeSHong Zhang   } while (sctx.chshift);
56278910aadSHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
56378910aadSHong Zhang 
564719d5645SBarry Smith   C->ops->solve                 = MatSolve_SeqSBAIJ_1_NaturalOrdering;
565719d5645SBarry Smith   C->ops->solvetranspose        = MatSolve_SeqSBAIJ_1_NaturalOrdering;
56678910aadSHong Zhang   C->assembled    = PETSC_TRUE;
56778910aadSHong Zhang   C->preallocated = PETSC_TRUE;
568d0f46423SBarry Smith   ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr);
5693cea4cbeSHong Zhang     if (sctx.nshift){
5703cea4cbeSHong Zhang     if (shiftnz) {
5711e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
5723cea4cbeSHong Zhang     } else if (shiftpd) {
5731e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
57478910aadSHong Zhang     }
57578910aadSHong Zhang   }
576c05c3958SHong Zhang   PetscFunctionReturn(0);
577c05c3958SHong Zhang }
578c05c3958SHong Zhang 
579c05c3958SHong Zhang #include "petscbt.h"
580c05c3958SHong Zhang #include "src/mat/utils/freespace.h"
581c05c3958SHong Zhang #undef __FUNCT__
582c05c3958SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqBAIJ"
583719d5645SBarry Smith PetscErrorCode MatICCFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS perm,MatFactorInfo *info)
584c05c3958SHong Zhang {
58578910aadSHong Zhang   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data;
58678910aadSHong Zhang   Mat_SeqSBAIJ       *b;
58778910aadSHong Zhang   Mat                B;
58878910aadSHong Zhang   PetscErrorCode     ierr;
58978910aadSHong Zhang   PetscTruth         perm_identity;
590*5d0c19d7SBarry Smith   PetscInt           reallocs=0,i,*ai=a->i,*aj=a->j,am=a->mbs,bs=A->rmap->bs,*ui;
591*5d0c19d7SBarry Smith   const PetscInt     *rip;
59278910aadSHong Zhang   PetscInt           jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow;
593cfdb8b8aSHong Zhang   PetscInt           nlnk,*lnk,*lnk_lvl=PETSC_NULL,ncols,ncols_upper,*cols,*cols_lvl,*uj,**uj_ptr,**uj_lvl_ptr;
59478910aadSHong Zhang   PetscReal          fill=info->fill,levels=info->levels;
595a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
596a1a86e44SBarry Smith   PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
59778910aadSHong Zhang   PetscBT            lnkbt;
59878910aadSHong Zhang 
599c05c3958SHong Zhang   PetscFunctionBegin;
6006ad2eaddSHong Zhang   if (bs > 1){
6016ad2eaddSHong Zhang     if (!a->sbaijMat){
602ceb03754SKris Buschelman       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
6036ad2eaddSHong Zhang     }
604719d5645SBarry Smith     (fact)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqSBAIJ;  /* undue the change made in MatGetFactor_seqbaij_petsc */
605719d5645SBarry Smith     ierr = MatICCFactorSymbolic(fact,a->sbaijMat,perm,info);CHKERRQ(ierr);
6066ad2eaddSHong Zhang     PetscFunctionReturn(0);
6076ad2eaddSHong Zhang   }
6086ad2eaddSHong Zhang 
60978910aadSHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
61078910aadSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
61178910aadSHong Zhang 
61278910aadSHong Zhang   /* special case that simply copies fill pattern */
61378910aadSHong Zhang   if (!levels && perm_identity) {
61478910aadSHong Zhang     ierr = MatMarkDiagonal_SeqBAIJ(A);CHKERRQ(ierr);
61578910aadSHong Zhang     ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
61678910aadSHong Zhang     for (i=0; i<am; i++) {
61778910aadSHong Zhang       ui[i] = ai[i+1] - a->diag[i]; /* ui: rowlengths - changes when !perm_identity */
61878910aadSHong Zhang     }
619719d5645SBarry Smith     B = fact;
62078910aadSHong Zhang     ierr = MatSeqSBAIJSetPreallocation(B,1,0,ui);CHKERRQ(ierr);
62178910aadSHong Zhang 
622b24902e0SBarry Smith 
62378910aadSHong Zhang     b  = (Mat_SeqSBAIJ*)B->data;
62478910aadSHong Zhang     uj = b->j;
62578910aadSHong Zhang     for (i=0; i<am; i++) {
62678910aadSHong Zhang       aj = a->j + a->diag[i];
62778910aadSHong Zhang       for (j=0; j<ui[i]; j++){
62878910aadSHong Zhang         *uj++ = *aj++;
62978910aadSHong Zhang       }
63078910aadSHong Zhang       b->ilen[i] = ui[i];
63178910aadSHong Zhang     }
63278910aadSHong Zhang     ierr = PetscFree(ui);CHKERRQ(ierr);
633719d5645SBarry Smith     B->factor = MAT_FACTOR_NONE;
63478910aadSHong Zhang     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
63578910aadSHong Zhang     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
636719d5645SBarry Smith     B->factor = MAT_FACTOR_ICC;
63778910aadSHong Zhang 
63878910aadSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
63978910aadSHong Zhang     PetscFunctionReturn(0);
64078910aadSHong Zhang   }
64178910aadSHong Zhang 
64278910aadSHong Zhang   /* initialization */
64378910aadSHong Zhang   ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
64478910aadSHong Zhang   ui[0] = 0;
64578910aadSHong Zhang   ierr  = PetscMalloc((2*am+1)*sizeof(PetscInt),&cols_lvl);CHKERRQ(ierr);
64678910aadSHong Zhang 
64778910aadSHong Zhang   /* jl: linked list for storing indices of the pivot rows
64878910aadSHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
64978910aadSHong Zhang   ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt*),&jl);CHKERRQ(ierr);
65078910aadSHong Zhang   il         = jl + am;
65178910aadSHong Zhang   uj_ptr     = (PetscInt**)(il + am);
65278910aadSHong Zhang   uj_lvl_ptr = (PetscInt**)(uj_ptr + am);
65378910aadSHong Zhang   for (i=0; i<am; i++){
65478910aadSHong Zhang     jl[i] = am; il[i] = 0;
65578910aadSHong Zhang   }
65678910aadSHong Zhang 
65778910aadSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
65878910aadSHong Zhang   nlnk = am + 1;
65978910aadSHong Zhang   ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
66078910aadSHong Zhang 
66178910aadSHong Zhang   /* initial FreeSpace size is fill*(ai[am]+1) */
662a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
66378910aadSHong Zhang   current_space = free_space;
664a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr);
66578910aadSHong Zhang   current_space_lvl = free_space_lvl;
66678910aadSHong Zhang 
66778910aadSHong Zhang   for (k=0; k<am; k++){  /* for each active row k */
66878910aadSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
66978910aadSHong Zhang     nzk   = 0;
67078910aadSHong Zhang     ncols = ai[rip[k]+1] - ai[rip[k]];
67178910aadSHong Zhang     ncols_upper = 0;
67278910aadSHong Zhang     cols        = cols_lvl + am;
67378910aadSHong Zhang     for (j=0; j<ncols; j++){
67478910aadSHong Zhang       i = rip[*(aj + ai[rip[k]] + j)];
67578910aadSHong Zhang       if (i >= k){ /* only take upper triangular entry */
67678910aadSHong Zhang         cols[ncols_upper] = i;
67778910aadSHong Zhang         cols_lvl[ncols_upper] = -1;  /* initialize level for nonzero entries */
67878910aadSHong Zhang         ncols_upper++;
67978910aadSHong Zhang       }
68078910aadSHong Zhang     }
68178910aadSHong Zhang     ierr = PetscIncompleteLLAdd(ncols_upper,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
68278910aadSHong Zhang     nzk += nlnk;
68378910aadSHong Zhang 
68478910aadSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
68578910aadSHong Zhang     prow = jl[k]; /* 1st pivot row */
68678910aadSHong Zhang 
68778910aadSHong Zhang     while (prow < k){
68878910aadSHong Zhang       nextprow = jl[prow];
68978910aadSHong Zhang 
69078910aadSHong Zhang       /* merge prow into k-th row */
69178910aadSHong Zhang       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
69278910aadSHong Zhang       jmax = ui[prow+1];
69378910aadSHong Zhang       ncols = jmax-jmin;
69478910aadSHong Zhang       i     = jmin - ui[prow];
69578910aadSHong Zhang       cols = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
69678910aadSHong Zhang       for (j=0; j<ncols; j++) cols_lvl[j] = *(uj_lvl_ptr[prow] + i + j);
6975a8e39fbSHong Zhang       ierr = PetscIncompleteLLAddSorted(ncols,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
69878910aadSHong Zhang       nzk += nlnk;
69978910aadSHong Zhang 
70078910aadSHong Zhang       /* update il and jl for prow */
70178910aadSHong Zhang       if (jmin < jmax){
70278910aadSHong Zhang         il[prow] = jmin;
70378910aadSHong Zhang         j = *cols; jl[prow] = jl[j]; jl[j] = prow;
70478910aadSHong Zhang       }
70578910aadSHong Zhang       prow = nextprow;
70678910aadSHong Zhang     }
70778910aadSHong Zhang 
70878910aadSHong Zhang     /* if free space is not available, make more free space */
70978910aadSHong Zhang     if (current_space->local_remaining<nzk) {
71078910aadSHong Zhang       i = am - k + 1; /* num of unfactored rows */
71178910aadSHong Zhang       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
712a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
713a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(i,&current_space_lvl);CHKERRQ(ierr);
71478910aadSHong Zhang       reallocs++;
71578910aadSHong Zhang     }
71678910aadSHong Zhang 
71778910aadSHong Zhang     /* copy data into free_space and free_space_lvl, then initialize lnk */
71878910aadSHong Zhang     ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
71978910aadSHong Zhang 
72078910aadSHong Zhang     /* add the k-th row into il and jl */
72178910aadSHong Zhang     if (nzk-1 > 0){
72278910aadSHong Zhang       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
72378910aadSHong Zhang       jl[k] = jl[i]; jl[i] = k;
72478910aadSHong Zhang       il[k] = ui[k] + 1;
72578910aadSHong Zhang     }
72678910aadSHong Zhang     uj_ptr[k]     = current_space->array;
72778910aadSHong Zhang     uj_lvl_ptr[k] = current_space_lvl->array;
72878910aadSHong Zhang 
72978910aadSHong Zhang     current_space->array           += nzk;
73078910aadSHong Zhang     current_space->local_used      += nzk;
73178910aadSHong Zhang     current_space->local_remaining -= nzk;
73278910aadSHong Zhang 
73378910aadSHong Zhang     current_space_lvl->array           += nzk;
73478910aadSHong Zhang     current_space_lvl->local_used      += nzk;
73578910aadSHong Zhang     current_space_lvl->local_remaining -= nzk;
73678910aadSHong Zhang 
73778910aadSHong Zhang     ui[k+1] = ui[k] + nzk;
73878910aadSHong Zhang   }
73978910aadSHong Zhang 
7406cf91177SBarry Smith #if defined(PETSC_USE_INFO)
74178910aadSHong Zhang   if (ai[am] != 0) {
74278910aadSHong Zhang     PetscReal af = ((PetscReal)(2*ui[am]-am))/((PetscReal)ai[am]);
743ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
744ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
745ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
74678910aadSHong Zhang   } else {
747ae15b995SBarry Smith     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
74878910aadSHong Zhang   }
74963ba0a88SBarry Smith #endif
75078910aadSHong Zhang 
75178910aadSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
75278910aadSHong Zhang   ierr = PetscFree(jl);CHKERRQ(ierr);
75378910aadSHong Zhang   ierr = PetscFree(cols_lvl);CHKERRQ(ierr);
75478910aadSHong Zhang 
75578910aadSHong Zhang   /* destroy list of free space and other temporary array(s) */
75678910aadSHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
757a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
75878910aadSHong Zhang   ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
759a1a86e44SBarry Smith   ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
76078910aadSHong Zhang 
76178910aadSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
762719d5645SBarry Smith   B = fact;
763ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation(B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
76478910aadSHong Zhang 
76578910aadSHong Zhang   b = (Mat_SeqSBAIJ*)B->data;
76678910aadSHong Zhang   b->singlemalloc = PETSC_FALSE;
767e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
768e6b907acSBarry Smith   b->free_ij       = PETSC_TRUE;
76978910aadSHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
77078910aadSHong Zhang   b->j    = uj;
77178910aadSHong Zhang   b->i    = ui;
77278910aadSHong Zhang   b->diag = 0;
77378910aadSHong Zhang   b->ilen = 0;
77478910aadSHong Zhang   b->imax = 0;
77578910aadSHong Zhang   b->row  = perm;
77678910aadSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
77778910aadSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
77878910aadSHong Zhang   b->icol = perm;
77978910aadSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
78078910aadSHong Zhang   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
78152e6d16bSBarry Smith   ierr    = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
78278910aadSHong Zhang   b->maxnz = b->nz = ui[am];
78378910aadSHong Zhang 
78478910aadSHong Zhang   B->info.factor_mallocs    = reallocs;
78578910aadSHong Zhang   B->info.fill_ratio_given  = fill;
78678910aadSHong Zhang   if (ai[am] != 0) {
78778910aadSHong Zhang     B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
78878910aadSHong Zhang   } else {
78978910aadSHong Zhang     B->info.fill_ratio_needed = 0.0;
79078910aadSHong Zhang   }
79178910aadSHong Zhang   if (perm_identity){
79278910aadSHong Zhang     B->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
79378910aadSHong Zhang     B->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
79478910aadSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
79578910aadSHong Zhang   } else {
796719d5645SBarry Smith     (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
79778910aadSHong Zhang   }
798c05c3958SHong Zhang   PetscFunctionReturn(0);
799c05c3958SHong Zhang }
800c05c3958SHong Zhang 
801c05c3958SHong Zhang #undef __FUNCT__
802c05c3958SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqBAIJ"
803719d5645SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS perm,MatFactorInfo *info)
804c05c3958SHong Zhang {
80578910aadSHong Zhang   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data;
80678910aadSHong Zhang   Mat_SeqSBAIJ       *b;
80778910aadSHong Zhang   Mat                B;
80878910aadSHong Zhang   PetscErrorCode     ierr;
80978910aadSHong Zhang   PetscTruth         perm_identity;
81078910aadSHong Zhang   PetscReal          fill = info->fill;
811*5d0c19d7SBarry Smith   const PetscInt     *rip;
812*5d0c19d7SBarry Smith   PetscInt           i,mbs=a->mbs,bs=A->rmap->bs,*ai=a->i,*aj=a->j,reallocs=0,prow;
81378910aadSHong Zhang   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
81478910aadSHong Zhang   PetscInt           nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr;
815a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
81678910aadSHong Zhang   PetscBT            lnkbt;
81778910aadSHong Zhang 
818c05c3958SHong Zhang   PetscFunctionBegin;
8196ad2eaddSHong Zhang   if (bs > 1) { /* convert to seqsbaij */
8206ad2eaddSHong Zhang     if (!a->sbaijMat){
821ceb03754SKris Buschelman       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
8226ad2eaddSHong Zhang     }
823719d5645SBarry Smith     (fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqSBAIJ; /* undue the change made in MatGetFactor_seqbaij_petsc */
824719d5645SBarry Smith     ierr = MatCholeskyFactorSymbolic(fact,a->sbaijMat,perm,info);CHKERRQ(ierr);
8256ad2eaddSHong Zhang     PetscFunctionReturn(0);
8266ad2eaddSHong Zhang   }
8276ad2eaddSHong Zhang 
82878910aadSHong Zhang   /* check whether perm is the identity mapping */
82978910aadSHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
830c84f5b01SHong Zhang   if (!perm_identity) SETERRQ(PETSC_ERR_SUP,"Matrix reordering is not supported");
83178910aadSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
83278910aadSHong Zhang 
83378910aadSHong Zhang   /* initialization */
83478910aadSHong Zhang   ierr  = PetscMalloc((mbs+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
83578910aadSHong Zhang   ui[0] = 0;
83678910aadSHong Zhang 
83778910aadSHong Zhang   /* jl: linked list for storing indices of the pivot rows
83878910aadSHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
83978910aadSHong Zhang   ierr = PetscMalloc((3*mbs+1)*sizeof(PetscInt)+mbs*sizeof(PetscInt*),&jl);CHKERRQ(ierr);
84078910aadSHong Zhang   il     = jl + mbs;
84178910aadSHong Zhang   cols   = il + mbs;
84278910aadSHong Zhang   ui_ptr = (PetscInt**)(cols + mbs);
84378910aadSHong Zhang   for (i=0; i<mbs; i++){
84478910aadSHong Zhang     jl[i] = mbs; il[i] = 0;
84578910aadSHong Zhang   }
84678910aadSHong Zhang 
84778910aadSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
84878910aadSHong Zhang   nlnk = mbs + 1;
84978910aadSHong Zhang   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
85078910aadSHong Zhang 
85178910aadSHong Zhang   /* initial FreeSpace size is fill*(ai[mbs]+1) */
852a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[mbs]+1)),&free_space);CHKERRQ(ierr);
85378910aadSHong Zhang   current_space = free_space;
85478910aadSHong Zhang 
85578910aadSHong Zhang   for (k=0; k<mbs; k++){  /* for each active row k */
85678910aadSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
85778910aadSHong Zhang     nzk   = 0;
85878910aadSHong Zhang     ncols = ai[rip[k]+1] - ai[rip[k]];
85978910aadSHong Zhang     ncols_upper = 0;
86078910aadSHong Zhang     for (j=0; j<ncols; j++){
86178910aadSHong Zhang       i = rip[*(aj + ai[rip[k]] + j)];
86278910aadSHong Zhang       if (i >= k){ /* only take upper triangular entry */
86378910aadSHong Zhang         cols[ncols_upper] = i;
86478910aadSHong Zhang         ncols_upper++;
86578910aadSHong Zhang       }
86678910aadSHong Zhang     }
86778910aadSHong Zhang     ierr = PetscLLAdd(ncols_upper,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
86878910aadSHong Zhang     nzk += nlnk;
86978910aadSHong Zhang 
87078910aadSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
87178910aadSHong Zhang     prow = jl[k]; /* 1st pivot row */
87278910aadSHong Zhang 
87378910aadSHong Zhang     while (prow < k){
87478910aadSHong Zhang       nextprow = jl[prow];
87578910aadSHong Zhang       /* merge prow into k-th row */
87678910aadSHong Zhang       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
87778910aadSHong Zhang       jmax = ui[prow+1];
87878910aadSHong Zhang       ncols = jmax-jmin;
87978910aadSHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
8805a8e39fbSHong Zhang       ierr = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
88178910aadSHong Zhang       nzk += nlnk;
88278910aadSHong Zhang 
88378910aadSHong Zhang       /* update il and jl for prow */
88478910aadSHong Zhang       if (jmin < jmax){
88578910aadSHong Zhang         il[prow] = jmin;
88678910aadSHong Zhang         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
88778910aadSHong Zhang       }
88878910aadSHong Zhang       prow = nextprow;
88978910aadSHong Zhang     }
89078910aadSHong Zhang 
89178910aadSHong Zhang     /* if free space is not available, make more free space */
89278910aadSHong Zhang     if (current_space->local_remaining<nzk) {
89378910aadSHong Zhang       i = mbs - k + 1; /* num of unfactored rows */
89478910aadSHong Zhang       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
895a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
89678910aadSHong Zhang       reallocs++;
89778910aadSHong Zhang     }
89878910aadSHong Zhang 
89978910aadSHong Zhang     /* copy data into free space, then initialize lnk */
90078910aadSHong Zhang     ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
90178910aadSHong Zhang 
90278910aadSHong Zhang     /* add the k-th row into il and jl */
90378910aadSHong Zhang     if (nzk-1 > 0){
90478910aadSHong Zhang       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
90578910aadSHong Zhang       jl[k] = jl[i]; jl[i] = k;
90678910aadSHong Zhang       il[k] = ui[k] + 1;
90778910aadSHong Zhang     }
90878910aadSHong Zhang     ui_ptr[k] = current_space->array;
90978910aadSHong Zhang     current_space->array           += nzk;
91078910aadSHong Zhang     current_space->local_used      += nzk;
91178910aadSHong Zhang     current_space->local_remaining -= nzk;
91278910aadSHong Zhang 
91378910aadSHong Zhang     ui[k+1] = ui[k] + nzk;
91478910aadSHong Zhang   }
91578910aadSHong Zhang 
9166cf91177SBarry Smith #if defined(PETSC_USE_INFO)
91778910aadSHong Zhang   if (ai[mbs] != 0) {
91878910aadSHong Zhang     PetscReal af = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
919ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
920ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
921ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
92278910aadSHong Zhang   } else {
923ae15b995SBarry Smith     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
92478910aadSHong Zhang   }
92563ba0a88SBarry Smith #endif
92678910aadSHong Zhang 
92778910aadSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
92878910aadSHong Zhang   ierr = PetscFree(jl);CHKERRQ(ierr);
92978910aadSHong Zhang 
93078910aadSHong Zhang   /* destroy list of free space and other temporary array(s) */
93178910aadSHong Zhang   ierr = PetscMalloc((ui[mbs]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
932a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
93378910aadSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
93478910aadSHong Zhang 
93578910aadSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
936719d5645SBarry Smith   B    = fact;
937ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
93878910aadSHong Zhang 
93978910aadSHong Zhang   b = (Mat_SeqSBAIJ*)B->data;
94078910aadSHong Zhang   b->singlemalloc = PETSC_FALSE;
941e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
942e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
94378910aadSHong Zhang   ierr = PetscMalloc((ui[mbs]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
94478910aadSHong Zhang   b->j    = uj;
94578910aadSHong Zhang   b->i    = ui;
94678910aadSHong Zhang   b->diag = 0;
94778910aadSHong Zhang   b->ilen = 0;
94878910aadSHong Zhang   b->imax = 0;
94978910aadSHong Zhang   b->row  = perm;
95078910aadSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
95178910aadSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
95278910aadSHong Zhang   b->icol = perm;
95378910aadSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
95478910aadSHong Zhang   ierr    = PetscMalloc((mbs+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
95552e6d16bSBarry Smith   ierr    = PetscLogObjectMemory(B,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
95678910aadSHong Zhang   b->maxnz = b->nz = ui[mbs];
95778910aadSHong Zhang 
95878910aadSHong Zhang   B->info.factor_mallocs    = reallocs;
95978910aadSHong Zhang   B->info.fill_ratio_given  = fill;
96078910aadSHong Zhang   if (ai[mbs] != 0) {
96178910aadSHong Zhang     B->info.fill_ratio_needed = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
96278910aadSHong Zhang   } else {
96378910aadSHong Zhang     B->info.fill_ratio_needed = 0.0;
96478910aadSHong Zhang   }
96578910aadSHong Zhang   if (perm_identity){
9666ad2eaddSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
96778910aadSHong Zhang   } else {
9686ad2eaddSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
96978910aadSHong Zhang   }
970c05c3958SHong Zhang   PetscFunctionReturn(0);
971c05c3958SHong Zhang }
972