xref: /petsc/src/mat/impls/sbaij/seq/sbaijfact.c (revision 910cf402deff9031709e244fa7c1e84ef2fdc3da)
149b5e25fSSatish Balay 
2c6db04a5SJed Brown #include <../src/mat/impls/baij/seq/baij.h>
3c6db04a5SJed Brown #include <../src/mat/impls/sbaij/seq/sbaij.h>
4af0996ceSBarry Smith #include <petsc/private/kernels/blockinvert.h>
5c6db04a5SJed Brown #include <petscis.h>
649b5e25fSSatish Balay 
74ff4e45dSHong Zhang PetscErrorCode MatGetInertia_SeqSBAIJ(Mat F,PetscInt *nneg,PetscInt *nzero,PetscInt *npos)
85f9f512dSHong Zhang {
94ff4e45dSHong Zhang   Mat_SeqSBAIJ *fact=(Mat_SeqSBAIJ*)F->data;
104ff4e45dSHong Zhang   MatScalar    *dd=fact->a;
114ff4e45dSHong Zhang   PetscInt     mbs=fact->mbs,bs=F->rmap->bs,i,nneg_tmp,npos_tmp,*fi=fact->diag;
125f9f512dSHong Zhang 
135f9f512dSHong Zhang   PetscFunctionBegin;
14e32f2f54SBarry Smith   if (bs != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for bs: %D >1 yet",bs);
157b6c816cSBarry Smith   if (F->factorerrortype==MAT_FACTOR_NUMERIC_ZEROPIVOT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatFactor fails with numeric zeropivot");
164ff4e45dSHong Zhang 
174ff4e45dSHong Zhang   nneg_tmp = 0; npos_tmp = 0;
18eeeff2ecSHong Zhang   for (i=0; i<mbs; i++) {
1926fbe8dcSKarl Rupp     if (PetscRealPart(dd[*fi]) > 0.0) npos_tmp++;
204ff4e45dSHong Zhang     else if (PetscRealPart(dd[*fi]) < 0.0) nneg_tmp++;
21eeeff2ecSHong Zhang     fi++;
223e0d88b5SBarry Smith   }
234ff4e45dSHong Zhang   if (nneg)  *nneg  = nneg_tmp;
24eeeff2ecSHong Zhang   if (npos)  *npos  = npos_tmp;
254ff4e45dSHong Zhang   if (nzero) *nzero = mbs - nneg_tmp - npos_tmp;
265f9f512dSHong Zhang   PetscFunctionReturn(0);
275f9f512dSHong Zhang }
285f9f512dSHong Zhang 
295f9f512dSHong Zhang /*
305f9f512dSHong Zhang   Symbolic U^T*D*U factorization for SBAIJ format. Modified from SSF of YSMP.
3110c27e3dSHong Zhang   Use Modified Sparse Row (MSR) storage for u and ju. See page 85, "Iterative Methods ..." by Saad.
325f9f512dSHong Zhang */
330481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ_MSR(Mat F,Mat A,IS perm,const MatFactorInfo *info)
3410c27e3dSHong Zhang {
3510c27e3dSHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b;
3610c27e3dSHong Zhang   PetscErrorCode ierr;
375d0c19d7SBarry Smith   const PetscInt *rip,*ai,*aj;
385d0c19d7SBarry Smith   PetscInt       i,mbs = a->mbs,*jutmp,bs = A->rmap->bs,bs2=a->bs2;
3910c27e3dSHong Zhang   PetscInt       m,reallocs = 0,prow;
4010c27e3dSHong Zhang   PetscInt       *jl,*q,jmin,jmax,juidx,nzk,qm,*iu,*ju,k,j,vj,umax,maxadd;
4110c27e3dSHong Zhang   PetscReal      f = info->fill;
42ace3abfcSBarry Smith   PetscBool      perm_identity;
4310c27e3dSHong Zhang 
4410c27e3dSHong Zhang   PetscFunctionBegin;
4510c27e3dSHong Zhang   /* check whether perm is the identity mapping */
4610c27e3dSHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
4710c27e3dSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
4810c27e3dSHong Zhang 
4910c27e3dSHong Zhang   if (perm_identity) { /* without permutation */
5010c27e3dSHong Zhang     a->permute = PETSC_FALSE;
5126fbe8dcSKarl Rupp 
5210c27e3dSHong Zhang     ai = a->i; aj = a->j;
5310c27e3dSHong Zhang   } else {            /* non-trivial permutation */
5410c27e3dSHong Zhang     a->permute = PETSC_TRUE;
5526fbe8dcSKarl Rupp 
5610c27e3dSHong Zhang     ierr = MatReorderingSeqSBAIJ(A,perm);CHKERRQ(ierr);
5726fbe8dcSKarl Rupp 
5810c27e3dSHong Zhang     ai = a->inew; aj = a->jnew;
5910c27e3dSHong Zhang   }
6010c27e3dSHong Zhang 
6110c27e3dSHong Zhang   /* initialization */
62854ce69bSBarry Smith   ierr  = PetscMalloc1(mbs+1,&iu);CHKERRQ(ierr);
6310c27e3dSHong Zhang   umax  = (PetscInt)(f*ai[mbs] + 1); umax += mbs + 1;
64785e854fSJed Brown   ierr  = PetscMalloc1(umax,&ju);CHKERRQ(ierr);
6510c27e3dSHong Zhang   iu[0] = mbs+1;
6610c27e3dSHong Zhang   juidx = mbs + 1; /* index for ju */
67d8c74875SBarry Smith   /* jl linked list for pivot row -- linked list for col index */
68dcca6d9dSJed Brown   ierr = PetscMalloc2(mbs,&jl,mbs,&q);CHKERRQ(ierr);
6910c27e3dSHong Zhang   for (i=0; i<mbs; i++) {
7010c27e3dSHong Zhang     jl[i] = mbs;
7110c27e3dSHong Zhang     q[i]  = 0;
7210c27e3dSHong Zhang   }
7310c27e3dSHong Zhang 
7410c27e3dSHong Zhang   /* for each row k */
7510c27e3dSHong Zhang   for (k=0; k<mbs; k++) {
7610c27e3dSHong Zhang     for (i=0; i<mbs; i++) q[i] = 0;  /* to be removed! */
7710c27e3dSHong Zhang     nzk  = 0; /* num. of nz blocks in k-th block row with diagonal block excluded */
7810c27e3dSHong Zhang     q[k] = mbs;
7910c27e3dSHong Zhang     /* initialize nonzero structure of k-th row to row rip[k] of A */
8010c27e3dSHong Zhang     jmin = ai[rip[k]] +1; /* exclude diag[k] */
8110c27e3dSHong Zhang     jmax = ai[rip[k]+1];
8210c27e3dSHong Zhang     for (j=jmin; j<jmax; j++) {
8310c27e3dSHong Zhang       vj = rip[aj[j]]; /* col. value */
8410c27e3dSHong Zhang       if (vj > k) {
8510c27e3dSHong Zhang         qm = k;
8610c27e3dSHong Zhang         do {
8710c27e3dSHong Zhang           m = qm; qm = q[m];
8810c27e3dSHong Zhang         } while (qm < vj);
89535b19f3SBarry Smith         if (qm == vj) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Duplicate entry in A\n");
9010c27e3dSHong Zhang         nzk++;
9110c27e3dSHong Zhang         q[m]  = vj;
9210c27e3dSHong Zhang         q[vj] = qm;
9310c27e3dSHong Zhang       } /* if (vj > k) */
9410c27e3dSHong Zhang     } /* for (j=jmin; j<jmax; j++) */
9510c27e3dSHong Zhang 
9610c27e3dSHong Zhang     /* modify nonzero structure of k-th row by computing fill-in
9710c27e3dSHong Zhang        for each row i to be merged in */
9810c27e3dSHong Zhang     prow = k;
9910c27e3dSHong Zhang     prow = jl[prow]; /* next pivot row (== mbs for symbolic factorization) */
10010c27e3dSHong Zhang 
10110c27e3dSHong Zhang     while (prow < k) {
10210c27e3dSHong Zhang       /* merge row prow into k-th row */
10310c27e3dSHong Zhang       jmin = iu[prow] + 1; jmax = iu[prow+1];
10410c27e3dSHong Zhang       qm   = k;
10510c27e3dSHong Zhang       for (j=jmin; j<jmax; j++) {
10610c27e3dSHong Zhang         vj = ju[j];
10710c27e3dSHong Zhang         do {
10810c27e3dSHong Zhang           m = qm; qm = q[m];
10910c27e3dSHong Zhang         } while (qm < vj);
11010c27e3dSHong Zhang         if (qm != vj) {
11110c27e3dSHong Zhang           nzk++; q[m] = vj; q[vj] = qm; qm = vj;
11210c27e3dSHong Zhang         }
11310c27e3dSHong Zhang       }
11410c27e3dSHong Zhang       prow = jl[prow]; /* next pivot row */
11510c27e3dSHong Zhang     }
11610c27e3dSHong Zhang 
11710c27e3dSHong Zhang     /* add k to row list for first nonzero element in k-th row */
11810c27e3dSHong Zhang     if (nzk > 0) {
11910c27e3dSHong Zhang       i     = q[k]; /* col value of first nonzero element in U(k, k+1:mbs-1) */
12010c27e3dSHong Zhang       jl[k] = jl[i]; jl[i] = k;
12110c27e3dSHong Zhang     }
12210c27e3dSHong Zhang     iu[k+1] = iu[k] + nzk;
12310c27e3dSHong Zhang 
12410c27e3dSHong Zhang     /* allocate more space to ju if needed */
12510c27e3dSHong Zhang     if (iu[k+1] > umax) {
12610c27e3dSHong Zhang       /* estimate how much additional space we will need */
12710c27e3dSHong Zhang       /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */
12810c27e3dSHong Zhang       /* just double the memory each time */
12910c27e3dSHong Zhang       maxadd = umax;
13010c27e3dSHong Zhang       if (maxadd < nzk) maxadd = (mbs-k)*(nzk+1)/2;
13110c27e3dSHong Zhang       umax += maxadd;
13210c27e3dSHong Zhang 
13310c27e3dSHong Zhang       /* allocate a longer ju */
134785e854fSJed Brown       ierr = PetscMalloc1(umax,&jutmp);CHKERRQ(ierr);
135580bdb30SBarry Smith       ierr = PetscArraycpy(jutmp,ju,iu[k]);CHKERRQ(ierr);
13610c27e3dSHong Zhang       ierr = PetscFree(ju);CHKERRQ(ierr);
13710c27e3dSHong Zhang       ju   = jutmp;
13810c27e3dSHong Zhang       reallocs++; /* count how many times we realloc */
13910c27e3dSHong Zhang     }
14010c27e3dSHong Zhang 
14110c27e3dSHong Zhang     /* save nonzero structure of k-th row in ju */
14210c27e3dSHong Zhang     i=k;
14310c27e3dSHong Zhang     while (nzk--) {
14410c27e3dSHong Zhang       i           = q[i];
14510c27e3dSHong Zhang       ju[juidx++] = i;
14610c27e3dSHong Zhang     }
14710c27e3dSHong Zhang   }
14810c27e3dSHong Zhang 
1496cf91177SBarry Smith #if defined(PETSC_USE_INFO)
15010c27e3dSHong Zhang   if (ai[mbs] != 0) {
15110c27e3dSHong Zhang     PetscReal af = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]);
15257622a8eSBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)f,(double)af);CHKERRQ(ierr);
15357622a8eSBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);CHKERRQ(ierr);
15457622a8eSBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%g);\n",(double)af);CHKERRQ(ierr);
155ae15b995SBarry Smith     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
15610c27e3dSHong Zhang   } else {
157dbf6bb8dSprj-     ierr = PetscInfo(A,"Empty matrix\n");CHKERRQ(ierr);
15810c27e3dSHong Zhang   }
15963ba0a88SBarry Smith #endif
16010c27e3dSHong Zhang 
16110c27e3dSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
162d8c74875SBarry Smith   ierr = PetscFree2(jl,q);CHKERRQ(ierr);
16310c27e3dSHong Zhang 
16410c27e3dSHong Zhang   /* put together the new matrix */
165367daffbSBarry Smith   ierr = MatSeqSBAIJSetPreallocation(F,bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr);
16610c27e3dSHong Zhang 
1673bb1ff40SBarry Smith   /* ierr = PetscLogObjectParent((PetscObject)B,(PetscObject)iperm);CHKERRQ(ierr); */
168719d5645SBarry Smith   b                = (Mat_SeqSBAIJ*)(F)->data;
16910c27e3dSHong Zhang   b->singlemalloc  = PETSC_FALSE;
170e6b907acSBarry Smith   b->free_a        = PETSC_TRUE;
171e6b907acSBarry Smith   b->free_ij       = PETSC_TRUE;
17226fbe8dcSKarl Rupp 
173854ce69bSBarry Smith   ierr    = PetscMalloc1((iu[mbs]+1)*bs2,&b->a);CHKERRQ(ierr);
17410c27e3dSHong Zhang   b->j    = ju;
17510c27e3dSHong Zhang   b->i    = iu;
17610c27e3dSHong Zhang   b->diag = 0;
17710c27e3dSHong Zhang   b->ilen = 0;
17810c27e3dSHong Zhang   b->imax = 0;
17910c27e3dSHong Zhang   b->row  = perm;
18026fbe8dcSKarl Rupp 
18110c27e3dSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
18226fbe8dcSKarl Rupp 
18310c27e3dSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
18426fbe8dcSKarl Rupp 
18510c27e3dSHong Zhang   b->icol = perm;
18610c27e3dSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
187854ce69bSBarry Smith   ierr    = PetscMalloc1(bs*mbs+bs,&b->solve_work);CHKERRQ(ierr);
18810c27e3dSHong Zhang   /* In b structure:  Free imax, ilen, old a, old j.
18910c27e3dSHong Zhang      Allocate idnew, solve_work, new a, new j */
1903bb1ff40SBarry Smith   ierr     = PetscLogObjectMemory((PetscObject)F,(iu[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
19110c27e3dSHong Zhang   b->maxnz = b->nz = iu[mbs];
19210c27e3dSHong Zhang 
193719d5645SBarry Smith   (F)->info.factor_mallocs   = reallocs;
194719d5645SBarry Smith   (F)->info.fill_ratio_given = f;
19510c27e3dSHong Zhang   if (ai[mbs] != 0) {
196719d5645SBarry Smith     (F)->info.fill_ratio_needed = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]);
19710c27e3dSHong Zhang   } else {
198719d5645SBarry Smith     (F)->info.fill_ratio_needed = 0.0;
19910c27e3dSHong Zhang   }
200d595f711SHong Zhang   ierr = MatSeqSBAIJSetNumericFactorization_inplace(F,perm_identity);CHKERRQ(ierr);
20110c27e3dSHong Zhang   PetscFunctionReturn(0);
20210c27e3dSHong Zhang }
20310c27e3dSHong Zhang /*
20410c27e3dSHong Zhang     Symbolic U^T*D*U factorization for SBAIJ format.
20580d3d5a7SHong Zhang     See MatICCFactorSymbolic_SeqAIJ() for description of its data structure.
20610c27e3dSHong Zhang */
207c6db04a5SJed Brown #include <petscbt.h>
208c6db04a5SJed Brown #include <../src/mat/utils/freespace.h>
2090481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
21049b5e25fSSatish Balay {
21198a8e78dSHong Zhang   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data;
21298a8e78dSHong Zhang   Mat_SeqSBAIJ       *b;
2136849ba73SBarry Smith   PetscErrorCode     ierr;
214ace3abfcSBarry Smith   PetscBool          perm_identity,missing;
21598a8e78dSHong Zhang   PetscReal          fill = info->fill;
2167b056e98SHong Zhang   const PetscInt     *rip,*ai=a->i,*aj=a->j;
2179186b105SHong Zhang   PetscInt           i,mbs=a->mbs,bs=A->rmap->bs,reallocs=0,prow;
21898a8e78dSHong Zhang   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
219c6d0d4f0SHong Zhang   PetscInt           nlnk,*lnk,ncols,*cols,*uj,**ui_ptr,*uj_ptr,*udiag;
2200298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
22198a8e78dSHong Zhang   PetscBT            lnkbt;
222d595f711SHong Zhang 
223d595f711SHong Zhang   PetscFunctionBegin;
2249186b105SHong Zhang   if (A->rmap->n != A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Must be square matrix, rows %D columns %D",A->rmap->n,A->cmap->n);
2259186b105SHong Zhang   ierr = MatMissingDiagonal(A,&missing,&i);CHKERRQ(ierr);
2269186b105SHong Zhang   if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i);
2276d07c2b1SHong Zhang   if (bs > 1) {
2286d07c2b1SHong Zhang     ierr = MatCholeskyFactorSymbolic_SeqSBAIJ_inplace(fact,A,perm,info);CHKERRQ(ierr);
2296d07c2b1SHong Zhang     PetscFunctionReturn(0);
2306d07c2b1SHong Zhang   }
231d595f711SHong Zhang 
232d595f711SHong Zhang   /* check whether perm is the identity mapping */
233d595f711SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
234e32f2f54SBarry Smith   if (!perm_identity) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix reordering is not supported for sbaij matrix. Use aij format");
235d595f711SHong Zhang   a->permute = PETSC_FALSE;
236d595f711SHong Zhang   ierr       = ISGetIndices(perm,&rip);CHKERRQ(ierr);
237d595f711SHong Zhang 
238d595f711SHong Zhang   /* initialization */
239854ce69bSBarry Smith   ierr  = PetscMalloc1(mbs+1,&ui);CHKERRQ(ierr);
240854ce69bSBarry Smith   ierr  = PetscMalloc1(mbs+1,&udiag);CHKERRQ(ierr);
241d595f711SHong Zhang   ui[0] = 0;
242d595f711SHong Zhang 
243d595f711SHong Zhang   /* jl: linked list for storing indices of the pivot rows
244d595f711SHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
245dcca6d9dSJed Brown   ierr = PetscMalloc4(mbs,&ui_ptr,mbs,&il,mbs,&jl,mbs,&cols);CHKERRQ(ierr);
246d595f711SHong Zhang   for (i=0; i<mbs; i++) {
247d595f711SHong Zhang     jl[i] = mbs; il[i] = 0;
248d595f711SHong Zhang   }
249d595f711SHong Zhang 
250d595f711SHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
251d595f711SHong Zhang   nlnk = mbs + 1;
252d595f711SHong Zhang   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
253d595f711SHong Zhang 
254d595f711SHong Zhang   /* initial FreeSpace size is fill*(ai[mbs]+1) */
255f91af8c7SBarry Smith   ierr          = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,ai[mbs]+1),&free_space);CHKERRQ(ierr);
256d595f711SHong Zhang   current_space = free_space;
257d595f711SHong Zhang 
258d595f711SHong Zhang   for (k=0; k<mbs; k++) {  /* for each active row k */
259d595f711SHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
260d595f711SHong Zhang     nzk   = 0;
261c6d0d4f0SHong Zhang     ncols = ai[k+1] - ai[k];
262e32f2f54SBarry Smith     if (!ncols) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_CH_ZRPVT,"Empty row %D in matrix ",k);
263d595f711SHong Zhang     for (j=0; j<ncols; j++) {
264c6d0d4f0SHong Zhang       i       = *(aj + ai[k] + j);
265c6d0d4f0SHong Zhang       cols[j] = i;
266d595f711SHong Zhang     }
267d595f711SHong Zhang     ierr = PetscLLAdd(ncols,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
268d595f711SHong Zhang     nzk += nlnk;
269d595f711SHong Zhang 
270d595f711SHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
271d595f711SHong Zhang     prow = jl[k]; /* 1st pivot row */
272d595f711SHong Zhang 
273d595f711SHong Zhang     while (prow < k) {
274d595f711SHong Zhang       nextprow = jl[prow];
275d595f711SHong Zhang       /* merge prow into k-th row */
276d595f711SHong Zhang       jmin   = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
277d595f711SHong Zhang       jmax   = ui[prow+1];
278d595f711SHong Zhang       ncols  = jmax-jmin;
279d595f711SHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
280d595f711SHong Zhang       ierr   = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
281d595f711SHong Zhang       nzk   += nlnk;
282d595f711SHong Zhang 
283d595f711SHong Zhang       /* update il and jl for prow */
284d595f711SHong Zhang       if (jmin < jmax) {
285d595f711SHong Zhang         il[prow] = jmin;
286d595f711SHong Zhang         j        = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
287d595f711SHong Zhang       }
288d595f711SHong Zhang       prow = nextprow;
289d595f711SHong Zhang     }
290d595f711SHong Zhang 
291d595f711SHong Zhang     /* if free space is not available, make more free space */
292d595f711SHong Zhang     if (current_space->local_remaining<nzk) {
293d595f711SHong Zhang       i    = mbs - k + 1; /* num of unfactored rows */
294f91af8c7SBarry Smith       i    = PetscIntMultTruncate(i,PetscMin(nzk, i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
295d595f711SHong Zhang       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
296d595f711SHong Zhang       reallocs++;
297d595f711SHong Zhang     }
298d595f711SHong Zhang 
299d595f711SHong Zhang     /* copy data into free space, then initialize lnk */
300d595f711SHong Zhang     ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
301d595f711SHong Zhang 
302d595f711SHong Zhang     /* add the k-th row into il and jl */
3037b056e98SHong Zhang     if (nzk > 1) {
304d595f711SHong Zhang       i     = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
305d595f711SHong Zhang       jl[k] = jl[i]; jl[i] = k;
306d595f711SHong Zhang       il[k] = ui[k] + 1;
307d595f711SHong Zhang     }
308d595f711SHong Zhang     ui_ptr[k] = current_space->array;
30926fbe8dcSKarl Rupp 
310d595f711SHong Zhang     current_space->array           += nzk;
311d595f711SHong Zhang     current_space->local_used      += nzk;
312d595f711SHong Zhang     current_space->local_remaining -= nzk;
313d595f711SHong Zhang 
314d595f711SHong Zhang     ui[k+1] = ui[k] + nzk;
315d595f711SHong Zhang   }
316d595f711SHong Zhang 
317d595f711SHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
318d595f711SHong Zhang   ierr = PetscFree4(ui_ptr,il,jl,cols);CHKERRQ(ierr);
319d595f711SHong Zhang 
320d595f711SHong Zhang   /* destroy list of free space and other temporary array(s) */
321854ce69bSBarry Smith   ierr = PetscMalloc1(ui[mbs]+1,&uj);CHKERRQ(ierr);
3227b056e98SHong Zhang   ierr = PetscFreeSpaceContiguous_Cholesky(&free_space,uj,mbs,ui,udiag);CHKERRQ(ierr); /* store matrix factor */
323d595f711SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
324d595f711SHong Zhang 
325d595f711SHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
326367daffbSBarry Smith   ierr = MatSeqSBAIJSetPreallocation(fact,bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr);
327d595f711SHong Zhang 
3287b056e98SHong Zhang   b               = (Mat_SeqSBAIJ*)fact->data;
329d595f711SHong Zhang   b->singlemalloc = PETSC_FALSE;
330d595f711SHong Zhang   b->free_a       = PETSC_TRUE;
331d595f711SHong Zhang   b->free_ij      = PETSC_TRUE;
33226fbe8dcSKarl Rupp 
333854ce69bSBarry Smith   ierr = PetscMalloc1(ui[mbs]+1,&b->a);CHKERRQ(ierr);
33426fbe8dcSKarl Rupp 
335d595f711SHong Zhang   b->j         = uj;
336d595f711SHong Zhang   b->i         = ui;
337c6d0d4f0SHong Zhang   b->diag      = udiag;
338c6d0d4f0SHong Zhang   b->free_diag = PETSC_TRUE;
339d595f711SHong Zhang   b->ilen      = 0;
340d595f711SHong Zhang   b->imax      = 0;
341d595f711SHong Zhang   b->row       = perm;
342d595f711SHong Zhang   b->icol      = perm;
34326fbe8dcSKarl Rupp 
344d595f711SHong Zhang   ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
3457b056e98SHong Zhang   ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
34626fbe8dcSKarl Rupp 
3477b056e98SHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
34826fbe8dcSKarl Rupp 
349854ce69bSBarry Smith   ierr = PetscMalloc1(mbs+1,&b->solve_work);CHKERRQ(ierr);
3503bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)fact,ui[mbs]*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
35126fbe8dcSKarl Rupp 
352d595f711SHong Zhang   b->maxnz = b->nz = ui[mbs];
353d595f711SHong Zhang 
35495b5ac22SHong Zhang   fact->info.factor_mallocs   = reallocs;
35595b5ac22SHong Zhang   fact->info.fill_ratio_given = fill;
356d595f711SHong Zhang   if (ai[mbs] != 0) {
35795b5ac22SHong Zhang     fact->info.fill_ratio_needed = ((PetscReal)ui[mbs])/ai[mbs];
358d595f711SHong Zhang   } else {
35995b5ac22SHong Zhang     fact->info.fill_ratio_needed = 0.0;
360d595f711SHong Zhang   }
36195b5ac22SHong Zhang #if defined(PETSC_USE_INFO)
36295b5ac22SHong Zhang   if (ai[mbs] != 0) {
36395b5ac22SHong Zhang     PetscReal af = fact->info.fill_ratio_needed;
36457622a8eSBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)fill,(double)af);CHKERRQ(ierr);
36557622a8eSBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);CHKERRQ(ierr);
36657622a8eSBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%g) for best performance.\n",(double)af);CHKERRQ(ierr);
36795b5ac22SHong Zhang   } else {
368dbf6bb8dSprj-     ierr = PetscInfo(A,"Empty matrix\n");CHKERRQ(ierr);
36995b5ac22SHong Zhang   }
37095b5ac22SHong Zhang #endif
371c6d0d4f0SHong Zhang   fact->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering;
372d595f711SHong Zhang   PetscFunctionReturn(0);
373d595f711SHong Zhang }
374d595f711SHong Zhang 
375d595f711SHong Zhang PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ_inplace(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
376d595f711SHong Zhang {
377d595f711SHong Zhang   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data;
378d595f711SHong Zhang   Mat_SeqSBAIJ       *b;
379d595f711SHong Zhang   PetscErrorCode     ierr;
380ace3abfcSBarry Smith   PetscBool          perm_identity,missing;
381d595f711SHong Zhang   PetscReal          fill = info->fill;
382d595f711SHong Zhang   const PetscInt     *rip,*ai,*aj;
383d595f711SHong Zhang   PetscInt           i,mbs=a->mbs,bs=A->rmap->bs,reallocs=0,prow,d;
384d595f711SHong Zhang   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
385d595f711SHong Zhang   PetscInt           nlnk,*lnk,ncols,*cols,*uj,**ui_ptr,*uj_ptr;
3860298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
387d595f711SHong Zhang   PetscBT            lnkbt;
38849b5e25fSSatish Balay 
38949b5e25fSSatish Balay   PetscFunctionBegin;
39058ebbce7SBarry Smith   ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr);
391e32f2f54SBarry Smith   if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
39258ebbce7SBarry Smith 
39310c27e3dSHong Zhang   /*
39410c27e3dSHong Zhang    This code originally uses Modified Sparse Row (MSR) storage
39510c27e3dSHong Zhang    (see page 85, "Iterative Methods ..." by Saad) for the output matrix B - bad choise!
39610c27e3dSHong Zhang    Then it is rewritten so the factor B takes seqsbaij format. However the associated
39710c27e3dSHong Zhang    MatCholeskyFactorNumeric_() have not been modified for the cases of bs>1 or !perm_identity,
39810c27e3dSHong Zhang    thus the original code in MSR format is still used for these cases.
39910c27e3dSHong Zhang    The code below should replace MatCholeskyFactorSymbolic_SeqSBAIJ_MSR() whenever
40010c27e3dSHong Zhang    MatCholeskyFactorNumeric_() is modified for using sbaij symbolic factor.
40110c27e3dSHong Zhang   */
402fff829cfSHong Zhang   if (bs > 1) {
403719d5645SBarry Smith     ierr = MatCholeskyFactorSymbolic_SeqSBAIJ_MSR(fact,A,perm,info);CHKERRQ(ierr);
404fff829cfSHong Zhang     PetscFunctionReturn(0);
405fff829cfSHong Zhang   }
40610c27e3dSHong Zhang 
40798a8e78dSHong Zhang   /* check whether perm is the identity mapping */
40898a8e78dSHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
40998a8e78dSHong Zhang 
410fff829cfSHong Zhang   if (perm_identity) {
411fff829cfSHong Zhang     a->permute = PETSC_FALSE;
41226fbe8dcSKarl Rupp 
413fff829cfSHong Zhang     ai = a->i; aj = a->j;
4146c4ed002SBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix reordering is not supported for sbaij matrix. Use aij format");
415fff829cfSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
416fff829cfSHong Zhang 
417fff829cfSHong Zhang   /* initialization */
418854ce69bSBarry Smith   ierr  = PetscMalloc1(mbs+1,&ui);CHKERRQ(ierr);
41998a8e78dSHong Zhang   ui[0] = 0;
42098a8e78dSHong Zhang 
42198a8e78dSHong Zhang   /* jl: linked list for storing indices of the pivot rows
4227625dc9aSHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
423dcca6d9dSJed Brown   ierr = PetscMalloc4(mbs,&ui_ptr,mbs,&il,mbs,&jl,mbs,&cols);CHKERRQ(ierr);
4247625dc9aSHong Zhang   for (i=0; i<mbs; i++) {
4257625dc9aSHong Zhang     jl[i] = mbs; il[i] = 0;
426fff829cfSHong Zhang   }
427fff829cfSHong Zhang 
42898a8e78dSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
4297625dc9aSHong Zhang   nlnk = mbs + 1;
4307625dc9aSHong Zhang   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
431fff829cfSHong Zhang 
4327625dc9aSHong Zhang   /* initial FreeSpace size is fill*(ai[mbs]+1) */
433f91af8c7SBarry Smith   ierr          = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,ai[mbs]+1),&free_space);CHKERRQ(ierr);
43498a8e78dSHong Zhang   current_space = free_space;
43598a8e78dSHong Zhang 
4367625dc9aSHong Zhang   for (k=0; k<mbs; k++) {  /* for each active row k */
43798a8e78dSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
43898a8e78dSHong Zhang     nzk   = 0;
43998a8e78dSHong Zhang     ncols = ai[rip[k]+1] - ai[rip[k]];
4407625dc9aSHong Zhang     for (j=0; j<ncols; j++) {
4417625dc9aSHong Zhang       i       = *(aj + ai[rip[k]] + j);
4427625dc9aSHong Zhang       cols[j] = rip[i];
4437625dc9aSHong Zhang     }
4447625dc9aSHong Zhang     ierr = PetscLLAdd(ncols,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
44598a8e78dSHong Zhang     nzk += nlnk;
44698a8e78dSHong Zhang 
44798a8e78dSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
44898a8e78dSHong Zhang     prow = jl[k]; /* 1st pivot row */
449fff829cfSHong Zhang 
450fff829cfSHong Zhang     while (prow < k) {
451fff829cfSHong Zhang       nextprow = jl[prow];
45298a8e78dSHong Zhang       /* merge prow into k-th row */
4537625dc9aSHong Zhang       jmin   = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
45498a8e78dSHong Zhang       jmax   = ui[prow+1];
45598a8e78dSHong Zhang       ncols  = jmax-jmin;
4567625dc9aSHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
4575a8e39fbSHong Zhang       ierr   = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
45898a8e78dSHong Zhang       nzk   += nlnk;
459fff829cfSHong Zhang 
46098a8e78dSHong Zhang       /* update il and jl for prow */
461fff829cfSHong Zhang       if (jmin < jmax) {
462fff829cfSHong Zhang         il[prow] = jmin;
46326fbe8dcSKarl Rupp 
4647625dc9aSHong Zhang         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
465fff829cfSHong Zhang       }
466fff829cfSHong Zhang       prow = nextprow;
467fff829cfSHong Zhang     }
468fff829cfSHong Zhang 
46998a8e78dSHong Zhang     /* if free space is not available, make more free space */
47098a8e78dSHong Zhang     if (current_space->local_remaining<nzk) {
4717625dc9aSHong Zhang       i    = mbs - k + 1; /* num of unfactored rows */
472f91af8c7SBarry Smith       i    = PetscMin(PetscIntMultTruncate(i,nzk), PetscIntMultTruncate(i,i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
473a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
47498a8e78dSHong Zhang       reallocs++;
47598a8e78dSHong Zhang     }
47698a8e78dSHong Zhang 
47798a8e78dSHong Zhang     /* copy data into free space, then initialize lnk */
4787625dc9aSHong Zhang     ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
47998a8e78dSHong Zhang 
48098a8e78dSHong Zhang     /* add the k-th row into il and jl */
48198a8e78dSHong Zhang     if (nzk-1 > 0) {
4827625dc9aSHong Zhang       i     = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
483fff829cfSHong Zhang       jl[k] = jl[i]; jl[i] = k;
48498a8e78dSHong Zhang       il[k] = ui[k] + 1;
485fff829cfSHong Zhang     }
4867625dc9aSHong Zhang     ui_ptr[k] = current_space->array;
48726fbe8dcSKarl Rupp 
48898a8e78dSHong Zhang     current_space->array           += nzk;
48998a8e78dSHong Zhang     current_space->local_used      += nzk;
49098a8e78dSHong Zhang     current_space->local_remaining -= nzk;
491fff829cfSHong Zhang 
49298a8e78dSHong Zhang     ui[k+1] = ui[k] + nzk;
493fff829cfSHong Zhang   }
494fff829cfSHong Zhang 
495fff829cfSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
496d8c74875SBarry Smith   ierr = PetscFree4(ui_ptr,il,jl,cols);CHKERRQ(ierr);
497fff829cfSHong Zhang 
49898a8e78dSHong Zhang   /* destroy list of free space and other temporary array(s) */
499854ce69bSBarry Smith   ierr = PetscMalloc1(ui[mbs]+1,&uj);CHKERRQ(ierr);
500a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
50198a8e78dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
502fff829cfSHong Zhang 
50398a8e78dSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
504367daffbSBarry Smith   ierr = MatSeqSBAIJSetPreallocation(fact,bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr);
50598a8e78dSHong Zhang 
50695b5ac22SHong Zhang   b               = (Mat_SeqSBAIJ*)fact->data;
507fff829cfSHong Zhang   b->singlemalloc = PETSC_FALSE;
508e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
509e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
51026fbe8dcSKarl Rupp 
511854ce69bSBarry Smith   ierr = PetscMalloc1(ui[mbs]+1,&b->a);CHKERRQ(ierr);
51226fbe8dcSKarl Rupp 
51398a8e78dSHong Zhang   b->j    = uj;
51498a8e78dSHong Zhang   b->i    = ui;
515fff829cfSHong Zhang   b->diag = 0;
516fff829cfSHong Zhang   b->ilen = 0;
517fff829cfSHong Zhang   b->imax = 0;
518fff829cfSHong Zhang   b->row  = perm;
51926fbe8dcSKarl Rupp 
520fff829cfSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
52126fbe8dcSKarl Rupp 
522fff829cfSHong Zhang   ierr     = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
523fff829cfSHong Zhang   b->icol  = perm;
524fff829cfSHong Zhang   ierr     = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
525854ce69bSBarry Smith   ierr     = PetscMalloc1(mbs+1,&b->solve_work);CHKERRQ(ierr);
5263bb1ff40SBarry Smith   ierr     = PetscLogObjectMemory((PetscObject)fact,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
5277625dc9aSHong Zhang   b->maxnz = b->nz = ui[mbs];
528fff829cfSHong Zhang 
52995b5ac22SHong Zhang   fact->info.factor_mallocs   = reallocs;
53095b5ac22SHong Zhang   fact->info.fill_ratio_given = fill;
5317625dc9aSHong Zhang   if (ai[mbs] != 0) {
53295b5ac22SHong Zhang     fact->info.fill_ratio_needed = ((PetscReal)ui[mbs])/ai[mbs];
533fff829cfSHong Zhang   } else {
53495b5ac22SHong Zhang     fact->info.fill_ratio_needed = 0.0;
535fff829cfSHong Zhang   }
53695b5ac22SHong Zhang #if defined(PETSC_USE_INFO)
53795b5ac22SHong Zhang   if (ai[mbs] != 0) {
53895b5ac22SHong Zhang     PetscReal af = fact->info.fill_ratio_needed;
53957622a8eSBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)fill,(double)af);CHKERRQ(ierr);
54057622a8eSBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);CHKERRQ(ierr);
54157622a8eSBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%g) for best performance.\n",(double)af);CHKERRQ(ierr);
54295b5ac22SHong Zhang   } else {
543dbf6bb8dSprj-     ierr = PetscInfo(A,"Empty matrix\n");CHKERRQ(ierr);
54495b5ac22SHong Zhang   }
54595b5ac22SHong Zhang #endif
546d595f711SHong Zhang   ierr = MatSeqSBAIJSetNumericFactorization_inplace(fact,perm_identity);CHKERRQ(ierr);
547064503c5SHong Zhang   PetscFunctionReturn(0);
548064503c5SHong Zhang }
549d595f711SHong Zhang 
5500481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat C,Mat A,const MatFactorInfo *info)
55149b5e25fSSatish Balay {
5524c16a6a6SHong Zhang   Mat_SeqSBAIJ   *a   = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data;
5534c16a6a6SHong Zhang   IS             perm = b->row;
5546849ba73SBarry Smith   PetscErrorCode ierr;
5555d0c19d7SBarry Smith   const PetscInt *ai,*aj,*perm_ptr,mbs=a->mbs,*bi=b->i,*bj=b->j;
5565d0c19d7SBarry Smith   PetscInt       i,j;
5575d0c19d7SBarry Smith   PetscInt       *a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
5583bc0b13bSBarry Smith   PetscInt       bs  =A->rmap->bs,bs2 = a->bs2;
5594c16a6a6SHong Zhang   MatScalar      *ba = b->a,*aa,*ap,*dk,*uik;
5604c16a6a6SHong Zhang   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
56128de702eSHong Zhang   MatScalar      *work;
56213f74950SBarry Smith   PetscInt       *pivots;
5635f8bbccaSHong Zhang   PetscBool      allowzeropivot,zeropivotdetected;
5644c16a6a6SHong Zhang 
5654c16a6a6SHong Zhang   PetscFunctionBegin;
5664c16a6a6SHong Zhang   /* initialization */
5671795a4d1SJed Brown   ierr = PetscCalloc1(bs2*mbs,&rtmp);CHKERRQ(ierr);
568dcca6d9dSJed Brown   ierr = PetscMalloc2(mbs,&il,mbs,&jl);CHKERRQ(ierr);
5695f8bbccaSHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
5705f8bbccaSHong Zhang 
5716df5ee2eSHong Zhang   il[0] = 0;
5726df5ee2eSHong Zhang   for (i=0; i<mbs; i++) jl[i] = mbs;
5736df5ee2eSHong Zhang 
574dcca6d9dSJed Brown   ierr = PetscMalloc3(bs2,&dk,bs2,&uik,bs,&work);CHKERRQ(ierr);
575785e854fSJed Brown   ierr = PetscMalloc1(bs,&pivots);CHKERRQ(ierr);
5764c16a6a6SHong Zhang 
5774c16a6a6SHong Zhang   ierr = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr);
5784c16a6a6SHong Zhang 
5794c16a6a6SHong Zhang   /* check permutation */
5804c16a6a6SHong Zhang   if (!a->permute) {
5814c16a6a6SHong Zhang     ai = a->i; aj = a->j; aa = a->a;
5824c16a6a6SHong Zhang   } else {
5834c16a6a6SHong Zhang     ai   = a->inew; aj = a->jnew;
584785e854fSJed Brown     ierr = PetscMalloc1(bs2*ai[mbs],&aa);CHKERRQ(ierr);
585580bdb30SBarry Smith     ierr = PetscArraycpy(aa,a->a,bs2*ai[mbs]);CHKERRQ(ierr);
586785e854fSJed Brown     ierr = PetscMalloc1(ai[mbs],&a2anew);CHKERRQ(ierr);
587580bdb30SBarry Smith     ierr = PetscArraycpy(a2anew,a->a2anew,ai[mbs]);CHKERRQ(ierr);
5884c16a6a6SHong Zhang 
5894c16a6a6SHong Zhang     for (i=0; i<mbs; i++) {
5904c16a6a6SHong Zhang       jmin = ai[i]; jmax = ai[i+1];
5914c16a6a6SHong Zhang       for (j=jmin; j<jmax; j++) {
5924c16a6a6SHong Zhang         while (a2anew[j] != j) {
5934c16a6a6SHong Zhang           k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
5944c16a6a6SHong Zhang           for (k1=0; k1<bs2; k1++) {
5954c16a6a6SHong Zhang             dk[k1]       = aa[k*bs2+k1];
5964c16a6a6SHong Zhang             aa[k*bs2+k1] = aa[j*bs2+k1];
5974c16a6a6SHong Zhang             aa[j*bs2+k1] = dk[k1];
5984c16a6a6SHong Zhang           }
5994c16a6a6SHong Zhang         }
6004c16a6a6SHong Zhang         /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
6014c16a6a6SHong Zhang         if (i > aj[j]) {
6024c16a6a6SHong Zhang           ap = aa + j*bs2;                     /* ptr to the beginning of j-th block of aa */
6034c16a6a6SHong Zhang           for (k=0; k<bs2; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */
6044c16a6a6SHong Zhang           for (k=0; k<bs; k++) {               /* j-th block of aa <- dk^T */
6054c16a6a6SHong Zhang             for (k1=0; k1<bs; k1++) *ap++ = dk[k + bs*k1];
6064c16a6a6SHong Zhang           }
6074c16a6a6SHong Zhang         }
6084c16a6a6SHong Zhang       }
6094c16a6a6SHong Zhang     }
610323b833fSBarry Smith     ierr = PetscFree(a2anew);CHKERRQ(ierr);
6114c16a6a6SHong Zhang   }
6124c16a6a6SHong Zhang 
6134c16a6a6SHong Zhang   /* for each row k */
6144c16a6a6SHong Zhang   for (k = 0; k<mbs; k++) {
6154c16a6a6SHong Zhang 
6164c16a6a6SHong Zhang     /*initialize k-th row with elements nonzero in row perm(k) of A */
6174c16a6a6SHong Zhang     jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];
618057f5ba7SHong Zhang 
6194c16a6a6SHong Zhang     ap = aa + jmin*bs2;
6204c16a6a6SHong Zhang     for (j = jmin; j < jmax; j++) {
6214c16a6a6SHong Zhang       vj       = perm_ptr[aj[j]];   /* block col. index */
6224c16a6a6SHong Zhang       rtmp_ptr = rtmp + vj*bs2;
6234c16a6a6SHong Zhang       for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++;
6244c16a6a6SHong Zhang     }
6254c16a6a6SHong Zhang 
6264c16a6a6SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
627580bdb30SBarry Smith     ierr = PetscArraycpy(dk,rtmp+k*bs2,bs2);CHKERRQ(ierr);
6284c16a6a6SHong Zhang     i    = jl[k]; /* first row to be added to k_th row  */
6294c16a6a6SHong Zhang 
630057f5ba7SHong Zhang     while (i < k) {
6314c16a6a6SHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
6324c16a6a6SHong Zhang 
6334c16a6a6SHong Zhang       /* compute multiplier */
6344c16a6a6SHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
6354c16a6a6SHong Zhang 
6364c16a6a6SHong Zhang       /* uik = -inv(Di)*U_bar(i,k) */
6374c16a6a6SHong Zhang       diag = ba + i*bs2;
6384c16a6a6SHong Zhang       u    = ba + ili*bs2;
639580bdb30SBarry Smith       ierr = PetscArrayzero(uik,bs2);CHKERRQ(ierr);
64096b95a6bSBarry Smith       PetscKernel_A_gets_A_minus_B_times_C(bs,uik,diag,u);
6414c16a6a6SHong Zhang 
6424c16a6a6SHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
64396b95a6bSBarry Smith       PetscKernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u);
6443bc0b13bSBarry Smith       ierr = PetscLogFlops(4.0*bs*bs2);CHKERRQ(ierr);
6454c16a6a6SHong Zhang 
6464c16a6a6SHong Zhang       /* update -U(i,k) */
647580bdb30SBarry Smith       ierr = PetscArraycpy(ba+ili*bs2,uik,bs2);CHKERRQ(ierr);
6484c16a6a6SHong Zhang 
6494c16a6a6SHong Zhang       /* add multiple of row i to k-th row ... */
6504c16a6a6SHong Zhang       jmin = ili + 1; jmax = bi[i+1];
6514c16a6a6SHong Zhang       if (jmin < jmax) {
6524c16a6a6SHong Zhang         for (j=jmin; j<jmax; j++) {
6534c16a6a6SHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j) */
6544c16a6a6SHong Zhang           rtmp_ptr = rtmp + bj[j]*bs2;
6554c16a6a6SHong Zhang           u        = ba + j*bs2;
65696b95a6bSBarry Smith           PetscKernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u);
6574c16a6a6SHong Zhang         }
6583bc0b13bSBarry Smith         ierr = PetscLogFlops(2.0*bs*bs2*(jmax-jmin));CHKERRQ(ierr);
6594c16a6a6SHong Zhang 
6604c16a6a6SHong Zhang         /* ... add i to row list for next nonzero entry */
6614c16a6a6SHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
6624c16a6a6SHong Zhang         j     = bj[jmin];
6634c16a6a6SHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
6644c16a6a6SHong Zhang       }
6654c16a6a6SHong Zhang       i = nexti;
6664c16a6a6SHong Zhang     }
6674c16a6a6SHong Zhang 
6684c16a6a6SHong Zhang     /* save nonzero entries in k-th row of U ... */
6694c16a6a6SHong Zhang 
6704c16a6a6SHong Zhang     /* invert diagonal block */
6714c16a6a6SHong Zhang     diag = ba+k*bs2;
672580bdb30SBarry Smith     ierr = PetscArraycpy(diag,dk,bs2);CHKERRQ(ierr);
6735f8bbccaSHong Zhang 
6745f8bbccaSHong Zhang     ierr = PetscKernel_A_gets_inverse_A(bs,diag,pivots,work,allowzeropivot,&zeropivotdetected);CHKERRQ(ierr);
6757b6c816cSBarry Smith     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
6764c16a6a6SHong Zhang 
6774c16a6a6SHong Zhang     jmin = bi[k]; jmax = bi[k+1];
6784c16a6a6SHong Zhang     if (jmin < jmax) {
6794c16a6a6SHong Zhang       for (j=jmin; j<jmax; j++) {
6804c16a6a6SHong Zhang         vj       = bj[j];      /* block col. index of U */
6814c16a6a6SHong Zhang         u        = ba + j*bs2;
6824c16a6a6SHong Zhang         rtmp_ptr = rtmp + vj*bs2;
6834c16a6a6SHong Zhang         for (k1=0; k1<bs2; k1++) {
6844c16a6a6SHong Zhang           *u++        = *rtmp_ptr;
6854c16a6a6SHong Zhang           *rtmp_ptr++ = 0.0;
6864c16a6a6SHong Zhang         }
6874c16a6a6SHong Zhang       }
6884c16a6a6SHong Zhang 
6894c16a6a6SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
6904c16a6a6SHong Zhang       il[k] = jmin;
6914c16a6a6SHong Zhang       i     = bj[jmin];
6924c16a6a6SHong Zhang       jl[k] = jl[i]; jl[i] = k;
6934c16a6a6SHong Zhang     }
6944c16a6a6SHong Zhang   }
6954c16a6a6SHong Zhang 
6964c16a6a6SHong Zhang   ierr = PetscFree(rtmp);CHKERRQ(ierr);
697d8c74875SBarry Smith   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
698d8c74875SBarry Smith   ierr = PetscFree3(dk,uik,work);CHKERRQ(ierr);
6994c16a6a6SHong Zhang   ierr = PetscFree(pivots);CHKERRQ(ierr);
7004c16a6a6SHong Zhang   if (a->permute) {
7014c16a6a6SHong Zhang     ierr = PetscFree(aa);CHKERRQ(ierr);
7024c16a6a6SHong Zhang   }
7034c16a6a6SHong Zhang 
7044c16a6a6SHong Zhang   ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr);
70526fbe8dcSKarl Rupp 
7064f79d315SHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_N_inplace;
7074f79d315SHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_N_inplace;
7084f79d315SHong Zhang   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_N_inplace;
7094f79d315SHong Zhang   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_N_inplace;
710db4efbfdSBarry Smith 
7114c16a6a6SHong Zhang   C->assembled    = PETSC_TRUE;
7124c16a6a6SHong Zhang   C->preallocated = PETSC_TRUE;
71326fbe8dcSKarl Rupp 
714efee365bSSatish Balay   ierr = PetscLogFlops(1.3333*bs*bs2*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
7154c16a6a6SHong Zhang   PetscFunctionReturn(0);
7164c16a6a6SHong Zhang }
717d4edadd4SHong Zhang 
7180481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info)
719671cb588SHong Zhang {
720671cb588SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data;
721dfbe8321SBarry Smith   PetscErrorCode ierr;
72213f74950SBarry Smith   PetscInt       i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
72313f74950SBarry Smith   PetscInt       *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
7243bc0b13bSBarry Smith   PetscInt       bs  =A->rmap->bs,bs2 = a->bs2;
725671cb588SHong Zhang   MatScalar      *ba = b->a,*aa,*ap,*dk,*uik;
726671cb588SHong Zhang   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
72728de702eSHong Zhang   MatScalar      *work;
72813f74950SBarry Smith   PetscInt       *pivots;
7295f8bbccaSHong Zhang   PetscBool      allowzeropivot,zeropivotdetected;
730671cb588SHong Zhang 
731671cb588SHong Zhang   PetscFunctionBegin;
7321795a4d1SJed Brown   ierr = PetscCalloc1(bs2*mbs,&rtmp);CHKERRQ(ierr);
733dcca6d9dSJed Brown   ierr = PetscMalloc2(mbs,&il,mbs,&jl);CHKERRQ(ierr);
7346df5ee2eSHong Zhang   il[0] = 0;
7356df5ee2eSHong Zhang   for (i=0; i<mbs; i++) jl[i] = mbs;
7366df5ee2eSHong Zhang 
737dcca6d9dSJed Brown   ierr = PetscMalloc3(bs2,&dk,bs2,&uik,bs,&work);CHKERRQ(ierr);
738785e854fSJed Brown   ierr = PetscMalloc1(bs,&pivots);CHKERRQ(ierr);
7395f8bbccaSHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
740671cb588SHong Zhang 
741671cb588SHong Zhang   ai = a->i; aj = a->j; aa = a->a;
742671cb588SHong Zhang 
743671cb588SHong Zhang   /* for each row k */
744671cb588SHong Zhang   for (k = 0; k<mbs; k++) {
745671cb588SHong Zhang 
746671cb588SHong Zhang     /*initialize k-th row with elements nonzero in row k of A */
747671cb588SHong Zhang     jmin = ai[k]; jmax = ai[k+1];
748671cb588SHong Zhang     ap   = aa + jmin*bs2;
749671cb588SHong Zhang     for (j = jmin; j < jmax; j++) {
750671cb588SHong Zhang       vj       = aj[j];   /* block col. index */
751671cb588SHong Zhang       rtmp_ptr = rtmp + vj*bs2;
752671cb588SHong Zhang       for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++;
753671cb588SHong Zhang     }
754671cb588SHong Zhang 
755671cb588SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
756580bdb30SBarry Smith     ierr = PetscArraycpy(dk,rtmp+k*bs2,bs2);CHKERRQ(ierr);
757671cb588SHong Zhang     i    = jl[k]; /* first row to be added to k_th row  */
758671cb588SHong Zhang 
759057f5ba7SHong Zhang     while (i < k) {
760671cb588SHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
761671cb588SHong Zhang 
762671cb588SHong Zhang       /* compute multiplier */
763671cb588SHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
764671cb588SHong Zhang 
765671cb588SHong Zhang       /* uik = -inv(Di)*U_bar(i,k) */
766671cb588SHong Zhang       diag = ba + i*bs2;
767671cb588SHong Zhang       u    = ba + ili*bs2;
768580bdb30SBarry Smith       ierr = PetscArrayzero(uik,bs2);CHKERRQ(ierr);
76996b95a6bSBarry Smith       PetscKernel_A_gets_A_minus_B_times_C(bs,uik,diag,u);
770671cb588SHong Zhang 
771671cb588SHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
77296b95a6bSBarry Smith       PetscKernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u);
7733bc0b13bSBarry Smith       ierr = PetscLogFlops(2.0*bs*bs2);CHKERRQ(ierr);
774671cb588SHong Zhang 
775671cb588SHong Zhang       /* update -U(i,k) */
776580bdb30SBarry Smith       ierr = PetscArraycpy(ba+ili*bs2,uik,bs2);CHKERRQ(ierr);
777671cb588SHong Zhang 
778671cb588SHong Zhang       /* add multiple of row i to k-th row ... */
779671cb588SHong Zhang       jmin = ili + 1; jmax = bi[i+1];
780671cb588SHong Zhang       if (jmin < jmax) {
781671cb588SHong Zhang         for (j=jmin; j<jmax; j++) {
782671cb588SHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j) */
783671cb588SHong Zhang           rtmp_ptr = rtmp + bj[j]*bs2;
784671cb588SHong Zhang           u        = ba + j*bs2;
78596b95a6bSBarry Smith           PetscKernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u);
786671cb588SHong Zhang         }
7873bc0b13bSBarry Smith         ierr = PetscLogFlops(2.0*bs*bs2*(jmax-jmin));CHKERRQ(ierr);
788671cb588SHong Zhang 
789671cb588SHong Zhang         /* ... add i to row list for next nonzero entry */
790671cb588SHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
791671cb588SHong Zhang         j     = bj[jmin];
792671cb588SHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
793671cb588SHong Zhang       }
794671cb588SHong Zhang       i = nexti;
795671cb588SHong Zhang     }
796671cb588SHong Zhang 
797671cb588SHong Zhang     /* save nonzero entries in k-th row of U ... */
798671cb588SHong Zhang 
799671cb588SHong Zhang     /* invert diagonal block */
800671cb588SHong Zhang     diag = ba+k*bs2;
801580bdb30SBarry Smith     ierr = PetscArraycpy(diag,dk,bs2);CHKERRQ(ierr);
8025f8bbccaSHong Zhang 
8035f8bbccaSHong Zhang     ierr = PetscKernel_A_gets_inverse_A(bs,diag,pivots,work,allowzeropivot,&zeropivotdetected);CHKERRQ(ierr);
8047b6c816cSBarry Smith     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
805671cb588SHong Zhang 
806671cb588SHong Zhang     jmin = bi[k]; jmax = bi[k+1];
807671cb588SHong Zhang     if (jmin < jmax) {
808671cb588SHong Zhang       for (j=jmin; j<jmax; j++) {
809671cb588SHong Zhang         vj       = bj[j];      /* block col. index of U */
810671cb588SHong Zhang         u        = ba + j*bs2;
811671cb588SHong Zhang         rtmp_ptr = rtmp + vj*bs2;
812671cb588SHong Zhang         for (k1=0; k1<bs2; k1++) {
813671cb588SHong Zhang           *u++        = *rtmp_ptr;
814671cb588SHong Zhang           *rtmp_ptr++ = 0.0;
815671cb588SHong Zhang         }
816671cb588SHong Zhang       }
817671cb588SHong Zhang 
818671cb588SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
819671cb588SHong Zhang       il[k] = jmin;
820671cb588SHong Zhang       i     = bj[jmin];
821671cb588SHong Zhang       jl[k] = jl[i]; jl[i] = k;
822671cb588SHong Zhang     }
823671cb588SHong Zhang   }
824671cb588SHong Zhang 
825671cb588SHong Zhang   ierr = PetscFree(rtmp);CHKERRQ(ierr);
826d8c74875SBarry Smith   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
827d8c74875SBarry Smith   ierr = PetscFree3(dk,uik,work);CHKERRQ(ierr);
828671cb588SHong Zhang   ierr = PetscFree(pivots);CHKERRQ(ierr);
829671cb588SHong Zhang 
8304f79d315SHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_N_NaturalOrdering_inplace;
8314f79d315SHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_N_NaturalOrdering_inplace;
8324f79d315SHong Zhang   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_N_NaturalOrdering_inplace;
8334f79d315SHong Zhang   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_N_NaturalOrdering_inplace;
834671cb588SHong Zhang   C->assembled           = PETSC_TRUE;
835671cb588SHong Zhang   C->preallocated        = PETSC_TRUE;
83626fbe8dcSKarl Rupp 
837efee365bSSatish Balay   ierr = PetscLogFlops(1.3333*bs*bs2*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
838671cb588SHong Zhang   PetscFunctionReturn(0);
839671cb588SHong Zhang }
840671cb588SHong Zhang 
84149b5e25fSSatish Balay /*
842fcf159c0SHong Zhang     Numeric U^T*D*U factorization for SBAIJ format. Modified from SNF of YSMP.
843cc0c071aSHong Zhang     Version for blocks 2 by 2.
84449b5e25fSSatish Balay */
8450481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat C,Mat A,const MatFactorInfo *info)
84649b5e25fSSatish Balay {
84791602c66SHong Zhang   Mat_SeqSBAIJ   *a   = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data;
848cc0c071aSHong Zhang   IS             perm = b->row;
8496849ba73SBarry Smith   PetscErrorCode ierr;
8505d0c19d7SBarry Smith   const PetscInt *ai,*aj,*perm_ptr;
8515d0c19d7SBarry Smith   PetscInt       i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
8525d0c19d7SBarry Smith   PetscInt       *a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
853d8c74875SBarry Smith   MatScalar      *ba = b->a,*aa,*ap;
854d8c74875SBarry Smith   MatScalar      *u,*diag,*rtmp,*rtmp_ptr,dk[4],uik[4];
85514d2772eSHong Zhang   PetscReal      shift = info->shiftamount;
856a455e926SHong Zhang   PetscBool      allowzeropivot,zeropivotdetected;
85749b5e25fSSatish Balay 
85849b5e25fSSatish Balay   PetscFunctionBegin;
8590164db54SHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
8600164db54SHong Zhang 
86191602c66SHong Zhang   /* initialization */
86291602c66SHong Zhang   /* il and jl record the first nonzero element in each row of the accessing
86391602c66SHong Zhang      window U(0:k, k:mbs-1).
86491602c66SHong Zhang      jl:    list of rows to be added to uneliminated rows
86591602c66SHong Zhang             i>= k: jl(i) is the first row to be added to row i
86691602c66SHong Zhang             i<  k: jl(i) is the row following row i in some list of rows
86791602c66SHong Zhang             jl(i) = mbs indicates the end of a list
86891602c66SHong Zhang      il(i): points to the first nonzero element in columns k,...,mbs-1 of
86991602c66SHong Zhang             row i of U */
8701795a4d1SJed Brown   ierr = PetscCalloc1(4*mbs,&rtmp);CHKERRQ(ierr);
871dcca6d9dSJed Brown   ierr = PetscMalloc2(mbs,&il,mbs,&jl);CHKERRQ(ierr);
8726df5ee2eSHong Zhang   il[0] = 0;
8736df5ee2eSHong Zhang   for (i=0; i<mbs; i++) jl[i] = mbs;
8746df5ee2eSHong Zhang 
875cc0c071aSHong Zhang   ierr = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr);
876cc0c071aSHong Zhang 
877cc0c071aSHong Zhang   /* check permutation */
878cc0c071aSHong Zhang   if (!a->permute) {
879cc0c071aSHong Zhang     ai = a->i; aj = a->j; aa = a->a;
880cc0c071aSHong Zhang   } else {
881cc0c071aSHong Zhang     ai   = a->inew; aj = a->jnew;
882785e854fSJed Brown     ierr = PetscMalloc1(4*ai[mbs],&aa);CHKERRQ(ierr);
883580bdb30SBarry Smith     ierr = PetscArraycpy(aa,a->a,4*ai[mbs]);CHKERRQ(ierr);
884785e854fSJed Brown     ierr = PetscMalloc1(ai[mbs],&a2anew);CHKERRQ(ierr);
885580bdb30SBarry Smith     ierr = PetscArraycpy(a2anew,a->a2anew,ai[mbs]);CHKERRQ(ierr);
886cc0c071aSHong Zhang 
887cc0c071aSHong Zhang     for (i=0; i<mbs; i++) {
888cc0c071aSHong Zhang       jmin = ai[i]; jmax = ai[i+1];
889cc0c071aSHong Zhang       for (j=jmin; j<jmax; j++) {
890cc0c071aSHong Zhang         while (a2anew[j] != j) {
891cc0c071aSHong Zhang           k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
892cc0c071aSHong Zhang           for (k1=0; k1<4; k1++) {
893cc0c071aSHong Zhang             dk[k1]     = aa[k*4+k1];
894cc0c071aSHong Zhang             aa[k*4+k1] = aa[j*4+k1];
895cc0c071aSHong Zhang             aa[j*4+k1] = dk[k1];
896cc0c071aSHong Zhang           }
897cc0c071aSHong Zhang         }
898cc0c071aSHong Zhang         /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
899cc0c071aSHong Zhang         if (i > aj[j]) {
900cc0c071aSHong Zhang           ap    = aa + j*4;  /* ptr to the beginning of the block */
901cc0c071aSHong Zhang           dk[1] = ap[1];     /* swap ap[1] and ap[2] */
902cc0c071aSHong Zhang           ap[1] = ap[2];
903cc0c071aSHong Zhang           ap[2] = dk[1];
904cc0c071aSHong Zhang         }
905cc0c071aSHong Zhang       }
906cc0c071aSHong Zhang     }
907ac355199SBarry Smith     ierr = PetscFree(a2anew);CHKERRQ(ierr);
908cc0c071aSHong Zhang   }
9093845f261SHong Zhang 
91091602c66SHong Zhang   /* for each row k */
91191602c66SHong Zhang   for (k = 0; k<mbs; k++) {
91291602c66SHong Zhang 
91391602c66SHong Zhang     /*initialize k-th row with elements nonzero in row perm(k) of A */
914cc0c071aSHong Zhang     jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];
915cc0c071aSHong Zhang     ap   = aa + jmin*4;
91691602c66SHong Zhang     for (j = jmin; j < jmax; j++) {
917cc0c071aSHong Zhang       vj       = perm_ptr[aj[j]];   /* block col. index */
918cc0c071aSHong Zhang       rtmp_ptr = rtmp + vj*4;
919cc0c071aSHong Zhang       for (i=0; i<4; i++) *rtmp_ptr++ = *ap++;
92091602c66SHong Zhang     }
92191602c66SHong Zhang 
92291602c66SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
923580bdb30SBarry Smith     ierr = PetscArraycpy(dk,rtmp+k*4,4);CHKERRQ(ierr);
92491602c66SHong Zhang     i    = jl[k]; /* first row to be added to k_th row  */
92591602c66SHong Zhang 
926057f5ba7SHong Zhang     while (i < k) {
92791602c66SHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
92891602c66SHong Zhang 
9293845f261SHong Zhang       /* compute multiplier */
93091602c66SHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
9313845f261SHong Zhang 
9323845f261SHong Zhang       /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */
933cc0c071aSHong Zhang       diag   = ba + i*4;
934cc0c071aSHong Zhang       u      = ba + ili*4;
935cc0c071aSHong Zhang       uik[0] = -(diag[0]*u[0] + diag[2]*u[1]);
936cc0c071aSHong Zhang       uik[1] = -(diag[1]*u[0] + diag[3]*u[1]);
937cc0c071aSHong Zhang       uik[2] = -(diag[0]*u[2] + diag[2]*u[3]);
938cc0c071aSHong Zhang       uik[3] = -(diag[1]*u[2] + diag[3]*u[3]);
9393845f261SHong Zhang 
9403845f261SHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */
941cc0c071aSHong Zhang       dk[0] += uik[0]*u[0] + uik[1]*u[1];
942cc0c071aSHong Zhang       dk[1] += uik[2]*u[0] + uik[3]*u[1];
943cc0c071aSHong Zhang       dk[2] += uik[0]*u[2] + uik[1]*u[3];
944cc0c071aSHong Zhang       dk[3] += uik[2]*u[2] + uik[3]*u[3];
9453845f261SHong Zhang 
946dc0b31edSSatish Balay       ierr = PetscLogFlops(16.0*2.0);CHKERRQ(ierr);
947187a9f4bSHong Zhang 
9483845f261SHong Zhang       /* update -U(i,k): ba[ili] = uik */
949580bdb30SBarry Smith       ierr = PetscArraycpy(ba+ili*4,uik,4);CHKERRQ(ierr);
95091602c66SHong Zhang 
95191602c66SHong Zhang       /* add multiple of row i to k-th row ... */
95291602c66SHong Zhang       jmin = ili + 1; jmax = bi[i+1];
95391602c66SHong Zhang       if (jmin < jmax) {
9543845f261SHong Zhang         for (j=jmin; j<jmax; j++) {
9553845f261SHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */
956cc0c071aSHong Zhang           rtmp_ptr     = rtmp + bj[j]*4;
957cc0c071aSHong Zhang           u            = ba + j*4;
958cc0c071aSHong Zhang           rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1];
959cc0c071aSHong Zhang           rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1];
960cc0c071aSHong Zhang           rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3];
961cc0c071aSHong Zhang           rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3];
9623845f261SHong Zhang         }
963dc0b31edSSatish Balay         ierr = PetscLogFlops(16.0*(jmax-jmin));CHKERRQ(ierr);
9643845f261SHong Zhang 
96591602c66SHong Zhang         /* ... add i to row list for next nonzero entry */
96691602c66SHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
96791602c66SHong Zhang         j     = bj[jmin];
96891602c66SHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
96991602c66SHong Zhang       }
970a1723e09SHong Zhang       i = nexti;
97191602c66SHong Zhang     }
972cc0c071aSHong Zhang 
97391602c66SHong Zhang     /* save nonzero entries in k-th row of U ... */
9743845f261SHong Zhang 
975cc0c071aSHong Zhang     /* invert diagonal block */
976cc0c071aSHong Zhang     diag = ba+k*4;
977580bdb30SBarry Smith     ierr = PetscArraycpy(diag,dk,4);CHKERRQ(ierr);
978a455e926SHong Zhang     ierr = PetscKernel_A_gets_inverse_A_2(diag,shift,allowzeropivot,&zeropivotdetected);CHKERRQ(ierr);
9797b6c816cSBarry Smith     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
9803845f261SHong Zhang 
98191602c66SHong Zhang     jmin = bi[k]; jmax = bi[k+1];
98291602c66SHong Zhang     if (jmin < jmax) {
98391602c66SHong Zhang       for (j=jmin; j<jmax; j++) {
984cc0c071aSHong Zhang         vj       = bj[j];      /* block col. index of U */
985cc0c071aSHong Zhang         u        = ba + j*4;
986cc0c071aSHong Zhang         rtmp_ptr = rtmp + vj*4;
987cc0c071aSHong Zhang         for (k1=0; k1<4; k1++) {
988cc0c071aSHong Zhang           *u++        = *rtmp_ptr;
989cc0c071aSHong Zhang           *rtmp_ptr++ = 0.0;
990cc0c071aSHong Zhang         }
991cc0c071aSHong Zhang       }
9923845f261SHong Zhang 
99391602c66SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
99491602c66SHong Zhang       il[k] = jmin;
99591602c66SHong Zhang       i     = bj[jmin];
99691602c66SHong Zhang       jl[k] = jl[i]; jl[i] = k;
99791602c66SHong Zhang     }
99891602c66SHong Zhang   }
9993845f261SHong Zhang 
100049b5e25fSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
1001d8c74875SBarry Smith   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
100291602c66SHong Zhang   if (a->permute) {
100391602c66SHong Zhang     ierr = PetscFree(aa);CHKERRQ(ierr);
100491602c66SHong Zhang   }
1005cc0c071aSHong Zhang   ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr);
100626fbe8dcSKarl Rupp 
10074f79d315SHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_2_inplace;
10084f79d315SHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_2_inplace;
100949b5e25fSSatish Balay   C->assembled           = PETSC_TRUE;
10105c0bcdfcSHong Zhang   C->preallocated        = PETSC_TRUE;
101126fbe8dcSKarl Rupp 
1012efee365bSSatish Balay   ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
101349b5e25fSSatish Balay   PetscFunctionReturn(0);
101449b5e25fSSatish Balay }
101591602c66SHong Zhang 
101649b5e25fSSatish Balay /*
101749b5e25fSSatish Balay       Version for when blocks are 2 by 2 Using natural ordering
101849b5e25fSSatish Balay */
10190481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info)
102049b5e25fSSatish Balay {
1021ab27746eSHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data;
1022dfbe8321SBarry Smith   PetscErrorCode ierr;
102313f74950SBarry Smith   PetscInt       i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
102413f74950SBarry Smith   PetscInt       *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
1025d8c74875SBarry Smith   MatScalar      *ba = b->a,*aa,*ap,dk[8],uik[8];
1026ab27746eSHong Zhang   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
102714d2772eSHong Zhang   PetscReal      shift = info->shiftamount;
1028a455e926SHong Zhang   PetscBool      allowzeropivot,zeropivotdetected;
102949b5e25fSSatish Balay 
103049b5e25fSSatish Balay   PetscFunctionBegin;
10310164db54SHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
10320164db54SHong Zhang 
1033ab27746eSHong Zhang   /* initialization */
1034ab27746eSHong Zhang   /* il and jl record the first nonzero element in each row of the accessing
1035ab27746eSHong Zhang      window U(0:k, k:mbs-1).
1036ab27746eSHong Zhang      jl:    list of rows to be added to uneliminated rows
1037ab27746eSHong Zhang             i>= k: jl(i) is the first row to be added to row i
1038ab27746eSHong Zhang             i<  k: jl(i) is the row following row i in some list of rows
1039ab27746eSHong Zhang             jl(i) = mbs indicates the end of a list
1040ab27746eSHong Zhang      il(i): points to the first nonzero element in columns k,...,mbs-1 of
1041ab27746eSHong Zhang             row i of U */
10421795a4d1SJed Brown   ierr = PetscCalloc1(4*mbs,&rtmp);CHKERRQ(ierr);
1043dcca6d9dSJed Brown   ierr = PetscMalloc2(mbs,&il,mbs,&jl);CHKERRQ(ierr);
10446df5ee2eSHong Zhang   il[0] = 0;
10456df5ee2eSHong Zhang   for (i=0; i<mbs; i++) jl[i] = mbs;
10466df5ee2eSHong Zhang 
1047ab27746eSHong Zhang   ai = a->i; aj = a->j; aa = a->a;
1048ab27746eSHong Zhang 
1049ab27746eSHong Zhang   /* for each row k */
1050ab27746eSHong Zhang   for (k = 0; k<mbs; k++) {
1051ab27746eSHong Zhang 
1052ab27746eSHong Zhang     /*initialize k-th row with elements nonzero in row k of A */
1053ab27746eSHong Zhang     jmin = ai[k]; jmax = ai[k+1];
1054ab27746eSHong Zhang     ap   = aa + jmin*4;
1055ab27746eSHong Zhang     for (j = jmin; j < jmax; j++) {
1056ab27746eSHong Zhang       vj       = aj[j];   /* block col. index */
1057ab27746eSHong Zhang       rtmp_ptr = rtmp + vj*4;
1058ab27746eSHong Zhang       for (i=0; i<4; i++) *rtmp_ptr++ = *ap++;
105949b5e25fSSatish Balay     }
1060ab27746eSHong Zhang 
1061ab27746eSHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
1062580bdb30SBarry Smith     ierr = PetscArraycpy(dk,rtmp+k*4,4);CHKERRQ(ierr);
1063ab27746eSHong Zhang     i    = jl[k]; /* first row to be added to k_th row  */
1064ab27746eSHong Zhang 
1065057f5ba7SHong Zhang     while (i < k) {
1066ab27746eSHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
1067ab27746eSHong Zhang 
1068ab27746eSHong Zhang       /* compute multiplier */
1069ab27746eSHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
1070ab27746eSHong Zhang 
1071ab27746eSHong Zhang       /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */
1072ab27746eSHong Zhang       diag   = ba + i*4;
1073ab27746eSHong Zhang       u      = ba + ili*4;
1074ab27746eSHong Zhang       uik[0] = -(diag[0]*u[0] + diag[2]*u[1]);
1075ab27746eSHong Zhang       uik[1] = -(diag[1]*u[0] + diag[3]*u[1]);
1076ab27746eSHong Zhang       uik[2] = -(diag[0]*u[2] + diag[2]*u[3]);
1077ab27746eSHong Zhang       uik[3] = -(diag[1]*u[2] + diag[3]*u[3]);
1078ab27746eSHong Zhang 
1079ab27746eSHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */
1080ab27746eSHong Zhang       dk[0] += uik[0]*u[0] + uik[1]*u[1];
1081ab27746eSHong Zhang       dk[1] += uik[2]*u[0] + uik[3]*u[1];
1082ab27746eSHong Zhang       dk[2] += uik[0]*u[2] + uik[1]*u[3];
1083ab27746eSHong Zhang       dk[3] += uik[2]*u[2] + uik[3]*u[3];
1084ab27746eSHong Zhang 
1085dc0b31edSSatish Balay       ierr = PetscLogFlops(16.0*2.0);CHKERRQ(ierr);
1086187a9f4bSHong Zhang 
1087ab27746eSHong Zhang       /* update -U(i,k): ba[ili] = uik */
1088580bdb30SBarry Smith       ierr = PetscArraycpy(ba+ili*4,uik,4);CHKERRQ(ierr);
1089ab27746eSHong Zhang 
1090ab27746eSHong Zhang       /* add multiple of row i to k-th row ... */
1091ab27746eSHong Zhang       jmin = ili + 1; jmax = bi[i+1];
1092ab27746eSHong Zhang       if (jmin < jmax) {
1093ab27746eSHong Zhang         for (j=jmin; j<jmax; j++) {
1094ab27746eSHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */
1095ab27746eSHong Zhang           rtmp_ptr     = rtmp + bj[j]*4;
1096ab27746eSHong Zhang           u            = ba + j*4;
1097ab27746eSHong Zhang           rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1];
1098ab27746eSHong Zhang           rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1];
1099ab27746eSHong Zhang           rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3];
1100ab27746eSHong Zhang           rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3];
110149b5e25fSSatish Balay         }
1102dc0b31edSSatish Balay         ierr = PetscLogFlops(16.0*(jmax-jmin));CHKERRQ(ierr);
1103ab27746eSHong Zhang 
1104ab27746eSHong Zhang         /* ... add i to row list for next nonzero entry */
1105ab27746eSHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
1106ab27746eSHong Zhang         j     = bj[jmin];
1107ab27746eSHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
110849b5e25fSSatish Balay       }
1109ab27746eSHong Zhang       i = nexti;
111049b5e25fSSatish Balay     }
1111ab27746eSHong Zhang 
1112ab27746eSHong Zhang     /* save nonzero entries in k-th row of U ... */
1113ab27746eSHong Zhang 
111449b5e25fSSatish Balay     /* invert diagonal block */
1115ab27746eSHong Zhang     diag = ba+k*4;
1116580bdb30SBarry Smith     ierr = PetscArraycpy(diag,dk,4);CHKERRQ(ierr);
1117a455e926SHong Zhang     ierr = PetscKernel_A_gets_inverse_A_2(diag,shift,allowzeropivot,&zeropivotdetected);CHKERRQ(ierr);
11187b6c816cSBarry Smith     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1119ab27746eSHong Zhang 
1120ab27746eSHong Zhang     jmin = bi[k]; jmax = bi[k+1];
1121ab27746eSHong Zhang     if (jmin < jmax) {
1122ab27746eSHong Zhang       for (j=jmin; j<jmax; j++) {
1123ab27746eSHong Zhang         vj       = bj[j];      /* block col. index of U */
1124ab27746eSHong Zhang         u        = ba + j*4;
1125ab27746eSHong Zhang         rtmp_ptr = rtmp + vj*4;
1126ab27746eSHong Zhang         for (k1=0; k1<4; k1++) {
1127ab27746eSHong Zhang           *u++        = *rtmp_ptr;
1128ab27746eSHong Zhang           *rtmp_ptr++ = 0.0;
1129ab27746eSHong Zhang         }
1130ab27746eSHong Zhang       }
1131ab27746eSHong Zhang 
1132ab27746eSHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
1133ab27746eSHong Zhang       il[k] = jmin;
1134ab27746eSHong Zhang       i     = bj[jmin];
1135ab27746eSHong Zhang       jl[k] = jl[i]; jl[i] = k;
1136ab27746eSHong Zhang     }
113749b5e25fSSatish Balay   }
113849b5e25fSSatish Balay 
113949b5e25fSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
1140d8c74875SBarry Smith   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
1141ab27746eSHong Zhang 
11424f79d315SHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_2_NaturalOrdering_inplace;
11434f79d315SHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_2_NaturalOrdering_inplace;
11444f79d315SHong Zhang   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_2_NaturalOrdering_inplace;
11454f79d315SHong Zhang   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_2_NaturalOrdering_inplace;
114649b5e25fSSatish Balay   C->assembled           = PETSC_TRUE;
11475c0bcdfcSHong Zhang   C->preallocated        = PETSC_TRUE;
114826fbe8dcSKarl Rupp 
1149efee365bSSatish Balay   ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
115049b5e25fSSatish Balay   PetscFunctionReturn(0);
115149b5e25fSSatish Balay }
115249b5e25fSSatish Balay 
115349b5e25fSSatish Balay /*
115498a8e78dSHong Zhang     Numeric U^T*D*U factorization for SBAIJ format.
115591602c66SHong Zhang     Version for blocks are 1 by 1.
115649b5e25fSSatish Balay */
1157d595f711SHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_inplace(Mat C,Mat A,const MatFactorInfo *info)
115849b5e25fSSatish Balay {
115949b5e25fSSatish Balay   Mat_SeqSBAIJ   *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ*)C->data;
116049b5e25fSSatish Balay   IS             ip=b->row;
11616849ba73SBarry Smith   PetscErrorCode ierr;
11625d0c19d7SBarry Smith   const PetscInt *ai,*aj,*rip;
11635d0c19d7SBarry Smith   PetscInt       *a2anew,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j,*bcol;
1164997a0949SHong Zhang   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
1165997a0949SHong Zhang   MatScalar      *rtmp,*ba=b->a,*bval,*aa,dk,uikdi;
11660e95ead3SHong Zhang   PetscReal      rs;
11670e95ead3SHong Zhang   FactorShiftCtx sctx;
116849b5e25fSSatish Balay 
116949b5e25fSSatish Balay   PetscFunctionBegin;
11700e95ead3SHong Zhang   /* MatPivotSetUp(): initialize shift context sctx */
11710e95ead3SHong Zhang   ierr = PetscMemzero(&sctx,sizeof(FactorShiftCtx));CHKERRQ(ierr);
11723cea4cbeSHong Zhang 
117349b5e25fSSatish Balay   ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr);
1174cb718733SHong Zhang   if (!a->permute) {
11752d007305SHong Zhang     ai = a->i; aj = a->j; aa = a->a;
11762d007305SHong Zhang   } else {
11772d007305SHong Zhang     ai     = a->inew; aj = a->jnew;
1178fff829cfSHong Zhang     nz     = ai[mbs];
1179785e854fSJed Brown     ierr   = PetscMalloc1(nz,&aa);CHKERRQ(ierr);
1180fff829cfSHong Zhang     a2anew = a->a2anew;
1181997a0949SHong Zhang     bval   = a->a;
1182fff829cfSHong Zhang     for (j=0; j<nz; j++) {
1183997a0949SHong Zhang       aa[a2anew[j]] = *(bval++);
11842d007305SHong Zhang     }
11852d007305SHong Zhang   }
11862d007305SHong Zhang 
118791602c66SHong Zhang   /* initialization */
118849b5e25fSSatish Balay   /* il and jl record the first nonzero element in each row of the accessing
118949b5e25fSSatish Balay      window U(0:k, k:mbs-1).
119049b5e25fSSatish Balay      jl:    list of rows to be added to uneliminated rows
119149b5e25fSSatish Balay             i>= k: jl(i) is the first row to be added to row i
119249b5e25fSSatish Balay             i<  k: jl(i) is the row following row i in some list of rows
119349b5e25fSSatish Balay             jl(i) = mbs indicates the end of a list
119449b5e25fSSatish Balay      il(i): points to the first nonzero element in columns k,...,mbs-1 of
119549b5e25fSSatish Balay             row i of U */
1196dcca6d9dSJed Brown   ierr = PetscMalloc3(mbs,&rtmp,mbs,&il,mbs,&jl);CHKERRQ(ierr);
1197997a0949SHong Zhang 
1198997a0949SHong Zhang   do {
119907b50cabSHong Zhang     sctx.newshift = PETSC_FALSE;
12006df5ee2eSHong Zhang     il[0] = 0;
120149b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
12026df5ee2eSHong Zhang       rtmp[i] = 0.0; jl[i] = mbs;
120349b5e25fSSatish Balay     }
1204997a0949SHong Zhang 
120549b5e25fSSatish Balay     for (k = 0; k<mbs; k++) {
1206997a0949SHong Zhang       /*initialize k-th row by the perm[k]-th row of A */
120749b5e25fSSatish Balay       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
12087701de02SHong Zhang       bval = ba + bi[k];
120949b5e25fSSatish Balay       for (j = jmin; j < jmax; j++) {
1210997a0949SHong Zhang         col       = rip[aj[j]];
1211997a0949SHong Zhang         rtmp[col] = aa[j];
12127701de02SHong Zhang         *bval++   = 0.0; /* for in-place factorization */
121349b5e25fSSatish Balay       }
12143cea4cbeSHong Zhang 
12153cea4cbeSHong Zhang       /* shift the diagonal of the matrix */
12163cea4cbeSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
121749b5e25fSSatish Balay 
121891602c66SHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
121949b5e25fSSatish Balay       dk = rtmp[k];
122049b5e25fSSatish Balay       i  = jl[k]; /* first row to be added to k_th row  */
122149b5e25fSSatish Balay 
1222057f5ba7SHong Zhang       while (i < k) {
122349b5e25fSSatish Balay         nexti = jl[i]; /* next row to be added to k_th row */
122449b5e25fSSatish Balay 
1225fff829cfSHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
122649b5e25fSSatish Balay         ili     = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
1227997a0949SHong Zhang         uikdi   = -ba[ili]*ba[bi[i]]; /* diagonal(k) */
122849b5e25fSSatish Balay         dk     += uikdi*ba[ili];
1229658e7b3eSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
123049b5e25fSSatish Balay 
1231997a0949SHong Zhang         /* add multiple of row i to k-th row */
123249b5e25fSSatish Balay         jmin = ili + 1; jmax = bi[i+1];
123349b5e25fSSatish Balay         if (jmin < jmax) {
123449b5e25fSSatish Balay           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
1235dc0b31edSSatish Balay           ierr = PetscLogFlops(2.0*(jmax-jmin));CHKERRQ(ierr);
1236187a9f4bSHong Zhang 
1237fff829cfSHong Zhang           /* update il and jl for row i */
1238fff829cfSHong Zhang           il[i] = jmin;
1239fff829cfSHong Zhang           j     = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
124049b5e25fSSatish Balay         }
1241ab27746eSHong Zhang         i = nexti;
124249b5e25fSSatish Balay       }
124349b5e25fSSatish Balay 
12443cea4cbeSHong Zhang       /* shift the diagonals when zero pivot is detected */
12453cea4cbeSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
12463cea4cbeSHong Zhang       rs   = 0.0;
1247997a0949SHong Zhang       jmin = bi[k]+1;
1248997a0949SHong Zhang       nz   = bi[k+1] - jmin;
1249997a0949SHong Zhang       if (nz) {
1250997a0949SHong Zhang         bcol = bj + jmin;
1251997a0949SHong Zhang         while (nz--) {
125289f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
125389f845c8SHong Zhang           bcol++;
1254997a0949SHong Zhang         }
1255997a0949SHong Zhang       }
12563cea4cbeSHong Zhang 
12573cea4cbeSHong Zhang       sctx.rs = rs;
12583cea4cbeSHong Zhang       sctx.pv = dk;
12594cccfbddSHong Zhang       ierr    = MatPivotCheck(C,A,info,&sctx,k);CHKERRQ(ierr);
126007b50cabSHong Zhang       if (sctx.newshift) break;    /* sctx.shift_amount is updated */
12610e95ead3SHong Zhang       dk = sctx.pv;
126249b5e25fSSatish Balay 
1263997a0949SHong Zhang       /* copy data into U(k,:) */
126498a8e78dSHong Zhang       ba[bi[k]] = 1.0/dk; /* U(k,k) */
1265fff829cfSHong Zhang       jmin      = bi[k]+1; jmax = bi[k+1];
126649b5e25fSSatish Balay       if (jmin < jmax) {
126749b5e25fSSatish Balay         for (j=jmin; j<jmax; j++) {
1268997a0949SHong Zhang           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
126949b5e25fSSatish Balay         }
127098a8e78dSHong Zhang         /* add the k-th row into il and jl */
127149b5e25fSSatish Balay         il[k] = jmin;
127298a8e78dSHong Zhang         i     = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
127349b5e25fSSatish Balay       }
127449b5e25fSSatish Balay     }
127507b50cabSHong Zhang   } while (sctx.newshift);
1276d8c74875SBarry Smith   ierr = PetscFree3(rtmp,il,jl);CHKERRQ(ierr);
127798a8e78dSHong Zhang   if (a->permute) {ierr = PetscFree(aa);CHKERRQ(ierr);}
127849b5e25fSSatish Balay 
127949b5e25fSSatish Balay   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
128026fbe8dcSKarl Rupp 
12810a3c351aSHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_1_inplace;
12824f79d315SHong Zhang   C->ops->solves         = MatSolves_SeqSBAIJ_1_inplace;
12830a3c351aSHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_1_inplace;
12840a3c351aSHong Zhang   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_1_inplace;
12850a3c351aSHong Zhang   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_1_inplace;
128649b5e25fSSatish Balay   C->assembled           = PETSC_TRUE;
12875c0bcdfcSHong Zhang   C->preallocated        = PETSC_TRUE;
128826fbe8dcSKarl Rupp 
1289d0f46423SBarry Smith   ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr);
12903cea4cbeSHong Zhang   if (sctx.nshift) {
1291f4db908eSBarry Smith     if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
129257622a8eSBarry Smith       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr);
1293f4db908eSBarry Smith     } else if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
129457622a8eSBarry Smith       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr);
1295997a0949SHong Zhang     }
1296997a0949SHong Zhang   }
129749b5e25fSSatish Balay   PetscFunctionReturn(0);
129849b5e25fSSatish Balay }
129949b5e25fSSatish Balay 
1300671cb588SHong Zhang /*
130180d3d5a7SHong Zhang   Version for when blocks are 1 by 1 Using natural ordering under new datastructure
130280d3d5a7SHong Zhang   Modified from MatCholeskyFactorNumeric_SeqAIJ()
1303671cb588SHong Zhang */
1304d595f711SHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info)
1305d595f711SHong Zhang {
1306d595f711SHong Zhang   Mat_SeqSBAIJ   *a=(Mat_SeqSBAIJ*)A->data;
13077b056e98SHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)B->data;
1308d595f711SHong Zhang   PetscErrorCode ierr;
1309d595f711SHong Zhang   PetscInt       i,j,mbs=A->rmap->n,*bi=b->i,*bj=b->j,*bdiag=b->diag,*bjtmp;
1310545dd064SHong Zhang   PetscInt       *ai=a->i,*aj=a->j,*ajtmp;
1311d595f711SHong Zhang   PetscInt       k,jmin,jmax,*c2r,*il,col,nexti,ili,nz;
1312d595f711SHong Zhang   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
1313d90ac83dSHong Zhang   FactorShiftCtx sctx;
1314d595f711SHong Zhang   PetscReal      rs;
1315d595f711SHong Zhang   MatScalar      d,*v;
1316d595f711SHong Zhang 
1317d595f711SHong Zhang   PetscFunctionBegin;
1318dcca6d9dSJed Brown   ierr = PetscMalloc3(mbs,&rtmp,mbs,&il,mbs,&c2r);CHKERRQ(ierr);
1319545dd064SHong Zhang 
1320d595f711SHong Zhang   /* MatPivotSetUp(): initialize shift context sctx */
1321d90ac83dSHong Zhang   ierr = PetscMemzero(&sctx,sizeof(FactorShiftCtx));CHKERRQ(ierr);
1322d595f711SHong Zhang 
1323f4db908eSBarry Smith   if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) { /* set sctx.shift_top=max{rs} */
1324d595f711SHong Zhang     sctx.shift_top = info->zeropivot;
132526fbe8dcSKarl Rupp 
1326580bdb30SBarry Smith     ierr = PetscArrayzero(rtmp,mbs);CHKERRQ(ierr);
132726fbe8dcSKarl Rupp 
1328d595f711SHong Zhang     for (i=0; i<mbs; i++) {
1329d595f711SHong Zhang       /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
1330d595f711SHong Zhang       d        = (aa)[a->diag[i]];
1331545dd064SHong Zhang       rtmp[i] += -PetscRealPart(d);  /* diagonal entry */
1332545dd064SHong Zhang       ajtmp    = aj + ai[i] + 1;     /* exclude diagonal */
1333545dd064SHong Zhang       v        = aa + ai[i] + 1;
1334545dd064SHong Zhang       nz       = ai[i+1] - ai[i] - 1;
1335545dd064SHong Zhang       for (j=0; j<nz; j++) {
1336545dd064SHong Zhang         rtmp[i]        += PetscAbsScalar(v[j]);
1337545dd064SHong Zhang         rtmp[ajtmp[j]] += PetscAbsScalar(v[j]);
1338545dd064SHong Zhang       }
13390838c725SBarry Smith       if (PetscRealPart(rtmp[i]) > sctx.shift_top) sctx.shift_top = PetscRealPart(rtmp[i]);
1340d595f711SHong Zhang     }
1341d595f711SHong Zhang     sctx.shift_top *= 1.1;
1342d595f711SHong Zhang     sctx.nshift_max = 5;
1343d595f711SHong Zhang     sctx.shift_lo   = 0.;
1344d595f711SHong Zhang     sctx.shift_hi   = 1.;
1345d595f711SHong Zhang   }
1346d595f711SHong Zhang 
1347d595f711SHong Zhang   /* allocate working arrays
1348d595f711SHong Zhang      c2r: linked list, keep track of pivot rows for a given column. c2r[col]: head of the list for a given col
1349d595f711SHong Zhang      il:  for active k row, il[i] gives the index of the 1st nonzero entry in U[i,k:n-1] in bj and ba arrays
1350d595f711SHong Zhang   */
1351d595f711SHong Zhang   do {
135207b50cabSHong Zhang     sctx.newshift = PETSC_FALSE;
1353d595f711SHong Zhang 
1354d595f711SHong Zhang     for (i=0; i<mbs; i++) c2r[i] = mbs;
13559e95ef84SSatish Balay     if (mbs) il[0] = 0;
1356d595f711SHong Zhang 
1357d595f711SHong Zhang     for (k = 0; k<mbs; k++) {
1358d595f711SHong Zhang       /* zero rtmp */
1359d595f711SHong Zhang       nz    = bi[k+1] - bi[k];
1360d595f711SHong Zhang       bjtmp = bj + bi[k];
13617b056e98SHong Zhang       for (j=0; j<nz; j++) rtmp[bjtmp[j]] = 0.0;
1362d595f711SHong Zhang 
1363d595f711SHong Zhang       /* load in initial unfactored row */
1364d595f711SHong Zhang       bval = ba + bi[k];
1365d595f711SHong Zhang       jmin = ai[k]; jmax = ai[k+1];
1366d595f711SHong Zhang       for (j = jmin; j < jmax; j++) {
1367d595f711SHong Zhang         col       = aj[j];
1368d595f711SHong Zhang         rtmp[col] = aa[j];
1369d595f711SHong Zhang         *bval++   = 0.0; /* for in-place factorization */
1370d595f711SHong Zhang       }
1371d595f711SHong Zhang       /* shift the diagonal of the matrix: ZeropivotApply() */
1372d595f711SHong Zhang       rtmp[k] += sctx.shift_amount;  /* shift the diagonal of the matrix */
1373d595f711SHong Zhang 
1374d595f711SHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
1375d595f711SHong Zhang       dk = rtmp[k];
1376d595f711SHong Zhang       i  = c2r[k]; /* first row to be added to k_th row  */
1377d595f711SHong Zhang 
1378d595f711SHong Zhang       while (i < k) {
1379d595f711SHong Zhang         nexti = c2r[i]; /* next row to be added to k_th row */
1380d595f711SHong Zhang 
1381d595f711SHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
1382d595f711SHong Zhang         ili     = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
1383d595f711SHong Zhang         uikdi   = -ba[ili]*ba[bdiag[i]]; /* diagonal(k) */
1384d595f711SHong Zhang         dk     += uikdi*ba[ili]; /* update diag[k] */
1385d595f711SHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
1386d595f711SHong Zhang 
1387d595f711SHong Zhang         /* add multiple of row i to k-th row */
1388d595f711SHong Zhang         jmin = ili + 1; jmax = bi[i+1];
1389d595f711SHong Zhang         if (jmin < jmax) {
1390d595f711SHong Zhang           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
1391d595f711SHong Zhang           /* update il and c2r for row i */
1392d595f711SHong Zhang           il[i] = jmin;
1393d595f711SHong Zhang           j     = bj[jmin]; c2r[i] = c2r[j]; c2r[j] = i;
1394d595f711SHong Zhang         }
1395d595f711SHong Zhang         i = nexti;
1396d595f711SHong Zhang       }
1397d595f711SHong Zhang 
1398d595f711SHong Zhang       /* copy data into U(k,:) */
1399d595f711SHong Zhang       rs   = 0.0;
1400d595f711SHong Zhang       jmin = bi[k]; jmax = bi[k+1]-1;
1401d595f711SHong Zhang       if (jmin < jmax) {
1402d595f711SHong Zhang         for (j=jmin; j<jmax; j++) {
1403d595f711SHong Zhang           col = bj[j]; ba[j] = rtmp[col]; rs += PetscAbsScalar(ba[j]);
1404d595f711SHong Zhang         }
1405d595f711SHong Zhang         /* add the k-th row into il and c2r */
1406d595f711SHong Zhang         il[k] = jmin;
1407d595f711SHong Zhang         i     = bj[jmin]; c2r[k] = c2r[i]; c2r[i] = k;
1408d595f711SHong Zhang       }
1409d595f711SHong Zhang 
1410d595f711SHong Zhang       sctx.rs = rs;
1411d595f711SHong Zhang       sctx.pv = dk;
14124cccfbddSHong Zhang       ierr    = MatPivotCheck(B,A,info,&sctx,k);CHKERRQ(ierr);
141307b50cabSHong Zhang       if (sctx.newshift) break;
1414d595f711SHong Zhang       dk = sctx.pv;
1415d595f711SHong Zhang 
1416d595f711SHong Zhang       ba[bdiag[k]] = 1.0/dk; /* U(k,k) */
1417d595f711SHong Zhang     }
141807b50cabSHong Zhang   } while (sctx.newshift);
1419d595f711SHong Zhang 
1420d595f711SHong Zhang   ierr = PetscFree3(rtmp,il,c2r);CHKERRQ(ierr);
1421d595f711SHong Zhang 
1422d595f711SHong Zhang   B->ops->solve          = MatSolve_SeqSBAIJ_1_NaturalOrdering;
14239dcdb97aSHong Zhang   B->ops->solves         = MatSolves_SeqSBAIJ_1;
1424d595f711SHong Zhang   B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1425*910cf402Sprj-   B->ops->matsolve       = MatMatSolve_SeqSBAIJ_1_NaturalOrdering;
142680d3d5a7SHong Zhang   B->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering;
142780d3d5a7SHong Zhang   B->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering;
1428d595f711SHong Zhang 
14297b056e98SHong Zhang   B->assembled    = PETSC_TRUE;
14307b056e98SHong Zhang   B->preallocated = PETSC_TRUE;
143126fbe8dcSKarl Rupp 
14327b056e98SHong Zhang   ierr = PetscLogFlops(B->rmap->n);CHKERRQ(ierr);
1433d595f711SHong Zhang 
1434d595f711SHong Zhang   /* MatPivotView() */
1435d595f711SHong Zhang   if (sctx.nshift) {
1436f4db908eSBarry Smith     if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
143757622a8eSBarry Smith       ierr = PetscInfo4(A,"number of shift_pd tries %D, shift_amount %g, diagonal shifted up by %e fraction top_value %e\n",sctx.nshift,(double)sctx.shift_amount,(double)sctx.shift_fraction,(double)sctx.shift_top);CHKERRQ(ierr);
1438f4db908eSBarry Smith     } else if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
143957622a8eSBarry Smith       ierr = PetscInfo2(A,"number of shift_nz tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr);
1440f4db908eSBarry Smith     } else if (info->shifttype == (PetscReal)MAT_SHIFT_INBLOCKS) {
144157622a8eSBarry Smith       ierr = PetscInfo2(A,"number of shift_inblocks applied %D, each shift_amount %g\n",sctx.nshift,(double)info->shiftamount);CHKERRQ(ierr);
1442d595f711SHong Zhang     }
1443d595f711SHong Zhang   }
1444d595f711SHong Zhang   PetscFunctionReturn(0);
1445d595f711SHong Zhang }
1446d595f711SHong Zhang 
1447d595f711SHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering_inplace(Mat C,Mat A,const MatFactorInfo *info)
1448671cb588SHong Zhang {
1449671cb588SHong Zhang   Mat_SeqSBAIJ   *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ*)C->data;
1450dfbe8321SBarry Smith   PetscErrorCode ierr;
145113f74950SBarry Smith   PetscInt       i,j,mbs = a->mbs;
145213f74950SBarry Smith   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
14533cea4cbeSHong Zhang   PetscInt       k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz;
1454653a6975SHong Zhang   MatScalar      *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval;
145514d2772eSHong Zhang   PetscReal      rs;
14560e95ead3SHong Zhang   FactorShiftCtx sctx;
1457653a6975SHong Zhang 
1458653a6975SHong Zhang   PetscFunctionBegin;
14590e95ead3SHong Zhang   /* MatPivotSetUp(): initialize shift context sctx */
14600e95ead3SHong Zhang   ierr = PetscMemzero(&sctx,sizeof(FactorShiftCtx));CHKERRQ(ierr);
14610e95ead3SHong Zhang 
1462653a6975SHong Zhang   /* initialization */
1463653a6975SHong Zhang   /* il and jl record the first nonzero element in each row of the accessing
1464653a6975SHong Zhang      window U(0:k, k:mbs-1).
1465653a6975SHong Zhang      jl:    list of rows to be added to uneliminated rows
1466653a6975SHong Zhang             i>= k: jl(i) is the first row to be added to row i
1467653a6975SHong Zhang             i<  k: jl(i) is the row following row i in some list of rows
1468653a6975SHong Zhang             jl(i) = mbs indicates the end of a list
1469653a6975SHong Zhang      il(i): points to the first nonzero element in U(i,k:mbs-1)
1470653a6975SHong Zhang   */
1471785e854fSJed Brown   ierr = PetscMalloc1(mbs,&rtmp);CHKERRQ(ierr);
1472dcca6d9dSJed Brown   ierr = PetscMalloc2(mbs,&il,mbs,&jl);CHKERRQ(ierr);
1473b00f7748SHong Zhang 
1474b00f7748SHong Zhang   do {
147507b50cabSHong Zhang     sctx.newshift = PETSC_FALSE;
14766df5ee2eSHong Zhang     il[0] = 0;
1477653a6975SHong Zhang     for (i=0; i<mbs; i++) {
14786df5ee2eSHong Zhang       rtmp[i] = 0.0; jl[i] = mbs;
1479653a6975SHong Zhang     }
1480653a6975SHong Zhang 
1481997a0949SHong Zhang     for (k = 0; k<mbs; k++) {
1482653a6975SHong Zhang       /*initialize k-th row with elements nonzero in row perm(k) of A */
1483653a6975SHong Zhang       nz   = ai[k+1] - ai[k];
1484653a6975SHong Zhang       acol = aj + ai[k];
1485653a6975SHong Zhang       aval = aa + ai[k];
1486653a6975SHong Zhang       bval = ba + bi[k];
1487653a6975SHong Zhang       while (nz--) {
1488653a6975SHong Zhang         rtmp[*acol++] = *aval++;
1489653a6975SHong Zhang         *bval++       = 0.0; /* for in-place factorization */
1490653a6975SHong Zhang       }
14913cea4cbeSHong Zhang 
14923cea4cbeSHong Zhang       /* shift the diagonal of the matrix */
14933cea4cbeSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
1494653a6975SHong Zhang 
1495653a6975SHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
1496653a6975SHong Zhang       dk = rtmp[k];
1497653a6975SHong Zhang       i  = jl[k]; /* first row to be added to k_th row  */
1498653a6975SHong Zhang 
1499653a6975SHong Zhang       while (i < k) {
1500653a6975SHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
1501653a6975SHong Zhang         /* compute multiplier, update D(k) and U(i,k) */
1502653a6975SHong Zhang         ili     = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
1503653a6975SHong Zhang         uikdi   = -ba[ili]*ba[bi[i]];
1504653a6975SHong Zhang         dk     += uikdi*ba[ili];
1505653a6975SHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
1506653a6975SHong Zhang 
1507653a6975SHong Zhang         /* add multiple of row i to k-th row ... */
1508653a6975SHong Zhang         jmin = ili + 1;
1509653a6975SHong Zhang         nz   = bi[i+1] - jmin;
1510653a6975SHong Zhang         if (nz > 0) {
1511653a6975SHong Zhang           bcol = bj + jmin;
1512653a6975SHong Zhang           bval = ba + jmin;
1513dc0b31edSSatish Balay           ierr = PetscLogFlops(2.0*nz);CHKERRQ(ierr);
1514653a6975SHong Zhang           while (nz--) rtmp[*bcol++] += uikdi*(*bval++);
1515187a9f4bSHong Zhang 
1516997a0949SHong Zhang           /* update il and jl for i-th row */
1517997a0949SHong Zhang           il[i] = jmin;
1518997a0949SHong Zhang           j     = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
1519653a6975SHong Zhang         }
1520653a6975SHong Zhang         i = nexti;
1521653a6975SHong Zhang       }
1522653a6975SHong Zhang 
15233cea4cbeSHong Zhang       /* shift the diagonals when zero pivot is detected */
15243cea4cbeSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
15253cea4cbeSHong Zhang       rs   = 0.0;
152621c26570Svictorle       jmin = bi[k]+1;
152721c26570Svictorle       nz   = bi[k+1] - jmin;
152821c26570Svictorle       if (nz) {
152921c26570Svictorle         bcol = bj + jmin;
153021c26570Svictorle         while (nz--) {
153189f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
153289f845c8SHong Zhang           bcol++;
153321c26570Svictorle         }
153421c26570Svictorle       }
15353cea4cbeSHong Zhang 
15363cea4cbeSHong Zhang       sctx.rs = rs;
15373cea4cbeSHong Zhang       sctx.pv = dk;
15384cccfbddSHong Zhang       ierr    = MatPivotCheck(C,A,info,&sctx,k);CHKERRQ(ierr);
153907b50cabSHong Zhang       if (sctx.newshift) break;    /* sctx.shift_amount is updated */
15400e95ead3SHong Zhang       dk = sctx.pv;
1541653a6975SHong Zhang 
1542997a0949SHong Zhang       /* copy data into U(k,:) */
1543653a6975SHong Zhang       ba[bi[k]] = 1.0/dk;
1544653a6975SHong Zhang       jmin      = bi[k]+1;
1545653a6975SHong Zhang       nz        = bi[k+1] - jmin;
1546653a6975SHong Zhang       if (nz) {
1547653a6975SHong Zhang         bcol = bj + jmin;
1548653a6975SHong Zhang         bval = ba + jmin;
1549653a6975SHong Zhang         while (nz--) {
1550653a6975SHong Zhang           *bval++       = rtmp[*bcol];
1551653a6975SHong Zhang           rtmp[*bcol++] = 0.0;
1552653a6975SHong Zhang         }
1553997a0949SHong Zhang         /* add k-th row into il and jl */
1554653a6975SHong Zhang         il[k] = jmin;
1555997a0949SHong Zhang         i     = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
1556653a6975SHong Zhang       }
1557b00f7748SHong Zhang     } /* end of for (k = 0; k<mbs; k++) */
155807b50cabSHong Zhang   } while (sctx.newshift);
1559653a6975SHong Zhang   ierr = PetscFree(rtmp);CHKERRQ(ierr);
1560d8c74875SBarry Smith   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
1561653a6975SHong Zhang 
15620a3c351aSHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
15634f79d315SHong Zhang   C->ops->solves         = MatSolves_SeqSBAIJ_1_inplace;
15640a3c351aSHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
15650a3c351aSHong Zhang   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
15660a3c351aSHong Zhang   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
1567db4efbfdSBarry Smith 
1568653a6975SHong Zhang   C->assembled    = PETSC_TRUE;
1569653a6975SHong Zhang   C->preallocated = PETSC_TRUE;
157026fbe8dcSKarl Rupp 
1571d0f46423SBarry Smith   ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr);
15723cea4cbeSHong Zhang   if (sctx.nshift) {
1573f4db908eSBarry Smith     if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
157457622a8eSBarry Smith       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr);
1575f4db908eSBarry Smith     } else if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
157657622a8eSBarry Smith       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr);
1577b00f7748SHong Zhang     }
1578fb10cecfSBarry Smith   }
1579653a6975SHong Zhang   PetscFunctionReturn(0);
1580653a6975SHong Zhang }
1581653a6975SHong Zhang 
15820481f469SBarry Smith PetscErrorCode MatCholeskyFactor_SeqSBAIJ(Mat A,IS perm,const MatFactorInfo *info)
158349b5e25fSSatish Balay {
1584dfbe8321SBarry Smith   PetscErrorCode ierr;
158549b5e25fSSatish Balay   Mat            C;
158649b5e25fSSatish Balay 
158749b5e25fSSatish Balay   PetscFunctionBegin;
1588719d5645SBarry Smith   ierr = MatGetFactor(A,"petsc",MAT_FACTOR_CHOLESKY,&C);CHKERRQ(ierr);
1589719d5645SBarry Smith   ierr = MatCholeskyFactorSymbolic(C,A,perm,info);CHKERRQ(ierr);
1590719d5645SBarry Smith   ierr = MatCholeskyFactorNumeric(C,A,info);CHKERRQ(ierr);
159126fbe8dcSKarl Rupp 
1592db4efbfdSBarry Smith   A->ops->solve          = C->ops->solve;
1593db4efbfdSBarry Smith   A->ops->solvetranspose = C->ops->solvetranspose;
159426fbe8dcSKarl Rupp 
159528be2f97SBarry Smith   ierr = MatHeaderMerge(A,&C);CHKERRQ(ierr);
159649b5e25fSSatish Balay   PetscFunctionReturn(0);
159749b5e25fSSatish Balay }
1598