xref: /petsc/src/mat/impls/baij/seq/baijfact.c (revision b588c5a2f37da9e74d7d91cdaec80d3d0decb4d9)
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"
7c60f0209SBarry Smith #include "../src/mat/blockinvert.h"
8ec3a8b7bSBarry Smith 
96d84be18SBarry Smith /* ------------------------------------------------------------*/
106d84be18SBarry Smith /*
114eeb42bcSBarry Smith       Version for when blocks are 2 by 2
124eeb42bcSBarry Smith */
13*b588c5a2SHong Zhang /* MatLUFactorNumeric_SeqBAIJ_2_newdatastruct -
14*b588c5a2SHong Zhang      copied from MatLUFactorNumeric_SeqBAIJ_N_newdatastruct() and re-implemented
15*b588c5a2SHong Zhang        Kernel_A_gets_A_times_B()
16*b588c5a2SHong Zhang        Kernel_A_gets_A_minus_B_times_C()
17*b588c5a2SHong Zhang        Kernel_A_gets_inverse_A()
18*b588c5a2SHong Zhang */
19*b588c5a2SHong Zhang #undef __FUNCT__
20*b588c5a2SHong Zhang #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2_newdatastruct"
21*b588c5a2SHong Zhang PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_newdatastruct(Mat B,Mat A,const MatFactorInfo *info)
22*b588c5a2SHong Zhang {
23*b588c5a2SHong Zhang   Mat            C=B;
24*b588c5a2SHong Zhang   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ *)C->data;
25*b588c5a2SHong Zhang   IS             isrow = b->row,isicol = b->icol;
26*b588c5a2SHong Zhang   PetscErrorCode ierr;
27*b588c5a2SHong Zhang   const PetscInt *r,*ic,*ics;
28*b588c5a2SHong Zhang   PetscInt       i,j,k,n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
29*b588c5a2SHong Zhang   PetscInt       *ajtmp,*bjtmp,nz,nzL,row,*bdiag=b->diag,*pj;
30*b588c5a2SHong Zhang   MatScalar      *rtmp,*pc,*mwork,*v,*pv,*aa=a->a;
31*b588c5a2SHong Zhang   PetscInt       bs=A->rmap->bs,bs2 = a->bs2,flg;
32*b588c5a2SHong Zhang   MatScalar      *x,m1,m2,m3,m4;
33*b588c5a2SHong Zhang   PetscReal      shift = info->shiftinblocks;
34*b588c5a2SHong Zhang 
35*b588c5a2SHong Zhang   PetscFunctionBegin;
36*b588c5a2SHong Zhang   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
37*b588c5a2SHong Zhang   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
38*b588c5a2SHong Zhang 
39*b588c5a2SHong Zhang   /* generate work space needed by LU factorization */
40*b588c5a2SHong Zhang   ierr = PetscMalloc((bs2*n+bs2+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
41*b588c5a2SHong Zhang   mwork = rtmp + bs2*n;
42*b588c5a2SHong Zhang   ierr = PetscMemzero(rtmp,(bs2*n+1)*sizeof(MatScalar));CHKERRQ(ierr);
43*b588c5a2SHong Zhang   ics  = ic;
44*b588c5a2SHong Zhang 
45*b588c5a2SHong Zhang   for (i=0; i<n; i++){
46*b588c5a2SHong Zhang     /* zero rtmp */
47*b588c5a2SHong Zhang     /* L part */
48*b588c5a2SHong Zhang     nz    = bi[i+1] - bi[i];
49*b588c5a2SHong Zhang     bjtmp = bj + bi[i];
50*b588c5a2SHong Zhang     for  (j=0; j<nz; j++){
51*b588c5a2SHong Zhang       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
52*b588c5a2SHong Zhang     }
53*b588c5a2SHong Zhang 
54*b588c5a2SHong Zhang     /* U part */
55*b588c5a2SHong Zhang     nz = bi[2*n-i+1] - bi[2*n-i];
56*b588c5a2SHong Zhang     bjtmp = bj + bi[2*n-i];
57*b588c5a2SHong Zhang     for  (j=0; j<nz; j++){
58*b588c5a2SHong Zhang       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
59*b588c5a2SHong Zhang     }
60*b588c5a2SHong Zhang 
61*b588c5a2SHong Zhang     /* load in initial (unfactored row) */
62*b588c5a2SHong Zhang     nz    = ai[r[i]+1] - ai[r[i]];
63*b588c5a2SHong Zhang     ajtmp = aj + ai[r[i]];
64*b588c5a2SHong Zhang     v     = aa + bs2*ai[r[i]];
65*b588c5a2SHong Zhang     for (j=0; j<nz; j++) {
66*b588c5a2SHong Zhang       ierr = PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr);
67*b588c5a2SHong Zhang     }
68*b588c5a2SHong Zhang 
69*b588c5a2SHong Zhang     /* elimination */
70*b588c5a2SHong Zhang     bjtmp = bj + bi[i];
71*b588c5a2SHong Zhang     row   = *bjtmp++;
72*b588c5a2SHong Zhang     nzL   = bi[i+1] - bi[i];
73*b588c5a2SHong Zhang     k   = 0;
74*b588c5a2SHong Zhang     while  (k < nzL) {
75*b588c5a2SHong Zhang       pc = rtmp + bs2*row;
76*b588c5a2SHong Zhang       for (flg=0,j=0; j<bs2; j++) { if (pc[j]!=0.0) { flg = 1; break; }}
77*b588c5a2SHong Zhang       if (flg) {
78*b588c5a2SHong Zhang         pv         = b->a + bs2*bdiag[row];
79*b588c5a2SHong Zhang         /* Kernel_A_gets_A_times_B(bs,pc,pv,mwork); */
80*b588c5a2SHong Zhang         /* ---Kernel_A_gets_A_times_B_2(): *pc = *pc * (*pv); ------- */
81*b588c5a2SHong Zhang         ierr = PetscMemcpy(mwork,pc,bs2*sizeof(MatScalar));CHKERRQ(ierr);
82*b588c5a2SHong Zhang         pc[0] = m1 = mwork[0]*pv[0] + mwork[2]*pv[1];
83*b588c5a2SHong Zhang         pc[1] = m2 = mwork[1]*pv[0] + mwork[3]*pv[1];
84*b588c5a2SHong Zhang         pc[2] = m3 = mwork[0]*pv[2] + mwork[2]*pv[3];
85*b588c5a2SHong Zhang         pc[3] = m4 = mwork[1]*pv[2] + mwork[3]*pv[3];
86*b588c5a2SHong Zhang         /* ------------ endof Kernel_A_gets_A_times_B_2() ----------- */
87*b588c5a2SHong Zhang 
88*b588c5a2SHong Zhang         pj         = b->j + bi[2*n-row]; /* begining of U(row,:) */
89*b588c5a2SHong Zhang         pv         = b->a + bs2*bi[2*n-row];
90*b588c5a2SHong Zhang         nz         = bi[2*n-row+1] - bi[2*n-row] - 1; /* num of entries inU(row,:), excluding diag */
91*b588c5a2SHong Zhang         for (j=0; j<nz; j++) {
92*b588c5a2SHong Zhang           /* Kernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
93*b588c5a2SHong Zhang           /* -- Kernel_A_gets_A_minus_B_times_C_2():  rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
94*b588c5a2SHong Zhang           x    = rtmp + 4*pj[j];
95*b588c5a2SHong Zhang           x[0] -= m1*pv[0] + m3*pv[1];
96*b588c5a2SHong Zhang           x[1] -= m2*pv[0] + m4*pv[1];
97*b588c5a2SHong Zhang           x[2] -= m1*pv[2] + m3*pv[3];
98*b588c5a2SHong Zhang           x[3] -= m2*pv[2] + m4*pv[3];
99*b588c5a2SHong Zhang           pv   += 4;
100*b588c5a2SHong Zhang           /*------------ endof Kernel_A_gets_A_minus_B_times_C_2() -- */
101*b588c5a2SHong Zhang         }
102*b588c5a2SHong Zhang         ierr = PetscLogFlops(16*nz+12);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
103*b588c5a2SHong Zhang       }
104*b588c5a2SHong Zhang       row = *bjtmp++; k++;
105*b588c5a2SHong Zhang     }
106*b588c5a2SHong Zhang 
107*b588c5a2SHong Zhang     /* finished row so stick it into b->a */
108*b588c5a2SHong Zhang     /* L part */
109*b588c5a2SHong Zhang     pv   = b->a + bs2*bi[i] ;
110*b588c5a2SHong Zhang     pj   = b->j + bi[i] ;
111*b588c5a2SHong Zhang     nz   = bi[i+1] - bi[i];
112*b588c5a2SHong Zhang     for (j=0; j<nz; j++) {
113*b588c5a2SHong Zhang       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
114*b588c5a2SHong Zhang     }
115*b588c5a2SHong Zhang 
116*b588c5a2SHong Zhang     /* Mark diagonal and invert diagonal for simplier triangular solves */
117*b588c5a2SHong Zhang     pv  = b->a + bs2*bdiag[i];
118*b588c5a2SHong Zhang     pj  = b->j + bdiag[i];
119*b588c5a2SHong Zhang     ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr);
120*b588c5a2SHong Zhang     // ierr = Kernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr);
121*b588c5a2SHong Zhang     ierr = Kernel_A_gets_inverse_A_2(pv,shift);CHKERRQ(ierr);
122*b588c5a2SHong Zhang 
123*b588c5a2SHong Zhang     /* U part */
124*b588c5a2SHong Zhang     pv = b->a + bs2*bi[2*n-i];
125*b588c5a2SHong Zhang     pj = b->j + bi[2*n-i];
126*b588c5a2SHong Zhang     nz = bi[2*n-i+1] - bi[2*n-i] - 1;
127*b588c5a2SHong Zhang     for (j=0; j<nz; j++){
128*b588c5a2SHong Zhang       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
129*b588c5a2SHong Zhang     }
130*b588c5a2SHong Zhang   }
131*b588c5a2SHong Zhang 
132*b588c5a2SHong Zhang   ierr = PetscFree(rtmp);CHKERRQ(ierr);
133*b588c5a2SHong Zhang   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
134*b588c5a2SHong Zhang   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
135*b588c5a2SHong Zhang 
136*b588c5a2SHong Zhang   C->assembled = PETSC_TRUE;
137*b588c5a2SHong Zhang   ierr = PetscLogFlops(1.3333*bs*bs2*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
138*b588c5a2SHong Zhang   PetscFunctionReturn(0);
139*b588c5a2SHong Zhang }
140*b588c5a2SHong Zhang 
1414a2ae208SSatish Balay #undef __FUNCT__
1424a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2"
1430481f469SBarry Smith PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2(Mat B,Mat A,const MatFactorInfo *info)
1444eeb42bcSBarry Smith {
145719d5645SBarry Smith   Mat            C = B;
1464eeb42bcSBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data;
1477cdf2dbbSSatish Balay   IS             isrow = b->row,isicol = b->icol;
1486849ba73SBarry Smith   PetscErrorCode ierr;
1495d0c19d7SBarry Smith   const PetscInt *r,*ic;
1505d0c19d7SBarry Smith   PetscInt       i,j,n = a->mbs,*bi = b->i,*bj = b->j;
151690b6cddSBarry Smith   PetscInt       *ajtmpold,*ajtmp,nz,row;
152690b6cddSBarry Smith   PetscInt       *diag_offset=b->diag,idx,*ai=a->i,*aj=a->j,*pj;
153329f5518SBarry Smith   MatScalar      *pv,*v,*rtmp,m1,m2,m3,m4,*pc,*w,*x,x1,x2,x3,x4;
1542515f8d2SSatish Balay   MatScalar      p1,p2,p3,p4;
1553f1db9ecSBarry Smith   MatScalar      *ba = b->a,*aa = a->a;
15662bba022SBarry Smith   PetscReal      shift = info->shiftinblocks;
1574eeb42bcSBarry Smith 
1583a40ed3dSBarry Smith   PetscFunctionBegin;
1594eeb42bcSBarry Smith   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
1604eeb42bcSBarry Smith   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
161b0a32e0cSBarry Smith   ierr  = PetscMalloc(4*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
1624eeb42bcSBarry Smith 
1634eeb42bcSBarry Smith   for (i=0; i<n; i++) {
1644078e994SBarry Smith     nz    = bi[i+1] - bi[i];
1654078e994SBarry Smith     ajtmp = bj + bi[i];
1664eeb42bcSBarry Smith     for  (j=0; j<nz; j++) {
1674eeb42bcSBarry Smith       x = rtmp+4*ajtmp[j]; x[0] = x[1] = x[2] = x[3] = 0.0;
1684eeb42bcSBarry Smith     }
1694eeb42bcSBarry Smith     /* load in initial (unfactored row) */
1704eeb42bcSBarry Smith     idx      = r[i];
1714078e994SBarry Smith     nz       = ai[idx+1] - ai[idx];
1724078e994SBarry Smith     ajtmpold = aj + ai[idx];
1734078e994SBarry Smith     v        = aa + 4*ai[idx];
1744eeb42bcSBarry Smith     for (j=0; j<nz; j++) {
1754eeb42bcSBarry Smith       x    = rtmp+4*ic[ajtmpold[j]];
1764eeb42bcSBarry Smith       x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3];
1774eeb42bcSBarry Smith       v    += 4;
1784eeb42bcSBarry Smith     }
1794eeb42bcSBarry Smith     row = *ajtmp++;
1804eeb42bcSBarry Smith     while (row < i) {
1814eeb42bcSBarry Smith       pc = rtmp + 4*row;
1824eeb42bcSBarry Smith       p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3];
18388685aaeSLois Curfman McInnes       if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0) {
1844078e994SBarry Smith         pv = ba + 4*diag_offset[row];
1854078e994SBarry Smith         pj = bj + diag_offset[row] + 1;
1864eeb42bcSBarry Smith         x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
1874eeb42bcSBarry Smith         pc[0] = m1 = p1*x1 + p3*x2;
1884eeb42bcSBarry Smith         pc[1] = m2 = p2*x1 + p4*x2;
1894eeb42bcSBarry Smith         pc[2] = m3 = p1*x3 + p3*x4;
1904eeb42bcSBarry Smith         pc[3] = m4 = p2*x3 + p4*x4;
1914078e994SBarry Smith         nz = bi[row+1] - diag_offset[row] - 1;
1924eeb42bcSBarry Smith         pv += 4;
1934eeb42bcSBarry Smith         for (j=0; j<nz; j++) {
1944eeb42bcSBarry Smith           x1   = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
1954eeb42bcSBarry Smith           x    = rtmp + 4*pj[j];
1964eeb42bcSBarry Smith           x[0] -= m1*x1 + m3*x2;
1974eeb42bcSBarry Smith           x[1] -= m2*x1 + m4*x2;
1984eeb42bcSBarry Smith           x[2] -= m1*x3 + m3*x4;
1994eeb42bcSBarry Smith           x[3] -= m2*x3 + m4*x4;
2004eeb42bcSBarry Smith           pv   += 4;
2014eeb42bcSBarry Smith         }
202dc0b31edSSatish Balay         ierr = PetscLogFlops(16.0*nz+12.0);CHKERRQ(ierr);
2034eeb42bcSBarry Smith       }
2044eeb42bcSBarry Smith       row = *ajtmp++;
2054eeb42bcSBarry Smith     }
2064eeb42bcSBarry Smith     /* finished row so stick it into b->a */
2074078e994SBarry Smith     pv = ba + 4*bi[i];
2084078e994SBarry Smith     pj = bj + bi[i];
2094078e994SBarry Smith     nz = bi[i+1] - bi[i];
2104eeb42bcSBarry Smith     for (j=0; j<nz; j++) {
2114eeb42bcSBarry Smith       x     = rtmp+4*pj[j];
2124eeb42bcSBarry Smith       pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3];
2134eeb42bcSBarry Smith       pv   += 4;
2144eeb42bcSBarry Smith     }
2154eeb42bcSBarry Smith     /* invert diagonal block */
2164078e994SBarry Smith     w = ba + 4*diag_offset[i];
21762bba022SBarry Smith     ierr = Kernel_A_gets_inverse_A_2(w,shift);CHKERRQ(ierr);
2184eeb42bcSBarry Smith   }
2194eeb42bcSBarry Smith 
220606d414cSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
2214eeb42bcSBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
2224eeb42bcSBarry Smith   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
223db4efbfdSBarry Smith   C->ops->solve          = MatSolve_SeqBAIJ_2;
224db4efbfdSBarry Smith   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2;
2254eeb42bcSBarry Smith   C->assembled = PETSC_TRUE;
226efee365bSSatish Balay   ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
2273a40ed3dSBarry Smith   PetscFunctionReturn(0);
2284eeb42bcSBarry Smith }
2294cd38560SBarry Smith /*
2304cd38560SBarry Smith       Version for when blocks are 2 by 2 Using natural ordering
2314cd38560SBarry Smith */
2324a2ae208SSatish Balay #undef __FUNCT__
2334a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering"
2340481f469SBarry Smith PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info)
2354cd38560SBarry Smith {
2364cd38560SBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data;
237dfbe8321SBarry Smith   PetscErrorCode ierr;
238690b6cddSBarry Smith   PetscInt       i,j,n = a->mbs,*bi = b->i,*bj = b->j;
239690b6cddSBarry Smith   PetscInt       *ajtmpold,*ajtmp,nz,row;
240690b6cddSBarry Smith   PetscInt       *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj;
241329f5518SBarry Smith   MatScalar      *pv,*v,*rtmp,*pc,*w,*x;
2422515f8d2SSatish Balay   MatScalar      p1,p2,p3,p4,m1,m2,m3,m4,x1,x2,x3,x4;
2434cd38560SBarry Smith   MatScalar      *ba = b->a,*aa = a->a;
24462bba022SBarry Smith   PetscReal      shift = info->shiftinblocks;
2454cd38560SBarry Smith 
2464cd38560SBarry Smith   PetscFunctionBegin;
247b0a32e0cSBarry Smith   ierr = PetscMalloc(4*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
2484cd38560SBarry Smith   for (i=0; i<n; i++) {
2494cd38560SBarry Smith     nz    = bi[i+1] - bi[i];
2504cd38560SBarry Smith     ajtmp = bj + bi[i];
2514cd38560SBarry Smith     for  (j=0; j<nz; j++) {
2524cd38560SBarry Smith       x = rtmp+4*ajtmp[j];
2534cd38560SBarry Smith       x[0]  = x[1]  = x[2]  = x[3]  = 0.0;
2544cd38560SBarry Smith     }
2554cd38560SBarry Smith     /* load in initial (unfactored row) */
2564cd38560SBarry Smith     nz       = ai[i+1] - ai[i];
2574cd38560SBarry Smith     ajtmpold = aj + ai[i];
2584cd38560SBarry Smith     v        = aa + 4*ai[i];
2594cd38560SBarry Smith     for (j=0; j<nz; j++) {
2604cd38560SBarry Smith       x    = rtmp+4*ajtmpold[j];
2614cd38560SBarry Smith       x[0]  = v[0];  x[1]  = v[1];  x[2]  = v[2];  x[3]  = v[3];
2624cd38560SBarry Smith       v    += 4;
2634cd38560SBarry Smith     }
2644cd38560SBarry Smith     row = *ajtmp++;
2654cd38560SBarry Smith     while (row < i) {
2664cd38560SBarry Smith       pc  = rtmp + 4*row;
2674cd38560SBarry Smith       p1  = pc[0];  p2  = pc[1];  p3  = pc[2];  p4  = pc[3];
2684cd38560SBarry Smith       if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0) {
2694cd38560SBarry Smith         pv = ba + 4*diag_offset[row];
2704cd38560SBarry Smith         pj = bj + diag_offset[row] + 1;
2714cd38560SBarry Smith         x1  = pv[0];  x2  = pv[1];  x3  = pv[2];  x4  = pv[3];
2724cd38560SBarry Smith         pc[0] = m1 = p1*x1 + p3*x2;
2734cd38560SBarry Smith         pc[1] = m2 = p2*x1 + p4*x2;
2744cd38560SBarry Smith         pc[2] = m3 = p1*x3 + p3*x4;
2754cd38560SBarry Smith         pc[3] = m4 = p2*x3 + p4*x4;
2764cd38560SBarry Smith         nz = bi[row+1] - diag_offset[row] - 1;
2774cd38560SBarry Smith         pv += 4;
2784cd38560SBarry Smith         for (j=0; j<nz; j++) {
2794cd38560SBarry Smith           x1   = pv[0];  x2  = pv[1];   x3 = pv[2];  x4  = pv[3];
2804cd38560SBarry Smith           x    = rtmp + 4*pj[j];
2814cd38560SBarry Smith           x[0] -= m1*x1 + m3*x2;
2824cd38560SBarry Smith           x[1] -= m2*x1 + m4*x2;
2834cd38560SBarry Smith           x[2] -= m1*x3 + m3*x4;
2844cd38560SBarry Smith           x[3] -= m2*x3 + m4*x4;
2854cd38560SBarry Smith           pv   += 4;
2864cd38560SBarry Smith         }
287dc0b31edSSatish Balay         ierr = PetscLogFlops(16.0*nz+12.0);CHKERRQ(ierr);
2884cd38560SBarry Smith       }
2894cd38560SBarry Smith       row = *ajtmp++;
2904cd38560SBarry Smith     }
2914cd38560SBarry Smith     /* finished row so stick it into b->a */
2924cd38560SBarry Smith     pv = ba + 4*bi[i];
2934cd38560SBarry Smith     pj = bj + bi[i];
2944cd38560SBarry Smith     nz = bi[i+1] - bi[i];
2954cd38560SBarry Smith     for (j=0; j<nz; j++) {
2964cd38560SBarry Smith       x      = rtmp+4*pj[j];
2974cd38560SBarry Smith       pv[0]  = x[0];  pv[1]  = x[1];  pv[2]  = x[2];  pv[3]  = x[3];
2982efa7f71SHong Zhang       /*
2992efa7f71SHong Zhang       printf(" col %d:",pj[j]);
3002efa7f71SHong Zhang       PetscInt j1;
3012efa7f71SHong Zhang       for (j1=0; j1<4; j1++) printf(" %g,",*(pv+j1));
3022efa7f71SHong Zhang       printf("\n");
3032efa7f71SHong Zhang       */
3044cd38560SBarry Smith       pv   += 4;
3054cd38560SBarry Smith     }
3064cd38560SBarry Smith     /* invert diagonal block */
3074cd38560SBarry Smith     w = ba + 4*diag_offset[i];
3082efa7f71SHong Zhang     /*
3092efa7f71SHong Zhang     printf(" \n%d -th: diag: ",i);
3102efa7f71SHong Zhang     for (j=0; j<4; j++){
3112efa7f71SHong Zhang       printf(" %g,",w[j]);
3122efa7f71SHong Zhang     }
3132efa7f71SHong Zhang     printf("\n----------------------------\n");
3142efa7f71SHong Zhang     */
31562bba022SBarry Smith     ierr = Kernel_A_gets_inverse_A_2(w,shift);CHKERRQ(ierr);
3164cd38560SBarry Smith   }
3174cd38560SBarry Smith 
318606d414cSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
319db4efbfdSBarry Smith   C->ops->solve          = MatSolve_SeqBAIJ_2_NaturalOrdering;
320db4efbfdSBarry Smith   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_NaturalOrdering;
3214cd38560SBarry Smith   C->assembled = PETSC_TRUE;
322efee365bSSatish Balay   ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
3234cd38560SBarry Smith   PetscFunctionReturn(0);
3244cd38560SBarry Smith }
3257fc0212eSBarry Smith 
3267fc0212eSBarry Smith /* ----------------------------------------------------------- */
3277fc0212eSBarry Smith /*
3287fc0212eSBarry Smith      Version for when blocks are 1 by 1.
3297fc0212eSBarry Smith */
3304a2ae208SSatish Balay #undef __FUNCT__
3314a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_1"
3320481f469SBarry Smith PetscErrorCode MatLUFactorNumeric_SeqBAIJ_1(Mat C,Mat A,const MatFactorInfo *info)
3337fc0212eSBarry Smith {
3347fc0212eSBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data;
3357cdf2dbbSSatish Balay   IS             isrow = b->row,isicol = b->icol;
3366849ba73SBarry Smith   PetscErrorCode ierr;
3375d0c19d7SBarry Smith   const PetscInt *r,*ic;
3385d0c19d7SBarry Smith   PetscInt       i,j,n = a->mbs,*bi = b->i,*bj = b->j;
339690b6cddSBarry Smith   PetscInt       *ajtmpold,*ajtmp,nz,row,*ai = a->i,*aj = a->j;
340690b6cddSBarry Smith   PetscInt       *diag_offset = b->diag,diag,*pj;
341329f5518SBarry Smith   MatScalar      *pv,*v,*rtmp,multiplier,*pc;
3423f1db9ecSBarry Smith   MatScalar      *ba = b->a,*aa = a->a;
343f58c8c74SBarry Smith   PetscTruth     row_identity, col_identity;
3447fc0212eSBarry Smith 
3453a40ed3dSBarry Smith   PetscFunctionBegin;
3467fc0212eSBarry Smith   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
3477fc0212eSBarry Smith   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
348b0a32e0cSBarry Smith   ierr  = PetscMalloc((n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
3497fc0212eSBarry Smith 
3507fc0212eSBarry Smith   for (i=0; i<n; i++) {
3514078e994SBarry Smith     nz    = bi[i+1] - bi[i];
3524078e994SBarry Smith     ajtmp = bj + bi[i];
3537fc0212eSBarry Smith     for  (j=0; j<nz; j++) rtmp[ajtmp[j]] = 0.0;
3547fc0212eSBarry Smith 
3557fc0212eSBarry Smith     /* load in initial (unfactored row) */
3564078e994SBarry Smith     nz       = ai[r[i]+1] - ai[r[i]];
3574078e994SBarry Smith     ajtmpold = aj + ai[r[i]];
3584078e994SBarry Smith     v        = aa + ai[r[i]];
3597fc0212eSBarry Smith     for (j=0; j<nz; j++) rtmp[ic[ajtmpold[j]]] =  v[j];
3607fc0212eSBarry Smith 
3617fc0212eSBarry Smith     row = *ajtmp++;
3627fc0212eSBarry Smith     while (row < i) {
3637fc0212eSBarry Smith       pc = rtmp + row;
3647fc0212eSBarry Smith       if (*pc != 0.0) {
3654078e994SBarry Smith         pv         = ba + diag_offset[row];
3664078e994SBarry Smith         pj         = bj + diag_offset[row] + 1;
3677fc0212eSBarry Smith         multiplier = *pc * *pv++;
3687fc0212eSBarry Smith         *pc        = multiplier;
3694078e994SBarry Smith         nz         = bi[row+1] - diag_offset[row] - 1;
3707fc0212eSBarry Smith         for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j];
371dc0b31edSSatish Balay         ierr = PetscLogFlops(1.0+2.0*nz);CHKERRQ(ierr);
3727fc0212eSBarry Smith       }
3737fc0212eSBarry Smith       row = *ajtmp++;
3747fc0212eSBarry Smith     }
3757fc0212eSBarry Smith     /* finished row so stick it into b->a */
3764078e994SBarry Smith     pv = ba + bi[i];
3774078e994SBarry Smith     pj = bj + bi[i];
3784078e994SBarry Smith     nz = bi[i+1] - bi[i];
3797fc0212eSBarry Smith     for (j=0; j<nz; j++) {pv[j] = rtmp[pj[j]];}
3804078e994SBarry Smith     diag = diag_offset[i] - bi[i];
3817fc0212eSBarry Smith     /* check pivot entry for current row */
382a73cf429SBarry Smith     if (pv[diag] == 0.0) {
3833b4a8b6dSBarry Smith       SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot: row in original ordering %D in permuted ordering %D",r[i],i);
3847fc0212eSBarry Smith     }
3857fc0212eSBarry Smith     pv[diag] = 1.0/pv[diag];
3867fc0212eSBarry Smith   }
3877fc0212eSBarry Smith 
388606d414cSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
3897fc0212eSBarry Smith   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
3907fc0212eSBarry Smith   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
391f58c8c74SBarry Smith   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
392f58c8c74SBarry Smith   ierr = ISIdentity(isicol,&col_identity);CHKERRQ(ierr);
393f58c8c74SBarry Smith   if (row_identity && col_identity) {
394f58c8c74SBarry Smith     C->ops->solve          = MatSolve_SeqBAIJ_1_NaturalOrdering;
395f58c8c74SBarry Smith     C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1_NaturalOrdering;
396f58c8c74SBarry Smith   } else {
397db4efbfdSBarry Smith     C->ops->solve          = MatSolve_SeqBAIJ_1;
398db4efbfdSBarry Smith     C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1;
399f58c8c74SBarry Smith   }
4007fc0212eSBarry Smith   C->assembled = PETSC_TRUE;
401d0f46423SBarry Smith   ierr = PetscLogFlops(C->cmap->n);CHKERRQ(ierr);
4023a40ed3dSBarry Smith   PetscFunctionReturn(0);
4037fc0212eSBarry Smith }
4047fc0212eSBarry Smith 
405e631078cSBarry Smith EXTERN_C_BEGIN
406b24902e0SBarry Smith #undef __FUNCT__
407b24902e0SBarry Smith #define __FUNCT__ "MatGetFactor_seqbaij_petsc"
408b24902e0SBarry Smith PetscErrorCode MatGetFactor_seqbaij_petsc(Mat A,MatFactorType ftype,Mat *B)
409b24902e0SBarry Smith {
410d0f46423SBarry Smith   PetscInt           n = A->rmap->n;
411b24902e0SBarry Smith   PetscErrorCode     ierr;
412b24902e0SBarry Smith 
413b24902e0SBarry Smith   PetscFunctionBegin;
414b24902e0SBarry Smith   ierr = MatCreate(((PetscObject)A)->comm,B);CHKERRQ(ierr);
415b24902e0SBarry Smith   ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr);
416c8342467SHong Zhang   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU || ftype == MAT_FACTOR_ILUDT) {
417b24902e0SBarry Smith     ierr = MatSetType(*B,MATSEQBAIJ);CHKERRQ(ierr);
418db4efbfdSBarry Smith     (*B)->ops->lufactorsymbolic  = MatLUFactorSymbolic_SeqBAIJ;
419db4efbfdSBarry Smith     (*B)->ops->ilufactorsymbolic = MatILUFactorSymbolic_SeqBAIJ;
420c8342467SHong Zhang     (*B)->ops->iludtfactor       = MatILUDTFactor_SeqBAIJ;
421b24902e0SBarry Smith   } else if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) {
422b24902e0SBarry Smith     ierr = MatSetType(*B,MATSEQSBAIJ);CHKERRQ(ierr);
4235c9eb25fSBarry Smith     ierr = MatSeqSBAIJSetPreallocation(*B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
4245c9eb25fSBarry Smith     (*B)->ops->iccfactorsymbolic      = MatICCFactorSymbolic_SeqBAIJ;
4255c9eb25fSBarry Smith     (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqBAIJ;
426b24902e0SBarry Smith   } else SETERRQ(PETSC_ERR_SUP,"Factor type not supported");
427719d5645SBarry Smith   (*B)->factor = ftype;
428b24902e0SBarry Smith   PetscFunctionReturn(0);
429b24902e0SBarry Smith }
430e631078cSBarry Smith EXTERN_C_END
431273d9f13SBarry Smith 
432db4efbfdSBarry Smith EXTERN_C_BEGIN
433db4efbfdSBarry Smith #undef __FUNCT__
434db4efbfdSBarry Smith #define __FUNCT__ "MatGetFactorAvailable_seqaij_petsc"
435db4efbfdSBarry Smith PetscErrorCode MatGetFactorAvailable_seqbaij_petsc(Mat A,MatFactorType ftype,PetscTruth *flg)
436db4efbfdSBarry Smith {
437db4efbfdSBarry Smith   PetscFunctionBegin;
438db4efbfdSBarry Smith   *flg = PETSC_TRUE;
439db4efbfdSBarry Smith   PetscFunctionReturn(0);
440db4efbfdSBarry Smith }
441db4efbfdSBarry Smith EXTERN_C_END
442db4efbfdSBarry Smith 
4435d517e7eSBarry Smith /* ----------------------------------------------------------- */
4444a2ae208SSatish Balay #undef __FUNCT__
4454a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqBAIJ"
4460481f469SBarry Smith PetscErrorCode MatLUFactor_SeqBAIJ(Mat A,IS row,IS col,const MatFactorInfo *info)
4475d517e7eSBarry Smith {
448dfbe8321SBarry Smith   PetscErrorCode ierr;
4495d517e7eSBarry Smith   Mat            C;
4505d517e7eSBarry Smith 
4513a40ed3dSBarry Smith   PetscFunctionBegin;
45243244d56SBarry Smith   ierr = MatGetFactor(A,MAT_SOLVER_PETSC,MAT_FACTOR_LU,&C);CHKERRQ(ierr);
453719d5645SBarry Smith   ierr = MatLUFactorSymbolic(C,A,row,col,info);CHKERRQ(ierr);
454719d5645SBarry Smith   ierr = MatLUFactorNumeric(C,A,info);CHKERRQ(ierr);
455db4efbfdSBarry Smith   A->ops->solve            = C->ops->solve;
456db4efbfdSBarry Smith   A->ops->solvetranspose   = C->ops->solvetranspose;
457273d9f13SBarry Smith   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
45852e6d16bSBarry Smith   ierr = PetscLogObjectParent(A,((Mat_SeqBAIJ*)(A->data))->icol);CHKERRQ(ierr);
4593a40ed3dSBarry Smith   PetscFunctionReturn(0);
4605d517e7eSBarry Smith }
4614cd38560SBarry Smith 
4627c4f633dSBarry Smith #include "../src/mat/impls/sbaij/seq/sbaij.h"
463c05c3958SHong Zhang #undef __FUNCT__
464c05c3958SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N"
4650481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N(Mat C,Mat A,const MatFactorInfo *info)
466c05c3958SHong Zhang {
467f3086b4bSHong Zhang   PetscErrorCode ierr;
46878910aadSHong Zhang   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data;
46978910aadSHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
47078910aadSHong Zhang   IS             ip=b->row;
4715d0c19d7SBarry Smith   const PetscInt *rip;
4725d0c19d7SBarry Smith   PetscInt       i,j,mbs=a->mbs,bs=A->rmap->bs,*bi=b->i,*bj=b->j,*bcol;
47378910aadSHong Zhang   PetscInt       *ai=a->i,*aj=a->j;
47478910aadSHong Zhang   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
47578910aadSHong Zhang   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
4763cea4cbeSHong Zhang   PetscReal      zeropivot,rs,shiftnz;
477fbf22428SSatish Balay   PetscReal      shiftpd;
4783cea4cbeSHong Zhang   ChShift_Ctx    sctx;
4793cea4cbeSHong Zhang   PetscInt       newshift;
48078910aadSHong Zhang 
481c05c3958SHong Zhang   PetscFunctionBegin;
4826ad2eaddSHong Zhang   if (bs > 1) {
4836ad2eaddSHong Zhang     if (!a->sbaijMat){
484ceb03754SKris Buschelman       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
4856ad2eaddSHong Zhang     }
486719d5645SBarry Smith     ierr = (a->sbaijMat)->ops->choleskyfactornumeric(C,a->sbaijMat,info);CHKERRQ(ierr);
4876ad2eaddSHong Zhang     ierr = MatDestroy(a->sbaijMat);CHKERRQ(ierr);
4886ad2eaddSHong Zhang     a->sbaijMat = PETSC_NULL;
4896ad2eaddSHong Zhang     PetscFunctionReturn(0);
4906ad2eaddSHong Zhang   }
49178910aadSHong Zhang 
49278910aadSHong Zhang   /* initialization */
4933cea4cbeSHong Zhang   shiftnz   = info->shiftnz;
4943cea4cbeSHong Zhang   shiftpd   = info->shiftpd;
4953cea4cbeSHong Zhang   zeropivot = info->zeropivot;
4963cea4cbeSHong Zhang 
4976ad2eaddSHong Zhang   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
49878910aadSHong Zhang   nz   = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar);
49978910aadSHong Zhang   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
50078910aadSHong Zhang   jl   = il + mbs;
50178910aadSHong Zhang   rtmp = (MatScalar*)(jl + mbs);
50278910aadSHong Zhang 
5033cea4cbeSHong Zhang   sctx.shift_amount = 0;
5043cea4cbeSHong Zhang   sctx.nshift       = 0;
50578910aadSHong Zhang   do {
5063cea4cbeSHong Zhang     sctx.chshift = PETSC_FALSE;
50778910aadSHong Zhang     for (i=0; i<mbs; i++) {
50878910aadSHong Zhang       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
50978910aadSHong Zhang     }
51078910aadSHong Zhang 
51178910aadSHong Zhang     for (k = 0; k<mbs; k++){
51278910aadSHong Zhang       bval = ba + bi[k];
51378910aadSHong Zhang       /* initialize k-th row by the perm[k]-th row of A */
51478910aadSHong Zhang       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
51578910aadSHong Zhang       for (j = jmin; j < jmax; j++){
51678910aadSHong Zhang         col = rip[aj[j]];
51778910aadSHong Zhang         if (col >= k){ /* only take upper triangular entry */
51878910aadSHong Zhang           rtmp[col] = aa[j];
51978910aadSHong Zhang           *bval++  = 0.0; /* for in-place factorization */
52078910aadSHong Zhang         }
52178910aadSHong Zhang       }
5223cea4cbeSHong Zhang 
5233cea4cbeSHong Zhang       /* shift the diagonal of the matrix */
5243cea4cbeSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
52578910aadSHong Zhang 
52678910aadSHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
52778910aadSHong Zhang       dk = rtmp[k];
52878910aadSHong Zhang       i = jl[k]; /* first row to be added to k_th row  */
52978910aadSHong Zhang 
53078910aadSHong Zhang       while (i < k){
53178910aadSHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
53278910aadSHong Zhang 
53378910aadSHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
53478910aadSHong Zhang         ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
53578910aadSHong Zhang         uikdi = - ba[ili]*ba[bi[i]];  /* diagonal(k) */
53678910aadSHong Zhang         dk += uikdi*ba[ili];
53778910aadSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
53878910aadSHong Zhang 
53978910aadSHong Zhang         /* add multiple of row i to k-th row */
54078910aadSHong Zhang         jmin = ili + 1; jmax = bi[i+1];
54178910aadSHong Zhang         if (jmin < jmax){
54278910aadSHong Zhang           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
54378910aadSHong Zhang           /* update il and jl for row i */
54478910aadSHong Zhang           il[i] = jmin;
54578910aadSHong Zhang           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
54678910aadSHong Zhang         }
54778910aadSHong Zhang         i = nexti;
54878910aadSHong Zhang       }
54978910aadSHong Zhang 
5503cea4cbeSHong Zhang       /* shift the diagonals when zero pivot is detected */
5513cea4cbeSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
5523cea4cbeSHong Zhang       rs   = 0.0;
55378910aadSHong Zhang       jmin = bi[k]+1;
55478910aadSHong Zhang       nz   = bi[k+1] - jmin;
55578910aadSHong Zhang       if (nz){
55678910aadSHong Zhang         bcol = bj + jmin;
55778910aadSHong Zhang         while (nz--){
55889f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
55989f845c8SHong Zhang           bcol++;
56078910aadSHong Zhang         }
56178910aadSHong Zhang       }
5623cea4cbeSHong Zhang 
5633cea4cbeSHong Zhang       sctx.rs = rs;
5643cea4cbeSHong Zhang       sctx.pv = dk;
56545938b79SHong Zhang       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
56645938b79SHong Zhang       if (newshift == 1) break;
56778910aadSHong Zhang 
56878910aadSHong Zhang       /* copy data into U(k,:) */
56978910aadSHong Zhang       ba[bi[k]] = 1.0/dk; /* U(k,k) */
57078910aadSHong Zhang       jmin = bi[k]+1; jmax = bi[k+1];
57178910aadSHong Zhang       if (jmin < jmax) {
57278910aadSHong Zhang         for (j=jmin; j<jmax; j++){
57378910aadSHong Zhang           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
57478910aadSHong Zhang         }
57578910aadSHong Zhang         /* add the k-th row into il and jl */
57678910aadSHong Zhang         il[k] = jmin;
57778910aadSHong Zhang         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
57878910aadSHong Zhang       }
57978910aadSHong Zhang     }
5803cea4cbeSHong Zhang   } while (sctx.chshift);
58178910aadSHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
58278910aadSHong Zhang 
58378910aadSHong Zhang   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
58478910aadSHong Zhang   C->assembled    = PETSC_TRUE;
58578910aadSHong Zhang   C->preallocated = PETSC_TRUE;
586d0f46423SBarry Smith   ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr);
5873cea4cbeSHong Zhang   if (sctx.nshift){
588e738e422SBarry Smith     if (shiftpd) {
5891e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
590e738e422SBarry Smith     } else if (shiftnz) {
591e738e422SBarry Smith       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
59278910aadSHong Zhang     }
59378910aadSHong Zhang   }
594c05c3958SHong Zhang   PetscFunctionReturn(0);
595c05c3958SHong Zhang }
5964cd38560SBarry Smith 
597c05c3958SHong Zhang #undef __FUNCT__
598c05c3958SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering"
5990481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info)
600c05c3958SHong Zhang {
60178910aadSHong Zhang   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data;
60278910aadSHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
60378910aadSHong Zhang   PetscErrorCode ierr;
6043cea4cbeSHong Zhang   PetscInt       i,j,am=a->mbs;
60578910aadSHong Zhang   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
6063cea4cbeSHong Zhang   PetscInt       k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz;
60778910aadSHong Zhang   MatScalar      *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval;
6083cea4cbeSHong Zhang   PetscReal      zeropivot,rs,shiftnz;
609fbf22428SSatish Balay   PetscReal      shiftpd;
6103cea4cbeSHong Zhang   ChShift_Ctx    sctx;
6113cea4cbeSHong Zhang   PetscInt       newshift;
61278910aadSHong Zhang 
613c05c3958SHong Zhang   PetscFunctionBegin;
61478910aadSHong Zhang   /* initialization */
6153cea4cbeSHong Zhang   shiftnz   = info->shiftnz;
6163cea4cbeSHong Zhang   shiftpd   = info->shiftpd;
6173cea4cbeSHong Zhang   zeropivot = info->zeropivot;
6183cea4cbeSHong Zhang 
61978910aadSHong Zhang   nz   = (2*am+1)*sizeof(PetscInt)+am*sizeof(MatScalar);
62078910aadSHong Zhang   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
62178910aadSHong Zhang   jl   = il + am;
62278910aadSHong Zhang   rtmp = (MatScalar*)(jl + am);
62378910aadSHong Zhang 
6243cea4cbeSHong Zhang   sctx.shift_amount = 0;
6253cea4cbeSHong Zhang   sctx.nshift       = 0;
62678910aadSHong Zhang   do {
6273cea4cbeSHong Zhang     sctx.chshift = PETSC_FALSE;
62878910aadSHong Zhang     for (i=0; i<am; i++) {
62978910aadSHong Zhang       rtmp[i] = 0.0; jl[i] = am; il[0] = 0;
63078910aadSHong Zhang     }
63178910aadSHong Zhang 
63278910aadSHong Zhang     for (k = 0; k<am; k++){
63378910aadSHong Zhang     /* initialize k-th row with elements nonzero in row perm(k) of A */
63478910aadSHong Zhang       nz   = ai[k+1] - ai[k];
63578910aadSHong Zhang       acol = aj + ai[k];
63678910aadSHong Zhang       aval = aa + ai[k];
63778910aadSHong Zhang       bval = ba + bi[k];
63878910aadSHong Zhang       while (nz -- ){
63978910aadSHong Zhang         if (*acol < k) { /* skip lower triangular entries */
64078910aadSHong Zhang           acol++; aval++;
64178910aadSHong Zhang         } else {
64278910aadSHong Zhang           rtmp[*acol++] = *aval++;
64378910aadSHong Zhang           *bval++       = 0.0; /* for in-place factorization */
64478910aadSHong Zhang         }
64578910aadSHong Zhang       }
6463cea4cbeSHong Zhang 
6473cea4cbeSHong Zhang       /* shift the diagonal of the matrix */
6483cea4cbeSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
64978910aadSHong Zhang 
65078910aadSHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
65178910aadSHong Zhang       dk = rtmp[k];
65278910aadSHong Zhang       i  = jl[k]; /* first row to be added to k_th row  */
65378910aadSHong Zhang 
65478910aadSHong Zhang       while (i < k){
65578910aadSHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
65678910aadSHong Zhang         /* compute multiplier, update D(k) and U(i,k) */
65778910aadSHong Zhang         ili   = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
65878910aadSHong Zhang         uikdi = - ba[ili]*ba[bi[i]];
65978910aadSHong Zhang         dk   += uikdi*ba[ili];
66078910aadSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
66178910aadSHong Zhang 
66278910aadSHong Zhang         /* add multiple of row i to k-th row ... */
66378910aadSHong Zhang         jmin = ili + 1;
66478910aadSHong Zhang         nz   = bi[i+1] - jmin;
66578910aadSHong Zhang         if (nz > 0){
66678910aadSHong Zhang           bcol = bj + jmin;
66778910aadSHong Zhang           bval = ba + jmin;
66878910aadSHong Zhang           while (nz --) rtmp[*bcol++] += uikdi*(*bval++);
66978910aadSHong Zhang           /* update il and jl for i-th row */
67078910aadSHong Zhang           il[i] = jmin;
67178910aadSHong Zhang           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
67278910aadSHong Zhang         }
67378910aadSHong Zhang         i = nexti;
67478910aadSHong Zhang       }
67578910aadSHong Zhang 
6763cea4cbeSHong Zhang       /* shift the diagonals when zero pivot is detected */
6773cea4cbeSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
6783cea4cbeSHong Zhang       rs   = 0.0;
67978910aadSHong Zhang       jmin = bi[k]+1;
68078910aadSHong Zhang       nz   = bi[k+1] - jmin;
68178910aadSHong Zhang       if (nz){
68278910aadSHong Zhang         bcol = bj + jmin;
68378910aadSHong Zhang         while (nz--){
68489f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
68589f845c8SHong Zhang           bcol++;
68678910aadSHong Zhang         }
68778910aadSHong Zhang       }
6883cea4cbeSHong Zhang 
6893cea4cbeSHong Zhang       sctx.rs = rs;
6903cea4cbeSHong Zhang       sctx.pv = dk;
69145938b79SHong Zhang       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
69245938b79SHong Zhang       if (newshift == 1) break;    /* sctx.shift_amount is updated */
69378910aadSHong Zhang 
69478910aadSHong Zhang       /* copy data into U(k,:) */
69578910aadSHong Zhang       ba[bi[k]] = 1.0/dk;
69678910aadSHong Zhang       jmin      = bi[k]+1;
69778910aadSHong Zhang       nz        = bi[k+1] - jmin;
69878910aadSHong Zhang       if (nz){
69978910aadSHong Zhang         bcol = bj + jmin;
70078910aadSHong Zhang         bval = ba + jmin;
70178910aadSHong Zhang         while (nz--){
70278910aadSHong Zhang           *bval++       = rtmp[*bcol];
70378910aadSHong Zhang           rtmp[*bcol++] = 0.0;
70478910aadSHong Zhang         }
70578910aadSHong Zhang         /* add k-th row into il and jl */
70678910aadSHong Zhang         il[k] = jmin;
70778910aadSHong Zhang         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
70878910aadSHong Zhang       }
70978910aadSHong Zhang     }
7103cea4cbeSHong Zhang   } while (sctx.chshift);
71178910aadSHong Zhang   ierr = PetscFree(il);CHKERRQ(ierr);
71278910aadSHong Zhang 
713719d5645SBarry Smith   C->ops->solve                 = MatSolve_SeqSBAIJ_1_NaturalOrdering;
714719d5645SBarry Smith   C->ops->solvetranspose        = MatSolve_SeqSBAIJ_1_NaturalOrdering;
71578910aadSHong Zhang   C->assembled    = PETSC_TRUE;
71678910aadSHong Zhang   C->preallocated = PETSC_TRUE;
717d0f46423SBarry Smith   ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr);
7183cea4cbeSHong Zhang     if (sctx.nshift){
7193cea4cbeSHong Zhang     if (shiftnz) {
7201e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
7213cea4cbeSHong Zhang     } else if (shiftpd) {
7221e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
72378910aadSHong Zhang     }
72478910aadSHong Zhang   }
725c05c3958SHong Zhang   PetscFunctionReturn(0);
726c05c3958SHong Zhang }
727c05c3958SHong Zhang 
728c05c3958SHong Zhang #include "petscbt.h"
7297c4f633dSBarry Smith #include "../src/mat/utils/freespace.h"
730c05c3958SHong Zhang #undef __FUNCT__
731c05c3958SHong Zhang #define __FUNCT__ "MatICCFactorSymbolic_SeqBAIJ"
7320481f469SBarry Smith PetscErrorCode MatICCFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
733c05c3958SHong Zhang {
73478910aadSHong Zhang   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data;
73578910aadSHong Zhang   Mat_SeqSBAIJ       *b;
73678910aadSHong Zhang   Mat                B;
73778910aadSHong Zhang   PetscErrorCode     ierr;
73878910aadSHong Zhang   PetscTruth         perm_identity;
7395d0c19d7SBarry Smith   PetscInt           reallocs=0,i,*ai=a->i,*aj=a->j,am=a->mbs,bs=A->rmap->bs,*ui;
7405d0c19d7SBarry Smith   const PetscInt     *rip;
74178910aadSHong Zhang   PetscInt           jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow;
742cfdb8b8aSHong Zhang   PetscInt           nlnk,*lnk,*lnk_lvl=PETSC_NULL,ncols,ncols_upper,*cols,*cols_lvl,*uj,**uj_ptr,**uj_lvl_ptr;
74378910aadSHong Zhang   PetscReal          fill=info->fill,levels=info->levels;
744a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
745a1a86e44SBarry Smith   PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
74678910aadSHong Zhang   PetscBT            lnkbt;
74778910aadSHong Zhang 
748c05c3958SHong Zhang   PetscFunctionBegin;
7496ad2eaddSHong Zhang   if (bs > 1){
7506ad2eaddSHong Zhang     if (!a->sbaijMat){
751ceb03754SKris Buschelman       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
7526ad2eaddSHong Zhang     }
753719d5645SBarry Smith     (fact)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqSBAIJ;  /* undue the change made in MatGetFactor_seqbaij_petsc */
754719d5645SBarry Smith     ierr = MatICCFactorSymbolic(fact,a->sbaijMat,perm,info);CHKERRQ(ierr);
7556ad2eaddSHong Zhang     PetscFunctionReturn(0);
7566ad2eaddSHong Zhang   }
7576ad2eaddSHong Zhang 
75878910aadSHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
75978910aadSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
76078910aadSHong Zhang 
76178910aadSHong Zhang   /* special case that simply copies fill pattern */
76278910aadSHong Zhang   if (!levels && perm_identity) {
76378910aadSHong Zhang     ierr = MatMarkDiagonal_SeqBAIJ(A);CHKERRQ(ierr);
76478910aadSHong Zhang     ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
76578910aadSHong Zhang     for (i=0; i<am; i++) {
76678910aadSHong Zhang       ui[i] = ai[i+1] - a->diag[i]; /* ui: rowlengths - changes when !perm_identity */
76778910aadSHong Zhang     }
768719d5645SBarry Smith     B = fact;
76978910aadSHong Zhang     ierr = MatSeqSBAIJSetPreallocation(B,1,0,ui);CHKERRQ(ierr);
77078910aadSHong Zhang 
771b24902e0SBarry Smith 
77278910aadSHong Zhang     b  = (Mat_SeqSBAIJ*)B->data;
77378910aadSHong Zhang     uj = b->j;
77478910aadSHong Zhang     for (i=0; i<am; i++) {
77578910aadSHong Zhang       aj = a->j + a->diag[i];
77678910aadSHong Zhang       for (j=0; j<ui[i]; j++){
77778910aadSHong Zhang         *uj++ = *aj++;
77878910aadSHong Zhang       }
77978910aadSHong Zhang       b->ilen[i] = ui[i];
78078910aadSHong Zhang     }
78178910aadSHong Zhang     ierr = PetscFree(ui);CHKERRQ(ierr);
782719d5645SBarry Smith     B->factor = MAT_FACTOR_NONE;
78378910aadSHong Zhang     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
78478910aadSHong Zhang     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
785719d5645SBarry Smith     B->factor = MAT_FACTOR_ICC;
78678910aadSHong Zhang 
78778910aadSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
78878910aadSHong Zhang     PetscFunctionReturn(0);
78978910aadSHong Zhang   }
79078910aadSHong Zhang 
79178910aadSHong Zhang   /* initialization */
79278910aadSHong Zhang   ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
79378910aadSHong Zhang   ui[0] = 0;
79478910aadSHong Zhang   ierr  = PetscMalloc((2*am+1)*sizeof(PetscInt),&cols_lvl);CHKERRQ(ierr);
79578910aadSHong Zhang 
79678910aadSHong Zhang   /* jl: linked list for storing indices of the pivot rows
79778910aadSHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
79878910aadSHong Zhang   ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt*),&jl);CHKERRQ(ierr);
79978910aadSHong Zhang   il         = jl + am;
80078910aadSHong Zhang   uj_ptr     = (PetscInt**)(il + am);
80178910aadSHong Zhang   uj_lvl_ptr = (PetscInt**)(uj_ptr + am);
80278910aadSHong Zhang   for (i=0; i<am; i++){
80378910aadSHong Zhang     jl[i] = am; il[i] = 0;
80478910aadSHong Zhang   }
80578910aadSHong Zhang 
80678910aadSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
80778910aadSHong Zhang   nlnk = am + 1;
80878910aadSHong Zhang   ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
80978910aadSHong Zhang 
81078910aadSHong Zhang   /* initial FreeSpace size is fill*(ai[am]+1) */
811a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
81278910aadSHong Zhang   current_space = free_space;
813a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr);
81478910aadSHong Zhang   current_space_lvl = free_space_lvl;
81578910aadSHong Zhang 
81678910aadSHong Zhang   for (k=0; k<am; k++){  /* for each active row k */
81778910aadSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
81878910aadSHong Zhang     nzk   = 0;
81978910aadSHong Zhang     ncols = ai[rip[k]+1] - ai[rip[k]];
82078910aadSHong Zhang     ncols_upper = 0;
82178910aadSHong Zhang     cols        = cols_lvl + am;
82278910aadSHong Zhang     for (j=0; j<ncols; j++){
82378910aadSHong Zhang       i = rip[*(aj + ai[rip[k]] + j)];
82478910aadSHong Zhang       if (i >= k){ /* only take upper triangular entry */
82578910aadSHong Zhang         cols[ncols_upper] = i;
82678910aadSHong Zhang         cols_lvl[ncols_upper] = -1;  /* initialize level for nonzero entries */
82778910aadSHong Zhang         ncols_upper++;
82878910aadSHong Zhang       }
82978910aadSHong Zhang     }
83078910aadSHong Zhang     ierr = PetscIncompleteLLAdd(ncols_upper,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
83178910aadSHong Zhang     nzk += nlnk;
83278910aadSHong Zhang 
83378910aadSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
83478910aadSHong Zhang     prow = jl[k]; /* 1st pivot row */
83578910aadSHong Zhang 
83678910aadSHong Zhang     while (prow < k){
83778910aadSHong Zhang       nextprow = jl[prow];
83878910aadSHong Zhang 
83978910aadSHong Zhang       /* merge prow into k-th row */
84078910aadSHong Zhang       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
84178910aadSHong Zhang       jmax = ui[prow+1];
84278910aadSHong Zhang       ncols = jmax-jmin;
84378910aadSHong Zhang       i     = jmin - ui[prow];
84478910aadSHong Zhang       cols = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
84578910aadSHong Zhang       for (j=0; j<ncols; j++) cols_lvl[j] = *(uj_lvl_ptr[prow] + i + j);
8465a8e39fbSHong Zhang       ierr = PetscIncompleteLLAddSorted(ncols,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
84778910aadSHong Zhang       nzk += nlnk;
84878910aadSHong Zhang 
84978910aadSHong Zhang       /* update il and jl for prow */
85078910aadSHong Zhang       if (jmin < jmax){
85178910aadSHong Zhang         il[prow] = jmin;
85278910aadSHong Zhang         j = *cols; jl[prow] = jl[j]; jl[j] = prow;
85378910aadSHong Zhang       }
85478910aadSHong Zhang       prow = nextprow;
85578910aadSHong Zhang     }
85678910aadSHong Zhang 
85778910aadSHong Zhang     /* if free space is not available, make more free space */
85878910aadSHong Zhang     if (current_space->local_remaining<nzk) {
85978910aadSHong Zhang       i = am - k + 1; /* num of unfactored rows */
86078910aadSHong Zhang       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
861a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
862a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(i,&current_space_lvl);CHKERRQ(ierr);
86378910aadSHong Zhang       reallocs++;
86478910aadSHong Zhang     }
86578910aadSHong Zhang 
86678910aadSHong Zhang     /* copy data into free_space and free_space_lvl, then initialize lnk */
86778910aadSHong Zhang     ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
86878910aadSHong Zhang 
86978910aadSHong Zhang     /* add the k-th row into il and jl */
87078910aadSHong Zhang     if (nzk-1 > 0){
87178910aadSHong Zhang       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
87278910aadSHong Zhang       jl[k] = jl[i]; jl[i] = k;
87378910aadSHong Zhang       il[k] = ui[k] + 1;
87478910aadSHong Zhang     }
87578910aadSHong Zhang     uj_ptr[k]     = current_space->array;
87678910aadSHong Zhang     uj_lvl_ptr[k] = current_space_lvl->array;
87778910aadSHong Zhang 
87878910aadSHong Zhang     current_space->array           += nzk;
87978910aadSHong Zhang     current_space->local_used      += nzk;
88078910aadSHong Zhang     current_space->local_remaining -= nzk;
88178910aadSHong Zhang 
88278910aadSHong Zhang     current_space_lvl->array           += nzk;
88378910aadSHong Zhang     current_space_lvl->local_used      += nzk;
88478910aadSHong Zhang     current_space_lvl->local_remaining -= nzk;
88578910aadSHong Zhang 
88678910aadSHong Zhang     ui[k+1] = ui[k] + nzk;
88778910aadSHong Zhang   }
88878910aadSHong Zhang 
8896cf91177SBarry Smith #if defined(PETSC_USE_INFO)
89078910aadSHong Zhang   if (ai[am] != 0) {
89178910aadSHong Zhang     PetscReal af = ((PetscReal)(2*ui[am]-am))/((PetscReal)ai[am]);
892ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
893ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
894ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
89578910aadSHong Zhang   } else {
896ae15b995SBarry Smith     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
89778910aadSHong Zhang   }
89863ba0a88SBarry Smith #endif
89978910aadSHong Zhang 
90078910aadSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
90178910aadSHong Zhang   ierr = PetscFree(jl);CHKERRQ(ierr);
90278910aadSHong Zhang   ierr = PetscFree(cols_lvl);CHKERRQ(ierr);
90378910aadSHong Zhang 
90478910aadSHong Zhang   /* destroy list of free space and other temporary array(s) */
90578910aadSHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
906a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
90778910aadSHong Zhang   ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
908a1a86e44SBarry Smith   ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
90978910aadSHong Zhang 
91078910aadSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
911719d5645SBarry Smith   B = fact;
912ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation(B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
91378910aadSHong Zhang 
91478910aadSHong Zhang   b = (Mat_SeqSBAIJ*)B->data;
91578910aadSHong Zhang   b->singlemalloc = PETSC_FALSE;
916e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
917e6b907acSBarry Smith   b->free_ij       = PETSC_TRUE;
91878910aadSHong Zhang   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
91978910aadSHong Zhang   b->j    = uj;
92078910aadSHong Zhang   b->i    = ui;
92178910aadSHong Zhang   b->diag = 0;
92278910aadSHong Zhang   b->ilen = 0;
92378910aadSHong Zhang   b->imax = 0;
92478910aadSHong Zhang   b->row  = perm;
92578910aadSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
92678910aadSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
92778910aadSHong Zhang   b->icol = perm;
92878910aadSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
92978910aadSHong Zhang   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
93052e6d16bSBarry Smith   ierr    = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
93178910aadSHong Zhang   b->maxnz = b->nz = ui[am];
93278910aadSHong Zhang 
93378910aadSHong Zhang   B->info.factor_mallocs    = reallocs;
93478910aadSHong Zhang   B->info.fill_ratio_given  = fill;
93578910aadSHong Zhang   if (ai[am] != 0) {
93678910aadSHong Zhang     B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
93778910aadSHong Zhang   } else {
93878910aadSHong Zhang     B->info.fill_ratio_needed = 0.0;
93978910aadSHong Zhang   }
94078910aadSHong Zhang   if (perm_identity){
94178910aadSHong Zhang     B->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
94278910aadSHong Zhang     B->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
94378910aadSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
94478910aadSHong Zhang   } else {
945719d5645SBarry Smith     (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
94678910aadSHong Zhang   }
947c05c3958SHong Zhang   PetscFunctionReturn(0);
948c05c3958SHong Zhang }
949c05c3958SHong Zhang 
950c05c3958SHong Zhang #undef __FUNCT__
951c05c3958SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqBAIJ"
9520481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
953c05c3958SHong Zhang {
95478910aadSHong Zhang   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data;
95578910aadSHong Zhang   Mat_SeqSBAIJ       *b;
95678910aadSHong Zhang   Mat                B;
95778910aadSHong Zhang   PetscErrorCode     ierr;
95878910aadSHong Zhang   PetscTruth         perm_identity;
95978910aadSHong Zhang   PetscReal          fill = info->fill;
9605d0c19d7SBarry Smith   const PetscInt     *rip;
9615d0c19d7SBarry Smith   PetscInt           i,mbs=a->mbs,bs=A->rmap->bs,*ai=a->i,*aj=a->j,reallocs=0,prow;
96278910aadSHong Zhang   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
96378910aadSHong Zhang   PetscInt           nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr;
964a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
96578910aadSHong Zhang   PetscBT            lnkbt;
96678910aadSHong Zhang 
967c05c3958SHong Zhang   PetscFunctionBegin;
9686ad2eaddSHong Zhang   if (bs > 1) { /* convert to seqsbaij */
9696ad2eaddSHong Zhang     if (!a->sbaijMat){
970ceb03754SKris Buschelman       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
9716ad2eaddSHong Zhang     }
972719d5645SBarry Smith     (fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqSBAIJ; /* undue the change made in MatGetFactor_seqbaij_petsc */
973719d5645SBarry Smith     ierr = MatCholeskyFactorSymbolic(fact,a->sbaijMat,perm,info);CHKERRQ(ierr);
9746ad2eaddSHong Zhang     PetscFunctionReturn(0);
9756ad2eaddSHong Zhang   }
9766ad2eaddSHong Zhang 
97778910aadSHong Zhang   /* check whether perm is the identity mapping */
97878910aadSHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
979c84f5b01SHong Zhang   if (!perm_identity) SETERRQ(PETSC_ERR_SUP,"Matrix reordering is not supported");
98078910aadSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
98178910aadSHong Zhang 
98278910aadSHong Zhang   /* initialization */
98378910aadSHong Zhang   ierr  = PetscMalloc((mbs+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
98478910aadSHong Zhang   ui[0] = 0;
98578910aadSHong Zhang 
98678910aadSHong Zhang   /* jl: linked list for storing indices of the pivot rows
98778910aadSHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
98878910aadSHong Zhang   ierr = PetscMalloc((3*mbs+1)*sizeof(PetscInt)+mbs*sizeof(PetscInt*),&jl);CHKERRQ(ierr);
98978910aadSHong Zhang   il     = jl + mbs;
99078910aadSHong Zhang   cols   = il + mbs;
99178910aadSHong Zhang   ui_ptr = (PetscInt**)(cols + mbs);
99278910aadSHong Zhang   for (i=0; i<mbs; i++){
99378910aadSHong Zhang     jl[i] = mbs; il[i] = 0;
99478910aadSHong Zhang   }
99578910aadSHong Zhang 
99678910aadSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
99778910aadSHong Zhang   nlnk = mbs + 1;
99878910aadSHong Zhang   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
99978910aadSHong Zhang 
100078910aadSHong Zhang   /* initial FreeSpace size is fill*(ai[mbs]+1) */
1001a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[mbs]+1)),&free_space);CHKERRQ(ierr);
100278910aadSHong Zhang   current_space = free_space;
100378910aadSHong Zhang 
100478910aadSHong Zhang   for (k=0; k<mbs; k++){  /* for each active row k */
100578910aadSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
100678910aadSHong Zhang     nzk   = 0;
100778910aadSHong Zhang     ncols = ai[rip[k]+1] - ai[rip[k]];
100878910aadSHong Zhang     ncols_upper = 0;
100978910aadSHong Zhang     for (j=0; j<ncols; j++){
101078910aadSHong Zhang       i = rip[*(aj + ai[rip[k]] + j)];
101178910aadSHong Zhang       if (i >= k){ /* only take upper triangular entry */
101278910aadSHong Zhang         cols[ncols_upper] = i;
101378910aadSHong Zhang         ncols_upper++;
101478910aadSHong Zhang       }
101578910aadSHong Zhang     }
101678910aadSHong Zhang     ierr = PetscLLAdd(ncols_upper,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
101778910aadSHong Zhang     nzk += nlnk;
101878910aadSHong Zhang 
101978910aadSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
102078910aadSHong Zhang     prow = jl[k]; /* 1st pivot row */
102178910aadSHong Zhang 
102278910aadSHong Zhang     while (prow < k){
102378910aadSHong Zhang       nextprow = jl[prow];
102478910aadSHong Zhang       /* merge prow into k-th row */
102578910aadSHong Zhang       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
102678910aadSHong Zhang       jmax = ui[prow+1];
102778910aadSHong Zhang       ncols = jmax-jmin;
102878910aadSHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
10295a8e39fbSHong Zhang       ierr = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
103078910aadSHong Zhang       nzk += nlnk;
103178910aadSHong Zhang 
103278910aadSHong Zhang       /* update il and jl for prow */
103378910aadSHong Zhang       if (jmin < jmax){
103478910aadSHong Zhang         il[prow] = jmin;
103578910aadSHong Zhang         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
103678910aadSHong Zhang       }
103778910aadSHong Zhang       prow = nextprow;
103878910aadSHong Zhang     }
103978910aadSHong Zhang 
104078910aadSHong Zhang     /* if free space is not available, make more free space */
104178910aadSHong Zhang     if (current_space->local_remaining<nzk) {
104278910aadSHong Zhang       i = mbs - k + 1; /* num of unfactored rows */
104378910aadSHong Zhang       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
1044a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
104578910aadSHong Zhang       reallocs++;
104678910aadSHong Zhang     }
104778910aadSHong Zhang 
104878910aadSHong Zhang     /* copy data into free space, then initialize lnk */
104978910aadSHong Zhang     ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
105078910aadSHong Zhang 
105178910aadSHong Zhang     /* add the k-th row into il and jl */
105278910aadSHong Zhang     if (nzk-1 > 0){
105378910aadSHong Zhang       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
105478910aadSHong Zhang       jl[k] = jl[i]; jl[i] = k;
105578910aadSHong Zhang       il[k] = ui[k] + 1;
105678910aadSHong Zhang     }
105778910aadSHong Zhang     ui_ptr[k] = current_space->array;
105878910aadSHong Zhang     current_space->array           += nzk;
105978910aadSHong Zhang     current_space->local_used      += nzk;
106078910aadSHong Zhang     current_space->local_remaining -= nzk;
106178910aadSHong Zhang 
106278910aadSHong Zhang     ui[k+1] = ui[k] + nzk;
106378910aadSHong Zhang   }
106478910aadSHong Zhang 
10656cf91177SBarry Smith #if defined(PETSC_USE_INFO)
106678910aadSHong Zhang   if (ai[mbs] != 0) {
106778910aadSHong Zhang     PetscReal af = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
1068ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
1069ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
1070ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
107178910aadSHong Zhang   } else {
1072ae15b995SBarry Smith     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
107378910aadSHong Zhang   }
107463ba0a88SBarry Smith #endif
107578910aadSHong Zhang 
107678910aadSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
107778910aadSHong Zhang   ierr = PetscFree(jl);CHKERRQ(ierr);
107878910aadSHong Zhang 
107978910aadSHong Zhang   /* destroy list of free space and other temporary array(s) */
108078910aadSHong Zhang   ierr = PetscMalloc((ui[mbs]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
1081a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
108278910aadSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
108378910aadSHong Zhang 
108478910aadSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
1085719d5645SBarry Smith   B    = fact;
1086ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
108778910aadSHong Zhang 
108878910aadSHong Zhang   b = (Mat_SeqSBAIJ*)B->data;
108978910aadSHong Zhang   b->singlemalloc = PETSC_FALSE;
1090e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
1091e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
109278910aadSHong Zhang   ierr = PetscMalloc((ui[mbs]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
109378910aadSHong Zhang   b->j    = uj;
109478910aadSHong Zhang   b->i    = ui;
109578910aadSHong Zhang   b->diag = 0;
109678910aadSHong Zhang   b->ilen = 0;
109778910aadSHong Zhang   b->imax = 0;
109878910aadSHong Zhang   b->row  = perm;
109978910aadSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
110078910aadSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
110178910aadSHong Zhang   b->icol = perm;
110278910aadSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
110378910aadSHong Zhang   ierr    = PetscMalloc((mbs+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
110452e6d16bSBarry Smith   ierr    = PetscLogObjectMemory(B,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
110578910aadSHong Zhang   b->maxnz = b->nz = ui[mbs];
110678910aadSHong Zhang 
110778910aadSHong Zhang   B->info.factor_mallocs    = reallocs;
110878910aadSHong Zhang   B->info.fill_ratio_given  = fill;
110978910aadSHong Zhang   if (ai[mbs] != 0) {
111078910aadSHong Zhang     B->info.fill_ratio_needed = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
111178910aadSHong Zhang   } else {
111278910aadSHong Zhang     B->info.fill_ratio_needed = 0.0;
111378910aadSHong Zhang   }
111478910aadSHong Zhang   if (perm_identity){
11156ad2eaddSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
111678910aadSHong Zhang   } else {
11176ad2eaddSHong Zhang     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
111878910aadSHong Zhang   }
1119c05c3958SHong Zhang   PetscFunctionReturn(0);
1120c05c3958SHong Zhang }
1121c8342467SHong Zhang 
11222efa7f71SHong Zhang /* --------------------------------------------------------- */
11232efa7f71SHong Zhang #undef __FUNCT__
112484a281e5SHong Zhang #define __FUNCT__ "MatSolve_SeqBAIJ_N_NaturalOrdering_newdatastruct"
112584a281e5SHong Zhang PetscErrorCode MatSolve_SeqBAIJ_N_NaturalOrdering_newdatastruct(Mat A,Vec bb,Vec xx)
11262efa7f71SHong Zhang {
11272efa7f71SHong Zhang   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ *)A->data;
11282efa7f71SHong Zhang   PetscErrorCode ierr;
11292efa7f71SHong Zhang   const PetscInt *ai=a->i,*aj=a->j,*vi;
11306464896eSShri Abhyankar   PetscInt       i,k,n=a->mbs;
11312efa7f71SHong Zhang   PetscInt       nz,bs=A->rmap->bs,bs2=a->bs2,*adiag=a->diag;
11322efa7f71SHong Zhang   MatScalar      *aa=a->a,*v;
11332efa7f71SHong Zhang   PetscScalar    *x,*b,*s,*t,*ls;
11342efa7f71SHong Zhang 
11352efa7f71SHong Zhang   PetscFunctionBegin;
11366bce7ff8SHong Zhang   /* printf("MatSolve_SeqBAIJ_NaturalOrdering_iludt..bs %d\n",bs); */
11372efa7f71SHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
11382efa7f71SHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
11392efa7f71SHong Zhang   t  = a->solve_work;
11402efa7f71SHong Zhang 
11412efa7f71SHong Zhang   /* forward solve the lower triangular */
11426da515a1SHong Zhang   ierr = PetscMemcpy(t,b,bs*sizeof(PetscScalar));CHKERRQ(ierr); /* copy 1st block of b to t */
11436bce7ff8SHong Zhang 
11442efa7f71SHong Zhang   for (i=1; i<n; i++) {
11452efa7f71SHong Zhang     v   = aa + bs2*ai[i];
11462efa7f71SHong Zhang     vi  = aj + ai[i];
11472efa7f71SHong Zhang     nz = ai[i+1] - ai[i];
11482efa7f71SHong Zhang     s = t + bs*i;
11496da515a1SHong Zhang     ierr = PetscMemcpy(s,b+bs*i,bs*sizeof(PetscScalar));CHKERRQ(ierr); /* copy i_th block of b to t */
11506464896eSShri Abhyankar     for(k=0;k<nz;k++){
11516464896eSShri Abhyankar       Kernel_v_gets_v_minus_A_times_w(bs,s,v,t+bs*vi[k]);
11522efa7f71SHong Zhang       v += bs2;
11532efa7f71SHong Zhang     }
11542efa7f71SHong Zhang   }
11552efa7f71SHong Zhang 
11562efa7f71SHong Zhang   /* backward solve the upper triangular */
11572efa7f71SHong Zhang   ls = a->solve_work + A->cmap->n;
11582efa7f71SHong Zhang   for (i=n-1; i>=0; i--){
11596bce7ff8SHong Zhang     v  = aa + bs2*ai[2*n-i];
11606bce7ff8SHong Zhang     vi = aj + ai[2*n-i];
11616bce7ff8SHong Zhang     nz = ai[2*n-i +1] - ai[2*n-i]-1;
11622efa7f71SHong Zhang     ierr = PetscMemcpy(ls,t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr);
11636464896eSShri Abhyankar     for(k=0;k<nz;k++){
11646464896eSShri Abhyankar       Kernel_v_gets_v_minus_A_times_w(bs,ls,v,t+bs*vi[k]);
11652efa7f71SHong Zhang       v += bs2;
11662efa7f71SHong Zhang     }
11672efa7f71SHong Zhang     Kernel_w_gets_A_times_v(bs,ls,aa+bs2*adiag[i],t+i*bs); /* *inv(diagonal[i]) */
11682efa7f71SHong Zhang     ierr = PetscMemcpy(x+i*bs,t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr);
11692efa7f71SHong Zhang   }
11706bce7ff8SHong Zhang 
11712efa7f71SHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
11722efa7f71SHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
11732efa7f71SHong Zhang   ierr = PetscLogFlops(2.0*(a->bs2)*(a->nz) - A->rmap->bs*A->cmap->n);CHKERRQ(ierr);
11742efa7f71SHong Zhang   PetscFunctionReturn(0);
11752efa7f71SHong Zhang }
11762efa7f71SHong Zhang 
11772efa7f71SHong Zhang #undef __FUNCT__
117884a281e5SHong Zhang #define __FUNCT__ "MatSolve_SeqBAIJ_N_newdatastruct"
117984a281e5SHong Zhang PetscErrorCode MatSolve_SeqBAIJ_N_newdatastruct(Mat A,Vec bb,Vec xx)
11802efa7f71SHong Zhang {
11812efa7f71SHong Zhang   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ *)A->data;
11822efa7f71SHong Zhang   IS             iscol=a->col,isrow=a->row;
11832efa7f71SHong Zhang   PetscErrorCode ierr;
11842efa7f71SHong Zhang   const PetscInt *r,*c,*rout,*cout,*ai=a->i,*aj=a->j,*vi;
118529b92fc1SShri Abhyankar   PetscInt       i,m,n=a->mbs;
11868f690400SShri Abhyankar   PetscInt       nz,bs=A->rmap->bs,bs2=a->bs2,k;
11872efa7f71SHong Zhang   MatScalar      *aa=a->a,*v;
11882efa7f71SHong Zhang   PetscScalar    *x,*b,*s,*t,*ls;
11892efa7f71SHong Zhang 
11902efa7f71SHong Zhang   PetscFunctionBegin;
11912efa7f71SHong Zhang   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
11922efa7f71SHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
11932efa7f71SHong Zhang   t  = a->solve_work;
11942efa7f71SHong Zhang 
11952efa7f71SHong Zhang   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
119629b92fc1SShri Abhyankar   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
11972efa7f71SHong Zhang 
11982efa7f71SHong Zhang   /* forward solve the lower triangular */
119929b92fc1SShri Abhyankar   ierr = PetscMemcpy(t,b+bs*r[0],bs*sizeof(PetscScalar));CHKERRQ(ierr);
12002efa7f71SHong Zhang   for (i=1; i<n; i++) {
12012efa7f71SHong Zhang     v   = aa + bs2*ai[i];
12022efa7f71SHong Zhang     vi  = aj + ai[i];
12032efa7f71SHong Zhang     nz = ai[i+1] - ai[i];
12042efa7f71SHong Zhang     s = t + bs*i;
120529b92fc1SShri Abhyankar     ierr = PetscMemcpy(s,b+bs*r[i],bs*sizeof(PetscScalar));CHKERRQ(ierr);
120629b92fc1SShri Abhyankar     for(m=0;m<nz;m++){
120729b92fc1SShri Abhyankar       Kernel_v_gets_v_minus_A_times_w(bs,s,v,t+bs*vi[m]);
12082efa7f71SHong Zhang       v += bs2;
12092efa7f71SHong Zhang     }
12102efa7f71SHong Zhang   }
12112efa7f71SHong Zhang 
12122efa7f71SHong Zhang   /* backward solve the upper triangular */
12132efa7f71SHong Zhang   ls = a->solve_work + A->cmap->n;
12142efa7f71SHong Zhang   for (i=n-1; i>=0; i--){
12158f690400SShri Abhyankar     k  = 2*n-i;
12168f690400SShri Abhyankar     v  = aa + bs2*ai[k];
12178f690400SShri Abhyankar     vi = aj + ai[k];
12188f690400SShri Abhyankar     nz = ai[k+1] - ai[k] - 1;
12192efa7f71SHong Zhang     ierr = PetscMemcpy(ls,t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr);
122029b92fc1SShri Abhyankar     for(m=0;m<nz;m++){
122129b92fc1SShri Abhyankar       Kernel_v_gets_v_minus_A_times_w(bs,ls,v,t+bs*vi[m]);
12222efa7f71SHong Zhang       v += bs2;
12232efa7f71SHong Zhang     }
12248f690400SShri Abhyankar     Kernel_w_gets_A_times_v(bs,ls,v,t+i*bs); /* *inv(diagonal[i]) */
122529b92fc1SShri Abhyankar     ierr = PetscMemcpy(x + bs*c[i],t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr);
12262efa7f71SHong Zhang   }
12272efa7f71SHong Zhang   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
12282efa7f71SHong Zhang   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
12292efa7f71SHong Zhang   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
12302efa7f71SHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
12312efa7f71SHong Zhang   ierr = PetscLogFlops(2.0*(a->bs2)*(a->nz) - A->rmap->bs*A->cmap->n);CHKERRQ(ierr);
12322efa7f71SHong Zhang   PetscFunctionReturn(0);
12332efa7f71SHong Zhang }
12342efa7f71SHong Zhang 
12352efa7f71SHong Zhang #undef __FUNCT__
12362efa7f71SHong Zhang #define __FUNCT__ "BlockAbs_privat"
123791d91c02SMatthew Knepley PetscErrorCode BlockAbs_private(PetscInt nbs,PetscInt bs2,PetscScalar *blockarray,PetscReal *absarray)
12382efa7f71SHong Zhang {
12392efa7f71SHong Zhang   PetscErrorCode     ierr;
12402efa7f71SHong Zhang   PetscInt           i,j;
12412efa7f71SHong Zhang   PetscFunctionBegin;
12422efa7f71SHong Zhang   ierr = PetscMemzero(absarray,(nbs+1)*sizeof(PetscReal));CHKERRQ(ierr);
12432efa7f71SHong Zhang   for (i=0; i<nbs; i++){
12442efa7f71SHong Zhang     for (j=0; j<bs2; j++){
12452efa7f71SHong Zhang       if (absarray[i] < PetscAbsScalar(blockarray[i*nbs+j])) absarray[i] = PetscAbsScalar(blockarray[i*nbs+j]);
12462efa7f71SHong Zhang     }
12472efa7f71SHong Zhang   }
12482efa7f71SHong Zhang   PetscFunctionReturn(0);
12492efa7f71SHong Zhang }
12502efa7f71SHong Zhang 
1251c8342467SHong Zhang #undef __FUNCT__
1252c8342467SHong Zhang #define __FUNCT__ "MatILUDTFactor_SeqBAIJ"
1253c8342467SHong Zhang PetscErrorCode MatILUDTFactor_SeqBAIJ(Mat A,IS isrow,IS iscol,const MatFactorInfo *info,Mat *fact)
1254c8342467SHong Zhang {
12552efa7f71SHong Zhang   Mat                B = *fact;
12562efa7f71SHong Zhang   Mat_SeqBAIJ        *a=(Mat_SeqBAIJ*)A->data,*b;
12572efa7f71SHong Zhang   IS                 isicol;
12582efa7f71SHong Zhang   PetscErrorCode     ierr;
12592efa7f71SHong Zhang   const PetscInt     *r,*ic;
12602efa7f71SHong Zhang   PetscInt           i,mbs=a->mbs,bs=A->rmap->bs,bs2=a->bs2,*ai=a->i,*aj=a->j,*ajtmp,*adiag;
12612efa7f71SHong Zhang   PetscInt           *bi,*bj,*bdiag;
12622efa7f71SHong Zhang 
12632efa7f71SHong Zhang   PetscInt           row,nzi,nzi_bl,nzi_bu,*im,dtcount,nzi_al,nzi_au;
12642efa7f71SHong Zhang   PetscInt           nlnk,*lnk;
12652efa7f71SHong Zhang   PetscBT            lnkbt;
12662efa7f71SHong Zhang   PetscTruth         row_identity,icol_identity,both_identity;
12672efa7f71SHong Zhang   MatScalar          *aatmp,*pv,*batmp,*ba,*rtmp,*pc,*multiplier,*vtmp;
12682efa7f71SHong Zhang   const PetscInt     *ics;
12692efa7f71SHong Zhang   PetscInt           j,nz,*pj,*bjtmp,k,ncut,*jtmp;
12702efa7f71SHong Zhang 
12716da515a1SHong Zhang   PetscReal          dt=info->dt; /* shift=info->shiftinblocks; */
12722efa7f71SHong Zhang   PetscInt           nnz_max;
12732efa7f71SHong Zhang   PetscTruth         missing;
1274021822bcSHong Zhang   PetscReal          *vtmp_abs;
127597ef94ebSSatish Balay   MatScalar          *v_work;
127697ef94ebSSatish Balay   PetscInt           *v_pivots;
12772efa7f71SHong Zhang 
1278c8342467SHong Zhang   PetscFunctionBegin;
12792efa7f71SHong Zhang   /* ------- symbolic factorization, can be reused ---------*/
12802efa7f71SHong Zhang   ierr = MatMissingDiagonal(A,&missing,&i);CHKERRQ(ierr);
12812efa7f71SHong Zhang   if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i);
12822efa7f71SHong Zhang   adiag=a->diag;
12832efa7f71SHong Zhang 
12842efa7f71SHong Zhang   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
12852efa7f71SHong Zhang 
12862efa7f71SHong Zhang   /* bdiag is location of diagonal in factor */
12872efa7f71SHong Zhang   ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
12882efa7f71SHong Zhang 
12892efa7f71SHong Zhang   /* allocate row pointers bi */
12906bce7ff8SHong Zhang   ierr = PetscMalloc((2*mbs+2)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
12912efa7f71SHong Zhang 
12922efa7f71SHong Zhang   /* allocate bj and ba; max num of nonzero entries is (ai[n]+2*n*dtcount+2) */
12932efa7f71SHong Zhang   dtcount = (PetscInt)info->dtcount;
12946bce7ff8SHong Zhang   if (dtcount > mbs-1) dtcount = mbs-1;
12956bce7ff8SHong Zhang   nnz_max  = ai[mbs]+2*mbs*dtcount +2;
12966da515a1SHong Zhang   /* printf("MatILUDTFactor_SeqBAIJ, bs %d, ai[mbs] %d, nnz_max  %d, dtcount %d\n",bs,ai[mbs],nnz_max,dtcount); */
12972efa7f71SHong Zhang   ierr = PetscMalloc(nnz_max*sizeof(PetscInt),&bj);CHKERRQ(ierr);
12986bce7ff8SHong Zhang   nnz_max = nnz_max*bs2;
12992efa7f71SHong Zhang   ierr = PetscMalloc(nnz_max*sizeof(MatScalar),&ba);CHKERRQ(ierr);
13002efa7f71SHong Zhang 
13012efa7f71SHong Zhang   /* put together the new matrix */
13022efa7f71SHong Zhang   ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
13032efa7f71SHong Zhang   ierr = PetscLogObjectParent(B,isicol);CHKERRQ(ierr);
13042efa7f71SHong Zhang   b    = (Mat_SeqBAIJ*)(B)->data;
13052efa7f71SHong Zhang   b->free_a       = PETSC_TRUE;
13062efa7f71SHong Zhang   b->free_ij      = PETSC_TRUE;
13072efa7f71SHong Zhang   b->singlemalloc = PETSC_FALSE;
13082efa7f71SHong Zhang   b->a          = ba;
13092efa7f71SHong Zhang   b->j          = bj;
13102efa7f71SHong Zhang   b->i          = bi;
13112efa7f71SHong Zhang   b->diag       = bdiag;
13122efa7f71SHong Zhang   b->ilen       = 0;
13132efa7f71SHong Zhang   b->imax       = 0;
13142efa7f71SHong Zhang   b->row        = isrow;
13152efa7f71SHong Zhang   b->col        = iscol;
13162efa7f71SHong Zhang   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
13172efa7f71SHong Zhang   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
13182efa7f71SHong Zhang   b->icol       = isicol;
13192efa7f71SHong Zhang   ierr = PetscMalloc((bs*(mbs+1))*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
13202efa7f71SHong Zhang 
13212efa7f71SHong Zhang   ierr = PetscLogObjectMemory(B,nnz_max*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
13222efa7f71SHong Zhang   b->maxnz = nnz_max/bs2;
13232efa7f71SHong Zhang 
13242efa7f71SHong Zhang   (B)->factor                = MAT_FACTOR_ILUDT;
13252efa7f71SHong Zhang   (B)->info.factor_mallocs   = 0;
13262efa7f71SHong Zhang   (B)->info.fill_ratio_given = ((PetscReal)nnz_max)/((PetscReal)(ai[mbs]*bs2));
13272efa7f71SHong Zhang   CHKMEMQ;
13282efa7f71SHong Zhang   /* ------- end of symbolic factorization ---------*/
13292efa7f71SHong Zhang   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
13302efa7f71SHong Zhang   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
13312efa7f71SHong Zhang   ics  = ic;
13322efa7f71SHong Zhang 
13332efa7f71SHong Zhang   /* linked list for storing column indices of the active row */
13342efa7f71SHong Zhang   nlnk = mbs + 1;
13352efa7f71SHong Zhang   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
13362efa7f71SHong Zhang 
13372efa7f71SHong Zhang   /* im: used by PetscLLAddSortedLU(); jtmp: working array for column indices of active row */
13382efa7f71SHong Zhang   ierr = PetscMalloc((2*mbs+1)*sizeof(PetscInt),&im);CHKERRQ(ierr);
13392efa7f71SHong Zhang   jtmp = im + mbs;
13402efa7f71SHong Zhang   /* rtmp, vtmp: working arrays for sparse and contiguous row entries of active row */
13412efa7f71SHong Zhang   ierr = PetscMalloc((2*mbs*bs2+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
13422efa7f71SHong Zhang   vtmp = rtmp + bs2*mbs;
1343021822bcSHong Zhang   ierr = PetscMalloc((mbs+1)*sizeof(PetscReal),&vtmp_abs);CHKERRQ(ierr);
13442efa7f71SHong Zhang 
13452efa7f71SHong Zhang   ierr       = PetscMalloc(bs*sizeof(PetscInt) + (bs+bs2)*sizeof(MatScalar),&v_work);CHKERRQ(ierr);
13462efa7f71SHong Zhang   multiplier = v_work + bs;
13472efa7f71SHong Zhang   v_pivots   = (PetscInt*)(multiplier + bs2);
13482efa7f71SHong Zhang 
13492efa7f71SHong Zhang   bi[0]    = 0;
13502efa7f71SHong Zhang   bdiag[0] = (nnz_max/bs2)-1; /* location of diagonal in factor B */
13516bce7ff8SHong Zhang   bi[2*mbs+1] = bdiag[0]+1; /* endof bj and ba array */
13522efa7f71SHong Zhang   for (i=0; i<mbs; i++) {
13532efa7f71SHong Zhang     /* copy initial fill into linked list */
13542efa7f71SHong Zhang     nzi = 0; /* nonzeros for active row i */
13552efa7f71SHong Zhang     nzi = ai[r[i]+1] - ai[r[i]];
13562efa7f71SHong Zhang     if (!nzi) SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i);
13572efa7f71SHong Zhang     nzi_al = adiag[r[i]] - ai[r[i]];
13582efa7f71SHong Zhang     nzi_au = ai[r[i]+1] - adiag[r[i]] -1;
13596da515a1SHong Zhang     /* printf("row %d, nzi_al/au %d %d\n",i,nzi_al,nzi_au); */
13602efa7f71SHong Zhang 
13612efa7f71SHong Zhang     /* load in initial unfactored row */
13622efa7f71SHong Zhang     ajtmp = aj + ai[r[i]];
13632efa7f71SHong Zhang     ierr = PetscLLAddPerm(nzi,ajtmp,ic,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
13642efa7f71SHong Zhang     ierr = PetscMemzero(rtmp,(mbs*bs2+1)*sizeof(PetscScalar));CHKERRQ(ierr);
13652efa7f71SHong Zhang     aatmp = a->a + bs2*ai[r[i]];
13662efa7f71SHong Zhang     for (j=0; j<nzi; j++) {
13672efa7f71SHong Zhang       ierr = PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],aatmp+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr);
13682efa7f71SHong Zhang     }
13692efa7f71SHong Zhang 
13702efa7f71SHong Zhang     /* add pivot rows into linked list */
13712efa7f71SHong Zhang     row = lnk[mbs];
13722efa7f71SHong Zhang     while (row < i) {
13732efa7f71SHong Zhang       nzi_bl = bi[row+1] - bi[row] + 1;
13742efa7f71SHong Zhang       bjtmp = bj + bdiag[row+1]+1; /* points to 1st column next to the diagonal in U */
13752efa7f71SHong Zhang       ierr  = PetscLLAddSortedLU(bjtmp,row,nlnk,lnk,lnkbt,i,nzi_bl,im);CHKERRQ(ierr);
13762efa7f71SHong Zhang       nzi  += nlnk;
13772efa7f71SHong Zhang       row   = lnk[row];
13782efa7f71SHong Zhang     }
13792efa7f71SHong Zhang 
13802efa7f71SHong Zhang     /* copy data from lnk into jtmp, then initialize lnk */
13812efa7f71SHong Zhang     ierr = PetscLLClean(mbs,mbs,nzi,lnk,jtmp,lnkbt);CHKERRQ(ierr);
13822efa7f71SHong Zhang 
13832efa7f71SHong Zhang     /* numerical factorization */
13842efa7f71SHong Zhang     bjtmp = jtmp;
13852efa7f71SHong Zhang     row   = *bjtmp++; /* 1st pivot row */
13862efa7f71SHong Zhang 
13872efa7f71SHong Zhang     while  (row < i) {
13882efa7f71SHong Zhang       pc = rtmp + bs2*row;
13892efa7f71SHong Zhang       pv = ba + bs2*bdiag[row]; /* inv(diag) of the pivot row */
13902efa7f71SHong Zhang       Kernel_A_gets_A_times_B(bs,pc,pv,multiplier); /* pc= multiplier = pc*inv(diag[row]) */
13912efa7f71SHong Zhang       ierr = BlockAbs_private(1,bs2,pc,vtmp_abs);CHKERRQ(ierr);
13922efa7f71SHong Zhang       if (vtmp_abs[0] > dt){ /* apply tolerance dropping rule */
13932efa7f71SHong Zhang         pj         = bj + bdiag[row+1] + 1; /* point to 1st entry of U(row,:) */
13942efa7f71SHong Zhang         pv         = ba + bs2*(bdiag[row+1] + 1);
13952efa7f71SHong Zhang         nz         = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:), excluding diagonal */
13962efa7f71SHong Zhang         for (j=0; j<nz; j++){
13972efa7f71SHong Zhang           Kernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j);
13982efa7f71SHong Zhang         }
13996da515a1SHong Zhang         /* ierr = PetscLogFlops(bslog*(nz+1.0)-bs);CHKERRQ(ierr); */
14002efa7f71SHong Zhang       }
14012efa7f71SHong Zhang       row = *bjtmp++;
14022efa7f71SHong Zhang     }
14032efa7f71SHong Zhang 
14042efa7f71SHong Zhang     /* copy sparse rtmp into contiguous vtmp; separate L and U part */
14052efa7f71SHong Zhang     nzi_bl = 0; j = 0;
14062efa7f71SHong Zhang     while (jtmp[j] < i){ /* L-part. Note: jtmp is sorted */
14072efa7f71SHong Zhang       ierr = PetscMemcpy(vtmp+bs2*j,rtmp+bs2*jtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
14082efa7f71SHong Zhang       nzi_bl++; j++;
14092efa7f71SHong Zhang     }
14102efa7f71SHong Zhang     nzi_bu = nzi - nzi_bl -1;
14116da515a1SHong Zhang     /* printf("nzi %d, nzi_bl %d, nzi_bu %d\n",nzi,nzi_bl,nzi_bu); */
14122efa7f71SHong Zhang 
14132efa7f71SHong Zhang     while (j < nzi){ /* U-part */
14142efa7f71SHong Zhang       ierr = PetscMemcpy(vtmp+bs2*j,rtmp+bs2*jtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
14152efa7f71SHong Zhang       /*
14162efa7f71SHong Zhang       printf(" col %d: ",jtmp[j]);
14172efa7f71SHong Zhang       for (j1=0; j1<bs2; j1++) printf(" %g",*(vtmp+bs2*j+j1));
14182efa7f71SHong Zhang       printf(" \n");
14192efa7f71SHong Zhang       */
14202efa7f71SHong Zhang       j++;
14212efa7f71SHong Zhang     }
14222efa7f71SHong Zhang 
14232efa7f71SHong Zhang     ierr = BlockAbs_private(nzi,bs2,vtmp,vtmp_abs);CHKERRQ(ierr);
14242efa7f71SHong Zhang     /*
14252efa7f71SHong Zhang     printf(" row %d, nzi %d, vtmp_abs\n",i,nzi);
14262efa7f71SHong Zhang     for (j1=0; j1<nzi; j1++) printf(" (%d %g),",jtmp[j1],vtmp_abs[j1]);
14272efa7f71SHong Zhang     printf(" \n");
14282efa7f71SHong Zhang     */
14292efa7f71SHong Zhang     bjtmp = bj + bi[i];
14302efa7f71SHong Zhang     batmp = ba + bs2*bi[i];
14312efa7f71SHong Zhang     /* apply level dropping rule to L part */
14322efa7f71SHong Zhang     ncut = nzi_al + dtcount;
14332efa7f71SHong Zhang     if (ncut < nzi_bl){
1434021822bcSHong Zhang       ierr = PetscSortSplitReal(ncut,nzi_bl,vtmp_abs,jtmp);CHKERRQ(ierr);
14352efa7f71SHong Zhang       ierr = PetscSortIntWithScalarArray(ncut,jtmp,vtmp);CHKERRQ(ierr);
14362efa7f71SHong Zhang     } else {
14372efa7f71SHong Zhang       ncut = nzi_bl;
14382efa7f71SHong Zhang     }
14392efa7f71SHong Zhang     for (j=0; j<ncut; j++){
14402efa7f71SHong Zhang       bjtmp[j] = jtmp[j];
14412efa7f71SHong Zhang       ierr = PetscMemcpy(batmp+bs2*j,rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
14422efa7f71SHong Zhang       /*
14432efa7f71SHong Zhang       printf(" col %d: ",bjtmp[j]);
14442efa7f71SHong Zhang       for (j1=0; j1<bs2; j1++) printf(" %g,",*(batmp+bs2*j+j1));
14452efa7f71SHong Zhang       printf("\n");
14462efa7f71SHong Zhang       */
14472efa7f71SHong Zhang     }
14482efa7f71SHong Zhang     bi[i+1] = bi[i] + ncut;
14492efa7f71SHong Zhang     nzi = ncut + 1;
14502efa7f71SHong Zhang 
14512efa7f71SHong Zhang     /* apply level dropping rule to U part */
14522efa7f71SHong Zhang     ncut = nzi_au + dtcount;
14532efa7f71SHong Zhang     if (ncut < nzi_bu){
1454021822bcSHong Zhang       ierr = PetscSortSplitReal(ncut,nzi_bu,vtmp_abs+nzi_bl+1,jtmp+nzi_bl+1);CHKERRQ(ierr);
14552efa7f71SHong Zhang       ierr = PetscSortIntWithScalarArray(ncut,jtmp+nzi_bl+1,vtmp+nzi_bl+1);CHKERRQ(ierr);
14562efa7f71SHong Zhang     } else {
14572efa7f71SHong Zhang       ncut = nzi_bu;
14582efa7f71SHong Zhang     }
14592efa7f71SHong Zhang     nzi += ncut;
14602efa7f71SHong Zhang 
14612efa7f71SHong Zhang     /* mark bdiagonal */
14622efa7f71SHong Zhang     bdiag[i+1]    = bdiag[i] - (ncut + 1);
14636bce7ff8SHong Zhang     bi[2*mbs - i] = bi[2*mbs - i +1] - (ncut + 1);
14646bce7ff8SHong Zhang 
14652efa7f71SHong Zhang     bjtmp = bj + bdiag[i];
14662efa7f71SHong Zhang     batmp = ba + bs2*bdiag[i];
14672efa7f71SHong Zhang     ierr = PetscMemcpy(batmp,rtmp+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
14682efa7f71SHong Zhang     *bjtmp = i;
14692efa7f71SHong Zhang     /*
14702efa7f71SHong Zhang     printf(" diag %d: ",*bjtmp);
14712efa7f71SHong Zhang     for (j=0; j<bs2; j++){
14722efa7f71SHong Zhang       printf(" %g,",batmp[j]);
14732efa7f71SHong Zhang     }
14742efa7f71SHong Zhang     printf("\n");
14752efa7f71SHong Zhang     */
14762efa7f71SHong Zhang     bjtmp = bj + bdiag[i+1]+1;
14772efa7f71SHong Zhang     batmp = ba + (bdiag[i+1]+1)*bs2;
14782efa7f71SHong Zhang 
14792efa7f71SHong Zhang     for (k=0; k<ncut; k++){
14802efa7f71SHong Zhang       bjtmp[k] = jtmp[nzi_bl+1+k];
14812efa7f71SHong Zhang       ierr = PetscMemcpy(batmp+bs2*k,rtmp+bs2*bjtmp[k],bs2*sizeof(MatScalar));CHKERRQ(ierr);
14822efa7f71SHong Zhang       /*
14832efa7f71SHong Zhang       printf(" col %d:",bjtmp[k]);
14842efa7f71SHong Zhang       for (j1=0; j1<bs2; j1++) printf(" %g,",*(batmp+bs2*k+j1));
14852efa7f71SHong Zhang       printf("\n");
14862efa7f71SHong Zhang       */
14872efa7f71SHong Zhang     }
14882efa7f71SHong Zhang 
14892efa7f71SHong Zhang     im[i] = nzi; /* used by PetscLLAddSortedLU() */
14902efa7f71SHong Zhang 
14912efa7f71SHong Zhang     /* invert diagonal block for simplier triangular solves - add shift??? */
14922efa7f71SHong Zhang     batmp = ba + bs2*bdiag[i];
14932efa7f71SHong Zhang     ierr = Kernel_A_gets_inverse_A(bs,batmp,v_pivots,v_work);CHKERRQ(ierr);
14942efa7f71SHong Zhang   } /* for (i=0; i<mbs; i++) */
14952efa7f71SHong Zhang 
14966da515a1SHong Zhang   /* printf("end of L %d, beginning of U %d\n",bi[mbs],bdiag[mbs]); */
14972efa7f71SHong Zhang   if (bi[mbs] >= bdiag[mbs]) SETERRQ2(PETSC_ERR_ARG_SIZ,"end of L array %d cannot >= the beginning of U array %d",bi[mbs],bdiag[mbs]);
14982efa7f71SHong Zhang 
14992efa7f71SHong Zhang   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
15002efa7f71SHong Zhang   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
15012efa7f71SHong Zhang 
15022efa7f71SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
15032efa7f71SHong Zhang 
15042efa7f71SHong Zhang   ierr = PetscFree(im);CHKERRQ(ierr);
15052efa7f71SHong Zhang   ierr = PetscFree(rtmp);CHKERRQ(ierr);
15062efa7f71SHong Zhang 
15072efa7f71SHong Zhang   ierr = PetscLogFlops(bs2*B->cmap->n);CHKERRQ(ierr);
15082efa7f71SHong Zhang   b->maxnz = b->nz = bi[mbs] + bdiag[0] - bdiag[mbs];
15092efa7f71SHong Zhang 
15102efa7f71SHong Zhang   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
15112efa7f71SHong Zhang   ierr = ISIdentity(isicol,&icol_identity);CHKERRQ(ierr);
15122efa7f71SHong Zhang   both_identity = (PetscTruth) (row_identity && icol_identity);
15132efa7f71SHong Zhang   if (row_identity && icol_identity) {
151484a281e5SHong Zhang     B->ops->solve = MatSolve_SeqBAIJ_N_NaturalOrdering_newdatastruct;
15152efa7f71SHong Zhang   } else {
151684a281e5SHong Zhang     B->ops->solve = MatSolve_SeqBAIJ_N_newdatastruct;
15172efa7f71SHong Zhang   }
15182efa7f71SHong Zhang 
15192efa7f71SHong Zhang   B->ops->solveadd          = 0;
15202efa7f71SHong Zhang   B->ops->solvetranspose    = 0;
15212efa7f71SHong Zhang   B->ops->solvetransposeadd = 0;
15222efa7f71SHong Zhang   B->ops->matsolve          = 0;
15232efa7f71SHong Zhang   B->assembled              = PETSC_TRUE;
15242efa7f71SHong Zhang   B->preallocated           = PETSC_TRUE;
1525c8342467SHong Zhang   PetscFunctionReturn(0);
1526c8342467SHong Zhang }
1527