xref: /petsc/src/mat/impls/baij/seq/baijfact.c (revision c8342467baa33ab03297e5da9c037bb665f9170c)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2be1d678aSKris Buschelman 
35d517e7eSBarry Smith /*
4ec3a8b7bSBarry Smith     Factorization code for BAIJ format.
55d517e7eSBarry Smith */
67c4f633dSBarry Smith #include "../src/mat/impls/baij/seq/baij.h"
77c4f633dSBarry 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"
150481f469SBarry Smith PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2(Mat B,Mat A,const 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;
215d0c19d7SBarry Smith   const PetscInt *r,*ic;
225d0c19d7SBarry 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         }
74dc0b31edSSatish Balay         ierr = PetscLogFlops(16.0*nz+12.0);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"
1060481f469SBarry Smith PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering(Mat C,Mat A,const 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         }
160dc0b31edSSatish Balay         ierr = PetscLogFlops(16.0*nz+12.0);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"
1920481f469SBarry Smith PetscErrorCode MatLUFactorNumeric_SeqBAIJ_1(Mat C,Mat A,const 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;
1975d0c19d7SBarry Smith   const PetscInt *r,*ic;
1985d0c19d7SBarry 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;
203f58c8c74SBarry Smith   PetscTruth     row_identity, col_identity;
2047fc0212eSBarry Smith 
2053a40ed3dSBarry Smith   PetscFunctionBegin;
2067fc0212eSBarry Smith   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
2077fc0212eSBarry Smith   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
208b0a32e0cSBarry Smith   ierr  = PetscMalloc((n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
2097fc0212eSBarry Smith 
2107fc0212eSBarry Smith   for (i=0; i<n; i++) {
2114078e994SBarry Smith     nz    = bi[i+1] - bi[i];
2124078e994SBarry Smith     ajtmp = bj + bi[i];
2137fc0212eSBarry Smith     for  (j=0; j<nz; j++) rtmp[ajtmp[j]] = 0.0;
2147fc0212eSBarry Smith 
2157fc0212eSBarry Smith     /* load in initial (unfactored row) */
2164078e994SBarry Smith     nz       = ai[r[i]+1] - ai[r[i]];
2174078e994SBarry Smith     ajtmpold = aj + ai[r[i]];
2184078e994SBarry Smith     v        = aa + ai[r[i]];
2197fc0212eSBarry Smith     for (j=0; j<nz; j++) rtmp[ic[ajtmpold[j]]] =  v[j];
2207fc0212eSBarry Smith 
2217fc0212eSBarry Smith     row = *ajtmp++;
2227fc0212eSBarry Smith     while (row < i) {
2237fc0212eSBarry Smith       pc = rtmp + row;
2247fc0212eSBarry Smith       if (*pc != 0.0) {
2254078e994SBarry Smith         pv         = ba + diag_offset[row];
2264078e994SBarry Smith         pj         = bj + diag_offset[row] + 1;
2277fc0212eSBarry Smith         multiplier = *pc * *pv++;
2287fc0212eSBarry Smith         *pc        = multiplier;
2294078e994SBarry Smith         nz         = bi[row+1] - diag_offset[row] - 1;
2307fc0212eSBarry Smith         for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j];
231dc0b31edSSatish Balay         ierr = PetscLogFlops(1.0+2.0*nz);CHKERRQ(ierr);
2327fc0212eSBarry Smith       }
2337fc0212eSBarry Smith       row = *ajtmp++;
2347fc0212eSBarry Smith     }
2357fc0212eSBarry Smith     /* finished row so stick it into b->a */
2364078e994SBarry Smith     pv = ba + bi[i];
2374078e994SBarry Smith     pj = bj + bi[i];
2384078e994SBarry Smith     nz = bi[i+1] - bi[i];
2397fc0212eSBarry Smith     for (j=0; j<nz; j++) {pv[j] = rtmp[pj[j]];}
2404078e994SBarry Smith     diag = diag_offset[i] - bi[i];
2417fc0212eSBarry Smith     /* check pivot entry for current row */
242a73cf429SBarry Smith     if (pv[diag] == 0.0) {
2433b4a8b6dSBarry Smith       SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot: row in original ordering %D in permuted ordering %D",r[i],i);
2447fc0212eSBarry Smith     }
2457fc0212eSBarry Smith     pv[diag] = 1.0/pv[diag];
2467fc0212eSBarry Smith   }
2477fc0212eSBarry Smith 
248606d414cSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
2497fc0212eSBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
2507fc0212eSBarry Smith   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
251f58c8c74SBarry Smith   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
252f58c8c74SBarry Smith   ierr = ISIdentity(isicol,&col_identity);CHKERRQ(ierr);
253f58c8c74SBarry Smith   if (row_identity && col_identity) {
254f58c8c74SBarry Smith     C->ops->solve          = MatSolve_SeqBAIJ_1_NaturalOrdering;
255f58c8c74SBarry Smith     C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1_NaturalOrdering;
256f58c8c74SBarry Smith   } else {
257db4efbfdSBarry Smith     C->ops->solve          = MatSolve_SeqBAIJ_1;
258db4efbfdSBarry Smith     C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1;
259f58c8c74SBarry Smith   }
2607fc0212eSBarry Smith   C->assembled = PETSC_TRUE;
261d0f46423SBarry Smith   ierr = PetscLogFlops(C->cmap->n);CHKERRQ(ierr);
2623a40ed3dSBarry Smith   PetscFunctionReturn(0);
2637fc0212eSBarry Smith }
2647fc0212eSBarry Smith 
265e631078cSBarry Smith EXTERN_C_BEGIN
266b24902e0SBarry Smith #undef __FUNCT__
267b24902e0SBarry Smith #define __FUNCT__ "MatGetFactor_seqbaij_petsc"
268b24902e0SBarry Smith PetscErrorCode MatGetFactor_seqbaij_petsc(Mat A,MatFactorType ftype,Mat *B)
269b24902e0SBarry Smith {
270d0f46423SBarry Smith   PetscInt           n = A->rmap->n;
271b24902e0SBarry Smith   PetscErrorCode     ierr;
272b24902e0SBarry Smith 
273b24902e0SBarry Smith   PetscFunctionBegin;
274b24902e0SBarry Smith   ierr = MatCreate(((PetscObject)A)->comm,B);CHKERRQ(ierr);
275b24902e0SBarry Smith   ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr);
276*c8342467SHong Zhang   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU || ftype == MAT_FACTOR_ILUDT) {
277b24902e0SBarry Smith     ierr = MatSetType(*B,MATSEQBAIJ);CHKERRQ(ierr);
278db4efbfdSBarry Smith     (*B)->ops->lufactorsymbolic  = MatLUFactorSymbolic_SeqBAIJ;
279db4efbfdSBarry Smith     (*B)->ops->ilufactorsymbolic = MatILUFactorSymbolic_SeqBAIJ;
280*c8342467SHong Zhang     (*B)->ops->iludtfactor       = MatILUDTFactor_SeqBAIJ;
281b24902e0SBarry Smith   } else if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) {
282b24902e0SBarry Smith     ierr = MatSetType(*B,MATSEQSBAIJ);CHKERRQ(ierr);
2835c9eb25fSBarry Smith     ierr = MatSeqSBAIJSetPreallocation(*B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
2845c9eb25fSBarry Smith     (*B)->ops->iccfactorsymbolic      = MatICCFactorSymbolic_SeqBAIJ;
2855c9eb25fSBarry Smith     (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqBAIJ;
286b24902e0SBarry Smith   } else SETERRQ(PETSC_ERR_SUP,"Factor type not supported");
287719d5645SBarry Smith   (*B)->factor = ftype;
288b24902e0SBarry Smith   PetscFunctionReturn(0);
289b24902e0SBarry Smith }
290e631078cSBarry Smith EXTERN_C_END
291273d9f13SBarry Smith 
292db4efbfdSBarry Smith EXTERN_C_BEGIN
293db4efbfdSBarry Smith #undef __FUNCT__
294db4efbfdSBarry Smith #define __FUNCT__ "MatGetFactorAvailable_seqaij_petsc"
295db4efbfdSBarry Smith PetscErrorCode MatGetFactorAvailable_seqbaij_petsc(Mat A,MatFactorType ftype,PetscTruth *flg)
296db4efbfdSBarry Smith {
297db4efbfdSBarry Smith   PetscFunctionBegin;
298db4efbfdSBarry Smith   *flg = PETSC_TRUE;
299db4efbfdSBarry Smith   PetscFunctionReturn(0);
300db4efbfdSBarry Smith }
301db4efbfdSBarry Smith EXTERN_C_END
302db4efbfdSBarry Smith 
3035d517e7eSBarry Smith /* ----------------------------------------------------------- */
3044a2ae208SSatish Balay #undef __FUNCT__
3054a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqBAIJ"
3060481f469SBarry Smith PetscErrorCode MatLUFactor_SeqBAIJ(Mat A,IS row,IS col,const MatFactorInfo *info)
3075d517e7eSBarry Smith {
308dfbe8321SBarry Smith   PetscErrorCode ierr;
3095d517e7eSBarry Smith   Mat            C;
3105d517e7eSBarry Smith 
3113a40ed3dSBarry Smith   PetscFunctionBegin;
31243244d56SBarry Smith   ierr = MatGetFactor(A,MAT_SOLVER_PETSC,MAT_FACTOR_LU,&C);CHKERRQ(ierr);
313719d5645SBarry Smith   ierr = MatLUFactorSymbolic(C,A,row,col,info);CHKERRQ(ierr);
314719d5645SBarry Smith   ierr = MatLUFactorNumeric(C,A,info);CHKERRQ(ierr);
315db4efbfdSBarry Smith   A->ops->solve            = C->ops->solve;
316db4efbfdSBarry Smith   A->ops->solvetranspose   = C->ops->solvetranspose;
317273d9f13SBarry Smith   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
31852e6d16bSBarry Smith   ierr = PetscLogObjectParent(A,((Mat_SeqBAIJ*)(A->data))->icol);CHKERRQ(ierr);
3193a40ed3dSBarry Smith   PetscFunctionReturn(0);
3205d517e7eSBarry Smith }
3214cd38560SBarry Smith 
3227c4f633dSBarry Smith #include "../src/mat/impls/sbaij/seq/sbaij.h"
323c05c3958SHong Zhang #undef __FUNCT__
324c05c3958SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N"
3250481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N(Mat C,Mat A,const MatFactorInfo *info)
326c05c3958SHong Zhang {
327f3086b4bSHong Zhang   PetscErrorCode ierr;
32878910aadSHong Zhang   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data;
32978910aadSHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
33078910aadSHong Zhang   IS             ip=b->row;
3315d0c19d7SBarry Smith   const PetscInt *rip;
3325d0c19d7SBarry Smith   PetscInt       i,j,mbs=a->mbs,bs=A->rmap->bs,*bi=b->i,*bj=b->j,*bcol;
33378910aadSHong Zhang   PetscInt       *ai=a->i,*aj=a->j;
33478910aadSHong Zhang   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
33578910aadSHong Zhang   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
3363cea4cbeSHong Zhang   PetscReal      zeropivot,rs,shiftnz;
337fbf22428SSatish Balay   PetscReal      shiftpd;
3383cea4cbeSHong Zhang   ChShift_Ctx    sctx;
3393cea4cbeSHong Zhang   PetscInt       newshift;
34078910aadSHong Zhang 
341c05c3958SHong Zhang   PetscFunctionBegin;
3426ad2eaddSHong Zhang   if (bs > 1) {
3436ad2eaddSHong Zhang     if (!a->sbaijMat){
344ceb03754SKris Buschelman       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
3456ad2eaddSHong Zhang     }
346719d5645SBarry Smith     ierr = (a->sbaijMat)->ops->choleskyfactornumeric(C,a->sbaijMat,info);CHKERRQ(ierr);
3476ad2eaddSHong Zhang     ierr = MatDestroy(a->sbaijMat);CHKERRQ(ierr);
3486ad2eaddSHong Zhang     a->sbaijMat = PETSC_NULL;
3496ad2eaddSHong Zhang     PetscFunctionReturn(0);
3506ad2eaddSHong Zhang   }
35178910aadSHong Zhang 
35278910aadSHong Zhang   /* initialization */
3533cea4cbeSHong Zhang   shiftnz   = info->shiftnz;
3543cea4cbeSHong Zhang   shiftpd   = info->shiftpd;
3553cea4cbeSHong Zhang   zeropivot = info->zeropivot;
3563cea4cbeSHong Zhang 
3576ad2eaddSHong Zhang   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
35878910aadSHong Zhang   nz   = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar);
35978910aadSHong Zhang   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
36078910aadSHong Zhang   jl   = il + mbs;
36178910aadSHong Zhang   rtmp = (MatScalar*)(jl + mbs);
36278910aadSHong Zhang 
3633cea4cbeSHong Zhang   sctx.shift_amount = 0;
3643cea4cbeSHong Zhang   sctx.nshift       = 0;
36578910aadSHong Zhang   do {
3663cea4cbeSHong Zhang     sctx.chshift = PETSC_FALSE;
36778910aadSHong Zhang     for (i=0; i<mbs; i++) {
36878910aadSHong Zhang       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
36978910aadSHong Zhang     }
37078910aadSHong Zhang 
37178910aadSHong Zhang     for (k = 0; k<mbs; k++){
37278910aadSHong Zhang       bval = ba + bi[k];
37378910aadSHong Zhang       /* initialize k-th row by the perm[k]-th row of A */
37478910aadSHong Zhang       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
37578910aadSHong Zhang       for (j = jmin; j < jmax; j++){
37678910aadSHong Zhang         col = rip[aj[j]];
37778910aadSHong Zhang         if (col >= k){ /* only take upper triangular entry */
37878910aadSHong Zhang           rtmp[col] = aa[j];
37978910aadSHong Zhang           *bval++  = 0.0; /* for in-place factorization */
38078910aadSHong Zhang         }
38178910aadSHong Zhang       }
3823cea4cbeSHong Zhang 
3833cea4cbeSHong Zhang       /* shift the diagonal of the matrix */
3843cea4cbeSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
38578910aadSHong Zhang 
38678910aadSHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
38778910aadSHong Zhang       dk = rtmp[k];
38878910aadSHong Zhang       i = jl[k]; /* first row to be added to k_th row  */
38978910aadSHong Zhang 
39078910aadSHong Zhang       while (i < k){
39178910aadSHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
39278910aadSHong Zhang 
39378910aadSHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
39478910aadSHong Zhang         ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
39578910aadSHong Zhang         uikdi = - ba[ili]*ba[bi[i]];  /* diagonal(k) */
39678910aadSHong Zhang         dk += uikdi*ba[ili];
39778910aadSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
39878910aadSHong Zhang 
39978910aadSHong Zhang         /* add multiple of row i to k-th row */
40078910aadSHong Zhang         jmin = ili + 1; jmax = bi[i+1];
40178910aadSHong Zhang         if (jmin < jmax){
40278910aadSHong Zhang           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
40378910aadSHong Zhang           /* update il and jl for row i */
40478910aadSHong Zhang           il[i] = jmin;
40578910aadSHong Zhang           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
40678910aadSHong Zhang         }
40778910aadSHong Zhang         i = nexti;
40878910aadSHong Zhang       }
40978910aadSHong Zhang 
4103cea4cbeSHong Zhang       /* shift the diagonals when zero pivot is detected */
4113cea4cbeSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
4123cea4cbeSHong Zhang       rs   = 0.0;
41378910aadSHong Zhang       jmin = bi[k]+1;
41478910aadSHong Zhang       nz   = bi[k+1] - jmin;
41578910aadSHong Zhang       if (nz){
41678910aadSHong Zhang         bcol = bj + jmin;
41778910aadSHong Zhang         while (nz--){
41889f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
41989f845c8SHong Zhang           bcol++;
42078910aadSHong Zhang         }
42178910aadSHong Zhang       }
4223cea4cbeSHong Zhang 
4233cea4cbeSHong Zhang       sctx.rs = rs;
4243cea4cbeSHong Zhang       sctx.pv = dk;
42545938b79SHong Zhang       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
42645938b79SHong Zhang       if (newshift == 1) break;
42778910aadSHong Zhang 
42878910aadSHong Zhang       /* copy data into U(k,:) */
42978910aadSHong Zhang       ba[bi[k]] = 1.0/dk; /* U(k,k) */
43078910aadSHong Zhang       jmin = bi[k]+1; jmax = bi[k+1];
43178910aadSHong Zhang       if (jmin < jmax) {
43278910aadSHong Zhang         for (j=jmin; j<jmax; j++){
43378910aadSHong Zhang           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
43478910aadSHong Zhang         }
43578910aadSHong Zhang         /* add the k-th row into il and jl */
43678910aadSHong Zhang         il[k] = jmin;
43778910aadSHong Zhang         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
43878910aadSHong Zhang       }
43978910aadSHong Zhang     }
4403cea4cbeSHong Zhang   } while (sctx.chshift);
44178910aadSHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
44278910aadSHong Zhang 
44378910aadSHong Zhang   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
44478910aadSHong Zhang   C->assembled    = PETSC_TRUE;
44578910aadSHong Zhang   C->preallocated = PETSC_TRUE;
446d0f46423SBarry Smith   ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr);
4473cea4cbeSHong Zhang   if (sctx.nshift){
448e738e422SBarry Smith     if (shiftpd) {
4491e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
450e738e422SBarry Smith     } else if (shiftnz) {
451e738e422SBarry Smith       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
45278910aadSHong Zhang     }
45378910aadSHong Zhang   }
454c05c3958SHong Zhang   PetscFunctionReturn(0);
455c05c3958SHong Zhang }
4564cd38560SBarry Smith 
457c05c3958SHong Zhang #undef __FUNCT__
458c05c3958SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering"
4590481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info)
460c05c3958SHong Zhang {
46178910aadSHong Zhang   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data;
46278910aadSHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
46378910aadSHong Zhang   PetscErrorCode ierr;
4643cea4cbeSHong Zhang   PetscInt       i,j,am=a->mbs;
46578910aadSHong Zhang   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
4663cea4cbeSHong Zhang   PetscInt       k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz;
46778910aadSHong Zhang   MatScalar      *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval;
4683cea4cbeSHong Zhang   PetscReal      zeropivot,rs,shiftnz;
469fbf22428SSatish Balay   PetscReal      shiftpd;
4703cea4cbeSHong Zhang   ChShift_Ctx    sctx;
4713cea4cbeSHong Zhang   PetscInt       newshift;
47278910aadSHong Zhang 
473c05c3958SHong Zhang   PetscFunctionBegin;
47478910aadSHong Zhang   /* initialization */
4753cea4cbeSHong Zhang   shiftnz   = info->shiftnz;
4763cea4cbeSHong Zhang   shiftpd   = info->shiftpd;
4773cea4cbeSHong Zhang   zeropivot = info->zeropivot;
4783cea4cbeSHong Zhang 
47978910aadSHong Zhang   nz   = (2*am+1)*sizeof(PetscInt)+am*sizeof(MatScalar);
48078910aadSHong Zhang   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
48178910aadSHong Zhang   jl   = il + am;
48278910aadSHong Zhang   rtmp = (MatScalar*)(jl + am);
48378910aadSHong Zhang 
4843cea4cbeSHong Zhang   sctx.shift_amount = 0;
4853cea4cbeSHong Zhang   sctx.nshift       = 0;
48678910aadSHong Zhang   do {
4873cea4cbeSHong Zhang     sctx.chshift = PETSC_FALSE;
48878910aadSHong Zhang     for (i=0; i<am; i++) {
48978910aadSHong Zhang       rtmp[i] = 0.0; jl[i] = am; il[0] = 0;
49078910aadSHong Zhang     }
49178910aadSHong Zhang 
49278910aadSHong Zhang     for (k = 0; k<am; k++){
49378910aadSHong Zhang     /* initialize k-th row with elements nonzero in row perm(k) of A */
49478910aadSHong Zhang       nz   = ai[k+1] - ai[k];
49578910aadSHong Zhang       acol = aj + ai[k];
49678910aadSHong Zhang       aval = aa + ai[k];
49778910aadSHong Zhang       bval = ba + bi[k];
49878910aadSHong Zhang       while (nz -- ){
49978910aadSHong Zhang         if (*acol < k) { /* skip lower triangular entries */
50078910aadSHong Zhang           acol++; aval++;
50178910aadSHong Zhang         } else {
50278910aadSHong Zhang           rtmp[*acol++] = *aval++;
50378910aadSHong Zhang           *bval++       = 0.0; /* for in-place factorization */
50478910aadSHong Zhang         }
50578910aadSHong Zhang       }
5063cea4cbeSHong Zhang 
5073cea4cbeSHong Zhang       /* shift the diagonal of the matrix */
5083cea4cbeSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
50978910aadSHong Zhang 
51078910aadSHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
51178910aadSHong Zhang       dk = rtmp[k];
51278910aadSHong Zhang       i  = jl[k]; /* first row to be added to k_th row  */
51378910aadSHong Zhang 
51478910aadSHong Zhang       while (i < k){
51578910aadSHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
51678910aadSHong Zhang         /* compute multiplier, update D(k) and U(i,k) */
51778910aadSHong Zhang         ili   = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
51878910aadSHong Zhang         uikdi = - ba[ili]*ba[bi[i]];
51978910aadSHong Zhang         dk   += uikdi*ba[ili];
52078910aadSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
52178910aadSHong Zhang 
52278910aadSHong Zhang         /* add multiple of row i to k-th row ... */
52378910aadSHong Zhang         jmin = ili + 1;
52478910aadSHong Zhang         nz   = bi[i+1] - jmin;
52578910aadSHong Zhang         if (nz > 0){
52678910aadSHong Zhang           bcol = bj + jmin;
52778910aadSHong Zhang           bval = ba + jmin;
52878910aadSHong Zhang           while (nz --) rtmp[*bcol++] += uikdi*(*bval++);
52978910aadSHong Zhang           /* update il and jl for i-th row */
53078910aadSHong Zhang           il[i] = jmin;
53178910aadSHong Zhang           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
53278910aadSHong Zhang         }
53378910aadSHong Zhang         i = nexti;
53478910aadSHong Zhang       }
53578910aadSHong Zhang 
5363cea4cbeSHong Zhang       /* shift the diagonals when zero pivot is detected */
5373cea4cbeSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
5383cea4cbeSHong Zhang       rs   = 0.0;
53978910aadSHong Zhang       jmin = bi[k]+1;
54078910aadSHong Zhang       nz   = bi[k+1] - jmin;
54178910aadSHong Zhang       if (nz){
54278910aadSHong Zhang         bcol = bj + jmin;
54378910aadSHong Zhang         while (nz--){
54489f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
54589f845c8SHong Zhang           bcol++;
54678910aadSHong Zhang         }
54778910aadSHong Zhang       }
5483cea4cbeSHong Zhang 
5493cea4cbeSHong Zhang       sctx.rs = rs;
5503cea4cbeSHong Zhang       sctx.pv = dk;
55145938b79SHong Zhang       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
55245938b79SHong Zhang       if (newshift == 1) break;    /* sctx.shift_amount is updated */
55378910aadSHong Zhang 
55478910aadSHong Zhang       /* copy data into U(k,:) */
55578910aadSHong Zhang       ba[bi[k]] = 1.0/dk;
55678910aadSHong Zhang       jmin      = bi[k]+1;
55778910aadSHong Zhang       nz        = bi[k+1] - jmin;
55878910aadSHong Zhang       if (nz){
55978910aadSHong Zhang         bcol = bj + jmin;
56078910aadSHong Zhang         bval = ba + jmin;
56178910aadSHong Zhang         while (nz--){
56278910aadSHong Zhang           *bval++       = rtmp[*bcol];
56378910aadSHong Zhang           rtmp[*bcol++] = 0.0;
56478910aadSHong Zhang         }
56578910aadSHong Zhang         /* add k-th row into il and jl */
56678910aadSHong Zhang         il[k] = jmin;
56778910aadSHong Zhang         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
56878910aadSHong Zhang       }
56978910aadSHong Zhang     }
5703cea4cbeSHong Zhang   } while (sctx.chshift);
57178910aadSHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
57278910aadSHong Zhang 
573719d5645SBarry Smith   C->ops->solve                 = MatSolve_SeqSBAIJ_1_NaturalOrdering;
574719d5645SBarry Smith   C->ops->solvetranspose        = MatSolve_SeqSBAIJ_1_NaturalOrdering;
57578910aadSHong Zhang   C->assembled    = PETSC_TRUE;
57678910aadSHong Zhang   C->preallocated = PETSC_TRUE;
577d0f46423SBarry Smith   ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr);
5783cea4cbeSHong Zhang     if (sctx.nshift){
5793cea4cbeSHong Zhang     if (shiftnz) {
5801e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
5813cea4cbeSHong Zhang     } else if (shiftpd) {
5821e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
58378910aadSHong Zhang     }
58478910aadSHong Zhang   }
585c05c3958SHong Zhang   PetscFunctionReturn(0);
586c05c3958SHong Zhang }
587c05c3958SHong Zhang 
588c05c3958SHong Zhang #include "petscbt.h"
5897c4f633dSBarry Smith #include "../src/mat/utils/freespace.h"
590c05c3958SHong Zhang #undef __FUNCT__
591c05c3958SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqBAIJ"
5920481f469SBarry Smith PetscErrorCode MatICCFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
593c05c3958SHong Zhang {
59478910aadSHong Zhang   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data;
59578910aadSHong Zhang   Mat_SeqSBAIJ       *b;
59678910aadSHong Zhang   Mat                B;
59778910aadSHong Zhang   PetscErrorCode     ierr;
59878910aadSHong Zhang   PetscTruth         perm_identity;
5995d0c19d7SBarry Smith   PetscInt           reallocs=0,i,*ai=a->i,*aj=a->j,am=a->mbs,bs=A->rmap->bs,*ui;
6005d0c19d7SBarry Smith   const PetscInt     *rip;
60178910aadSHong Zhang   PetscInt           jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow;
602cfdb8b8aSHong Zhang   PetscInt           nlnk,*lnk,*lnk_lvl=PETSC_NULL,ncols,ncols_upper,*cols,*cols_lvl,*uj,**uj_ptr,**uj_lvl_ptr;
60378910aadSHong Zhang   PetscReal          fill=info->fill,levels=info->levels;
604a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
605a1a86e44SBarry Smith   PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
60678910aadSHong Zhang   PetscBT            lnkbt;
60778910aadSHong Zhang 
608c05c3958SHong Zhang   PetscFunctionBegin;
6096ad2eaddSHong Zhang   if (bs > 1){
6106ad2eaddSHong Zhang     if (!a->sbaijMat){
611ceb03754SKris Buschelman       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
6126ad2eaddSHong Zhang     }
613719d5645SBarry Smith     (fact)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqSBAIJ;  /* undue the change made in MatGetFactor_seqbaij_petsc */
614719d5645SBarry Smith     ierr = MatICCFactorSymbolic(fact,a->sbaijMat,perm,info);CHKERRQ(ierr);
6156ad2eaddSHong Zhang     PetscFunctionReturn(0);
6166ad2eaddSHong Zhang   }
6176ad2eaddSHong Zhang 
61878910aadSHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
61978910aadSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
62078910aadSHong Zhang 
62178910aadSHong Zhang   /* special case that simply copies fill pattern */
62278910aadSHong Zhang   if (!levels && perm_identity) {
62378910aadSHong Zhang     ierr = MatMarkDiagonal_SeqBAIJ(A);CHKERRQ(ierr);
62478910aadSHong Zhang     ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
62578910aadSHong Zhang     for (i=0; i<am; i++) {
62678910aadSHong Zhang       ui[i] = ai[i+1] - a->diag[i]; /* ui: rowlengths - changes when !perm_identity */
62778910aadSHong Zhang     }
628719d5645SBarry Smith     B = fact;
62978910aadSHong Zhang     ierr = MatSeqSBAIJSetPreallocation(B,1,0,ui);CHKERRQ(ierr);
63078910aadSHong Zhang 
631b24902e0SBarry Smith 
63278910aadSHong Zhang     b  = (Mat_SeqSBAIJ*)B->data;
63378910aadSHong Zhang     uj = b->j;
63478910aadSHong Zhang     for (i=0; i<am; i++) {
63578910aadSHong Zhang       aj = a->j + a->diag[i];
63678910aadSHong Zhang       for (j=0; j<ui[i]; j++){
63778910aadSHong Zhang         *uj++ = *aj++;
63878910aadSHong Zhang       }
63978910aadSHong Zhang       b->ilen[i] = ui[i];
64078910aadSHong Zhang     }
64178910aadSHong Zhang     ierr = PetscFree(ui);CHKERRQ(ierr);
642719d5645SBarry Smith     B->factor = MAT_FACTOR_NONE;
64378910aadSHong Zhang     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
64478910aadSHong Zhang     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
645719d5645SBarry Smith     B->factor = MAT_FACTOR_ICC;
64678910aadSHong Zhang 
64778910aadSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
64878910aadSHong Zhang     PetscFunctionReturn(0);
64978910aadSHong Zhang   }
65078910aadSHong Zhang 
65178910aadSHong Zhang   /* initialization */
65278910aadSHong Zhang   ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
65378910aadSHong Zhang   ui[0] = 0;
65478910aadSHong Zhang   ierr  = PetscMalloc((2*am+1)*sizeof(PetscInt),&cols_lvl);CHKERRQ(ierr);
65578910aadSHong Zhang 
65678910aadSHong Zhang   /* jl: linked list for storing indices of the pivot rows
65778910aadSHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
65878910aadSHong Zhang   ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt*),&jl);CHKERRQ(ierr);
65978910aadSHong Zhang   il         = jl + am;
66078910aadSHong Zhang   uj_ptr     = (PetscInt**)(il + am);
66178910aadSHong Zhang   uj_lvl_ptr = (PetscInt**)(uj_ptr + am);
66278910aadSHong Zhang   for (i=0; i<am; i++){
66378910aadSHong Zhang     jl[i] = am; il[i] = 0;
66478910aadSHong Zhang   }
66578910aadSHong Zhang 
66678910aadSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
66778910aadSHong Zhang   nlnk = am + 1;
66878910aadSHong Zhang   ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
66978910aadSHong Zhang 
67078910aadSHong Zhang   /* initial FreeSpace size is fill*(ai[am]+1) */
671a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
67278910aadSHong Zhang   current_space = free_space;
673a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr);
67478910aadSHong Zhang   current_space_lvl = free_space_lvl;
67578910aadSHong Zhang 
67678910aadSHong Zhang   for (k=0; k<am; k++){  /* for each active row k */
67778910aadSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
67878910aadSHong Zhang     nzk   = 0;
67978910aadSHong Zhang     ncols = ai[rip[k]+1] - ai[rip[k]];
68078910aadSHong Zhang     ncols_upper = 0;
68178910aadSHong Zhang     cols        = cols_lvl + am;
68278910aadSHong Zhang     for (j=0; j<ncols; j++){
68378910aadSHong Zhang       i = rip[*(aj + ai[rip[k]] + j)];
68478910aadSHong Zhang       if (i >= k){ /* only take upper triangular entry */
68578910aadSHong Zhang         cols[ncols_upper] = i;
68678910aadSHong Zhang         cols_lvl[ncols_upper] = -1;  /* initialize level for nonzero entries */
68778910aadSHong Zhang         ncols_upper++;
68878910aadSHong Zhang       }
68978910aadSHong Zhang     }
69078910aadSHong Zhang     ierr = PetscIncompleteLLAdd(ncols_upper,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
69178910aadSHong Zhang     nzk += nlnk;
69278910aadSHong Zhang 
69378910aadSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
69478910aadSHong Zhang     prow = jl[k]; /* 1st pivot row */
69578910aadSHong Zhang 
69678910aadSHong Zhang     while (prow < k){
69778910aadSHong Zhang       nextprow = jl[prow];
69878910aadSHong Zhang 
69978910aadSHong Zhang       /* merge prow into k-th row */
70078910aadSHong Zhang       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
70178910aadSHong Zhang       jmax = ui[prow+1];
70278910aadSHong Zhang       ncols = jmax-jmin;
70378910aadSHong Zhang       i     = jmin - ui[prow];
70478910aadSHong Zhang       cols = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
70578910aadSHong Zhang       for (j=0; j<ncols; j++) cols_lvl[j] = *(uj_lvl_ptr[prow] + i + j);
7065a8e39fbSHong Zhang       ierr = PetscIncompleteLLAddSorted(ncols,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
70778910aadSHong Zhang       nzk += nlnk;
70878910aadSHong Zhang 
70978910aadSHong Zhang       /* update il and jl for prow */
71078910aadSHong Zhang       if (jmin < jmax){
71178910aadSHong Zhang         il[prow] = jmin;
71278910aadSHong Zhang         j = *cols; jl[prow] = jl[j]; jl[j] = prow;
71378910aadSHong Zhang       }
71478910aadSHong Zhang       prow = nextprow;
71578910aadSHong Zhang     }
71678910aadSHong Zhang 
71778910aadSHong Zhang     /* if free space is not available, make more free space */
71878910aadSHong Zhang     if (current_space->local_remaining<nzk) {
71978910aadSHong Zhang       i = am - k + 1; /* num of unfactored rows */
72078910aadSHong Zhang       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
721a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
722a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(i,&current_space_lvl);CHKERRQ(ierr);
72378910aadSHong Zhang       reallocs++;
72478910aadSHong Zhang     }
72578910aadSHong Zhang 
72678910aadSHong Zhang     /* copy data into free_space and free_space_lvl, then initialize lnk */
72778910aadSHong Zhang     ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
72878910aadSHong Zhang 
72978910aadSHong Zhang     /* add the k-th row into il and jl */
73078910aadSHong Zhang     if (nzk-1 > 0){
73178910aadSHong Zhang       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
73278910aadSHong Zhang       jl[k] = jl[i]; jl[i] = k;
73378910aadSHong Zhang       il[k] = ui[k] + 1;
73478910aadSHong Zhang     }
73578910aadSHong Zhang     uj_ptr[k]     = current_space->array;
73678910aadSHong Zhang     uj_lvl_ptr[k] = current_space_lvl->array;
73778910aadSHong Zhang 
73878910aadSHong Zhang     current_space->array           += nzk;
73978910aadSHong Zhang     current_space->local_used      += nzk;
74078910aadSHong Zhang     current_space->local_remaining -= nzk;
74178910aadSHong Zhang 
74278910aadSHong Zhang     current_space_lvl->array           += nzk;
74378910aadSHong Zhang     current_space_lvl->local_used      += nzk;
74478910aadSHong Zhang     current_space_lvl->local_remaining -= nzk;
74578910aadSHong Zhang 
74678910aadSHong Zhang     ui[k+1] = ui[k] + nzk;
74778910aadSHong Zhang   }
74878910aadSHong Zhang 
7496cf91177SBarry Smith #if defined(PETSC_USE_INFO)
75078910aadSHong Zhang   if (ai[am] != 0) {
75178910aadSHong Zhang     PetscReal af = ((PetscReal)(2*ui[am]-am))/((PetscReal)ai[am]);
752ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
753ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
754ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
75578910aadSHong Zhang   } else {
756ae15b995SBarry Smith     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
75778910aadSHong Zhang   }
75863ba0a88SBarry Smith #endif
75978910aadSHong Zhang 
76078910aadSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
76178910aadSHong Zhang   ierr = PetscFree(jl);CHKERRQ(ierr);
76278910aadSHong Zhang   ierr = PetscFree(cols_lvl);CHKERRQ(ierr);
76378910aadSHong Zhang 
76478910aadSHong Zhang   /* destroy list of free space and other temporary array(s) */
76578910aadSHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
766a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
76778910aadSHong Zhang   ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
768a1a86e44SBarry Smith   ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
76978910aadSHong Zhang 
77078910aadSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
771719d5645SBarry Smith   B = fact;
772ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation(B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
77378910aadSHong Zhang 
77478910aadSHong Zhang   b = (Mat_SeqSBAIJ*)B->data;
77578910aadSHong Zhang   b->singlemalloc = PETSC_FALSE;
776e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
777e6b907acSBarry Smith   b->free_ij       = PETSC_TRUE;
77878910aadSHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
77978910aadSHong Zhang   b->j    = uj;
78078910aadSHong Zhang   b->i    = ui;
78178910aadSHong Zhang   b->diag = 0;
78278910aadSHong Zhang   b->ilen = 0;
78378910aadSHong Zhang   b->imax = 0;
78478910aadSHong Zhang   b->row  = perm;
78578910aadSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
78678910aadSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
78778910aadSHong Zhang   b->icol = perm;
78878910aadSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
78978910aadSHong Zhang   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
79052e6d16bSBarry Smith   ierr    = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
79178910aadSHong Zhang   b->maxnz = b->nz = ui[am];
79278910aadSHong Zhang 
79378910aadSHong Zhang   B->info.factor_mallocs    = reallocs;
79478910aadSHong Zhang   B->info.fill_ratio_given  = fill;
79578910aadSHong Zhang   if (ai[am] != 0) {
79678910aadSHong Zhang     B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
79778910aadSHong Zhang   } else {
79878910aadSHong Zhang     B->info.fill_ratio_needed = 0.0;
79978910aadSHong Zhang   }
80078910aadSHong Zhang   if (perm_identity){
80178910aadSHong Zhang     B->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
80278910aadSHong Zhang     B->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
80378910aadSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
80478910aadSHong Zhang   } else {
805719d5645SBarry Smith     (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
80678910aadSHong Zhang   }
807c05c3958SHong Zhang   PetscFunctionReturn(0);
808c05c3958SHong Zhang }
809c05c3958SHong Zhang 
810c05c3958SHong Zhang #undef __FUNCT__
811c05c3958SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqBAIJ"
8120481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
813c05c3958SHong Zhang {
81478910aadSHong Zhang   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data;
81578910aadSHong Zhang   Mat_SeqSBAIJ       *b;
81678910aadSHong Zhang   Mat                B;
81778910aadSHong Zhang   PetscErrorCode     ierr;
81878910aadSHong Zhang   PetscTruth         perm_identity;
81978910aadSHong Zhang   PetscReal          fill = info->fill;
8205d0c19d7SBarry Smith   const PetscInt     *rip;
8215d0c19d7SBarry Smith   PetscInt           i,mbs=a->mbs,bs=A->rmap->bs,*ai=a->i,*aj=a->j,reallocs=0,prow;
82278910aadSHong Zhang   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
82378910aadSHong Zhang   PetscInt           nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr;
824a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
82578910aadSHong Zhang   PetscBT            lnkbt;
82678910aadSHong Zhang 
827c05c3958SHong Zhang   PetscFunctionBegin;
8286ad2eaddSHong Zhang   if (bs > 1) { /* convert to seqsbaij */
8296ad2eaddSHong Zhang     if (!a->sbaijMat){
830ceb03754SKris Buschelman       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
8316ad2eaddSHong Zhang     }
832719d5645SBarry Smith     (fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqSBAIJ; /* undue the change made in MatGetFactor_seqbaij_petsc */
833719d5645SBarry Smith     ierr = MatCholeskyFactorSymbolic(fact,a->sbaijMat,perm,info);CHKERRQ(ierr);
8346ad2eaddSHong Zhang     PetscFunctionReturn(0);
8356ad2eaddSHong Zhang   }
8366ad2eaddSHong Zhang 
83778910aadSHong Zhang   /* check whether perm is the identity mapping */
83878910aadSHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
839c84f5b01SHong Zhang   if (!perm_identity) SETERRQ(PETSC_ERR_SUP,"Matrix reordering is not supported");
84078910aadSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
84178910aadSHong Zhang 
84278910aadSHong Zhang   /* initialization */
84378910aadSHong Zhang   ierr  = PetscMalloc((mbs+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
84478910aadSHong Zhang   ui[0] = 0;
84578910aadSHong Zhang 
84678910aadSHong Zhang   /* jl: linked list for storing indices of the pivot rows
84778910aadSHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
84878910aadSHong Zhang   ierr = PetscMalloc((3*mbs+1)*sizeof(PetscInt)+mbs*sizeof(PetscInt*),&jl);CHKERRQ(ierr);
84978910aadSHong Zhang   il     = jl + mbs;
85078910aadSHong Zhang   cols   = il + mbs;
85178910aadSHong Zhang   ui_ptr = (PetscInt**)(cols + mbs);
85278910aadSHong Zhang   for (i=0; i<mbs; i++){
85378910aadSHong Zhang     jl[i] = mbs; il[i] = 0;
85478910aadSHong Zhang   }
85578910aadSHong Zhang 
85678910aadSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
85778910aadSHong Zhang   nlnk = mbs + 1;
85878910aadSHong Zhang   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
85978910aadSHong Zhang 
86078910aadSHong Zhang   /* initial FreeSpace size is fill*(ai[mbs]+1) */
861a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[mbs]+1)),&free_space);CHKERRQ(ierr);
86278910aadSHong Zhang   current_space = free_space;
86378910aadSHong Zhang 
86478910aadSHong Zhang   for (k=0; k<mbs; k++){  /* for each active row k */
86578910aadSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
86678910aadSHong Zhang     nzk   = 0;
86778910aadSHong Zhang     ncols = ai[rip[k]+1] - ai[rip[k]];
86878910aadSHong Zhang     ncols_upper = 0;
86978910aadSHong Zhang     for (j=0; j<ncols; j++){
87078910aadSHong Zhang       i = rip[*(aj + ai[rip[k]] + j)];
87178910aadSHong Zhang       if (i >= k){ /* only take upper triangular entry */
87278910aadSHong Zhang         cols[ncols_upper] = i;
87378910aadSHong Zhang         ncols_upper++;
87478910aadSHong Zhang       }
87578910aadSHong Zhang     }
87678910aadSHong Zhang     ierr = PetscLLAdd(ncols_upper,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
87778910aadSHong Zhang     nzk += nlnk;
87878910aadSHong Zhang 
87978910aadSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
88078910aadSHong Zhang     prow = jl[k]; /* 1st pivot row */
88178910aadSHong Zhang 
88278910aadSHong Zhang     while (prow < k){
88378910aadSHong Zhang       nextprow = jl[prow];
88478910aadSHong Zhang       /* merge prow into k-th row */
88578910aadSHong Zhang       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
88678910aadSHong Zhang       jmax = ui[prow+1];
88778910aadSHong Zhang       ncols = jmax-jmin;
88878910aadSHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
8895a8e39fbSHong Zhang       ierr = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
89078910aadSHong Zhang       nzk += nlnk;
89178910aadSHong Zhang 
89278910aadSHong Zhang       /* update il and jl for prow */
89378910aadSHong Zhang       if (jmin < jmax){
89478910aadSHong Zhang         il[prow] = jmin;
89578910aadSHong Zhang         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
89678910aadSHong Zhang       }
89778910aadSHong Zhang       prow = nextprow;
89878910aadSHong Zhang     }
89978910aadSHong Zhang 
90078910aadSHong Zhang     /* if free space is not available, make more free space */
90178910aadSHong Zhang     if (current_space->local_remaining<nzk) {
90278910aadSHong Zhang       i = mbs - k + 1; /* num of unfactored rows */
90378910aadSHong Zhang       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
904a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
90578910aadSHong Zhang       reallocs++;
90678910aadSHong Zhang     }
90778910aadSHong Zhang 
90878910aadSHong Zhang     /* copy data into free space, then initialize lnk */
90978910aadSHong Zhang     ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
91078910aadSHong Zhang 
91178910aadSHong Zhang     /* add the k-th row into il and jl */
91278910aadSHong Zhang     if (nzk-1 > 0){
91378910aadSHong Zhang       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
91478910aadSHong Zhang       jl[k] = jl[i]; jl[i] = k;
91578910aadSHong Zhang       il[k] = ui[k] + 1;
91678910aadSHong Zhang     }
91778910aadSHong Zhang     ui_ptr[k] = current_space->array;
91878910aadSHong Zhang     current_space->array           += nzk;
91978910aadSHong Zhang     current_space->local_used      += nzk;
92078910aadSHong Zhang     current_space->local_remaining -= nzk;
92178910aadSHong Zhang 
92278910aadSHong Zhang     ui[k+1] = ui[k] + nzk;
92378910aadSHong Zhang   }
92478910aadSHong Zhang 
9256cf91177SBarry Smith #if defined(PETSC_USE_INFO)
92678910aadSHong Zhang   if (ai[mbs] != 0) {
92778910aadSHong Zhang     PetscReal af = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
928ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
929ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
930ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
93178910aadSHong Zhang   } else {
932ae15b995SBarry Smith     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
93378910aadSHong Zhang   }
93463ba0a88SBarry Smith #endif
93578910aadSHong Zhang 
93678910aadSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
93778910aadSHong Zhang   ierr = PetscFree(jl);CHKERRQ(ierr);
93878910aadSHong Zhang 
93978910aadSHong Zhang   /* destroy list of free space and other temporary array(s) */
94078910aadSHong Zhang   ierr = PetscMalloc((ui[mbs]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
941a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
94278910aadSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
94378910aadSHong Zhang 
94478910aadSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
945719d5645SBarry Smith   B    = fact;
946ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
94778910aadSHong Zhang 
94878910aadSHong Zhang   b = (Mat_SeqSBAIJ*)B->data;
94978910aadSHong Zhang   b->singlemalloc = PETSC_FALSE;
950e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
951e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
95278910aadSHong Zhang   ierr = PetscMalloc((ui[mbs]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
95378910aadSHong Zhang   b->j    = uj;
95478910aadSHong Zhang   b->i    = ui;
95578910aadSHong Zhang   b->diag = 0;
95678910aadSHong Zhang   b->ilen = 0;
95778910aadSHong Zhang   b->imax = 0;
95878910aadSHong Zhang   b->row  = perm;
95978910aadSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
96078910aadSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
96178910aadSHong Zhang   b->icol = perm;
96278910aadSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
96378910aadSHong Zhang   ierr    = PetscMalloc((mbs+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
96452e6d16bSBarry Smith   ierr    = PetscLogObjectMemory(B,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
96578910aadSHong Zhang   b->maxnz = b->nz = ui[mbs];
96678910aadSHong Zhang 
96778910aadSHong Zhang   B->info.factor_mallocs    = reallocs;
96878910aadSHong Zhang   B->info.fill_ratio_given  = fill;
96978910aadSHong Zhang   if (ai[mbs] != 0) {
97078910aadSHong Zhang     B->info.fill_ratio_needed = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
97178910aadSHong Zhang   } else {
97278910aadSHong Zhang     B->info.fill_ratio_needed = 0.0;
97378910aadSHong Zhang   }
97478910aadSHong Zhang   if (perm_identity){
9756ad2eaddSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
97678910aadSHong Zhang   } else {
9776ad2eaddSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
97878910aadSHong Zhang   }
979c05c3958SHong Zhang   PetscFunctionReturn(0);
980c05c3958SHong Zhang }
981*c8342467SHong Zhang 
982*c8342467SHong Zhang #undef __FUNCT__
983*c8342467SHong Zhang #define __FUNCT__ "MatILUDTFactor_SeqBAIJ"
984*c8342467SHong Zhang PetscErrorCode MatILUDTFactor_SeqBAIJ(Mat A,IS isrow,IS iscol,const MatFactorInfo *info,Mat *fact)
985*c8342467SHong Zhang {
986*c8342467SHong Zhang     PetscFunctionBegin;
987*c8342467SHong Zhang     PetscFunctionReturn(0);
988*c8342467SHong Zhang }
989