xref: /petsc/src/mat/impls/sbaij/seq/sbaijfact.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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;
14*5f80ce2aSJacob Faibussowitsch   PetscCheck(bs == 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for bs: %" PetscInt_FMT " >1 yet",bs);
154ff4e45dSHong Zhang 
164ff4e45dSHong Zhang   nneg_tmp = 0; npos_tmp = 0;
17eeeff2ecSHong Zhang   for (i=0; i<mbs; i++) {
1826fbe8dcSKarl Rupp     if (PetscRealPart(dd[*fi]) > 0.0) npos_tmp++;
194ff4e45dSHong Zhang     else if (PetscRealPart(dd[*fi]) < 0.0) nneg_tmp++;
20eeeff2ecSHong Zhang     fi++;
213e0d88b5SBarry Smith   }
224ff4e45dSHong Zhang   if (nneg)  *nneg  = nneg_tmp;
23eeeff2ecSHong Zhang   if (npos)  *npos  = npos_tmp;
244ff4e45dSHong Zhang   if (nzero) *nzero = mbs - nneg_tmp - npos_tmp;
255f9f512dSHong Zhang   PetscFunctionReturn(0);
265f9f512dSHong Zhang }
275f9f512dSHong Zhang 
285f9f512dSHong Zhang /*
295f9f512dSHong Zhang   Symbolic U^T*D*U factorization for SBAIJ format. Modified from SSF of YSMP.
3010c27e3dSHong Zhang   Use Modified Sparse Row (MSR) storage for u and ju. See page 85, "Iterative Methods ..." by Saad.
315f9f512dSHong Zhang */
320481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ_MSR(Mat F,Mat A,IS perm,const MatFactorInfo *info)
3310c27e3dSHong Zhang {
3410c27e3dSHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b;
355d0c19d7SBarry Smith   const PetscInt *rip,*ai,*aj;
365d0c19d7SBarry Smith   PetscInt       i,mbs = a->mbs,*jutmp,bs = A->rmap->bs,bs2=a->bs2;
3710c27e3dSHong Zhang   PetscInt       m,reallocs = 0,prow;
3810c27e3dSHong Zhang   PetscInt       *jl,*q,jmin,jmax,juidx,nzk,qm,*iu,*ju,k,j,vj,umax,maxadd;
3910c27e3dSHong Zhang   PetscReal      f = info->fill;
40ace3abfcSBarry Smith   PetscBool      perm_identity;
4110c27e3dSHong Zhang 
4210c27e3dSHong Zhang   PetscFunctionBegin;
4310c27e3dSHong Zhang   /* check whether perm is the identity mapping */
44*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISIdentity(perm,&perm_identity));
45*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetIndices(perm,&rip));
4610c27e3dSHong Zhang 
4710c27e3dSHong Zhang   if (perm_identity) { /* without permutation */
4810c27e3dSHong Zhang     a->permute = PETSC_FALSE;
4926fbe8dcSKarl Rupp 
5010c27e3dSHong Zhang     ai = a->i; aj = a->j;
5110c27e3dSHong Zhang   } else {            /* non-trivial permutation */
5210c27e3dSHong Zhang     a->permute = PETSC_TRUE;
5326fbe8dcSKarl Rupp 
54*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatReorderingSeqSBAIJ(A,perm));
5526fbe8dcSKarl Rupp 
5610c27e3dSHong Zhang     ai = a->inew; aj = a->jnew;
5710c27e3dSHong Zhang   }
5810c27e3dSHong Zhang 
5910c27e3dSHong Zhang   /* initialization */
60*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(mbs+1,&iu));
6110c27e3dSHong Zhang   umax  = (PetscInt)(f*ai[mbs] + 1); umax += mbs + 1;
62*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(umax,&ju));
6310c27e3dSHong Zhang   iu[0] = mbs+1;
6410c27e3dSHong Zhang   juidx = mbs + 1; /* index for ju */
65d8c74875SBarry Smith   /* jl linked list for pivot row -- linked list for col index */
66*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc2(mbs,&jl,mbs,&q));
6710c27e3dSHong Zhang   for (i=0; i<mbs; i++) {
6810c27e3dSHong Zhang     jl[i] = mbs;
6910c27e3dSHong Zhang     q[i]  = 0;
7010c27e3dSHong Zhang   }
7110c27e3dSHong Zhang 
7210c27e3dSHong Zhang   /* for each row k */
7310c27e3dSHong Zhang   for (k=0; k<mbs; k++) {
7410c27e3dSHong Zhang     for (i=0; i<mbs; i++) q[i] = 0;  /* to be removed! */
7510c27e3dSHong Zhang     nzk  = 0; /* num. of nz blocks in k-th block row with diagonal block excluded */
7610c27e3dSHong Zhang     q[k] = mbs;
7710c27e3dSHong Zhang     /* initialize nonzero structure of k-th row to row rip[k] of A */
7810c27e3dSHong Zhang     jmin = ai[rip[k]] +1; /* exclude diag[k] */
7910c27e3dSHong Zhang     jmax = ai[rip[k]+1];
8010c27e3dSHong Zhang     for (j=jmin; j<jmax; j++) {
8110c27e3dSHong Zhang       vj = rip[aj[j]]; /* col. value */
8210c27e3dSHong Zhang       if (vj > k) {
8310c27e3dSHong Zhang         qm = k;
8410c27e3dSHong Zhang         do {
8510c27e3dSHong Zhang           m = qm; qm = q[m];
8610c27e3dSHong Zhang         } while (qm < vj);
87*5f80ce2aSJacob Faibussowitsch         PetscCheck(qm != vj,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Duplicate entry in A");
8810c27e3dSHong Zhang         nzk++;
8910c27e3dSHong Zhang         q[m]  = vj;
9010c27e3dSHong Zhang         q[vj] = qm;
9110c27e3dSHong Zhang       } /* if (vj > k) */
9210c27e3dSHong Zhang     } /* for (j=jmin; j<jmax; j++) */
9310c27e3dSHong Zhang 
9410c27e3dSHong Zhang     /* modify nonzero structure of k-th row by computing fill-in
9510c27e3dSHong Zhang        for each row i to be merged in */
9610c27e3dSHong Zhang     prow = k;
9710c27e3dSHong Zhang     prow = jl[prow]; /* next pivot row (== mbs for symbolic factorization) */
9810c27e3dSHong Zhang 
9910c27e3dSHong Zhang     while (prow < k) {
10010c27e3dSHong Zhang       /* merge row prow into k-th row */
10110c27e3dSHong Zhang       jmin = iu[prow] + 1; jmax = iu[prow+1];
10210c27e3dSHong Zhang       qm   = k;
10310c27e3dSHong Zhang       for (j=jmin; j<jmax; j++) {
10410c27e3dSHong Zhang         vj = ju[j];
10510c27e3dSHong Zhang         do {
10610c27e3dSHong Zhang           m = qm; qm = q[m];
10710c27e3dSHong Zhang         } while (qm < vj);
10810c27e3dSHong Zhang         if (qm != vj) {
10910c27e3dSHong Zhang           nzk++; q[m] = vj; q[vj] = qm; qm = vj;
11010c27e3dSHong Zhang         }
11110c27e3dSHong Zhang       }
11210c27e3dSHong Zhang       prow = jl[prow]; /* next pivot row */
11310c27e3dSHong Zhang     }
11410c27e3dSHong Zhang 
11510c27e3dSHong Zhang     /* add k to row list for first nonzero element in k-th row */
11610c27e3dSHong Zhang     if (nzk > 0) {
11710c27e3dSHong Zhang       i     = q[k]; /* col value of first nonzero element in U(k, k+1:mbs-1) */
11810c27e3dSHong Zhang       jl[k] = jl[i]; jl[i] = k;
11910c27e3dSHong Zhang     }
12010c27e3dSHong Zhang     iu[k+1] = iu[k] + nzk;
12110c27e3dSHong Zhang 
12210c27e3dSHong Zhang     /* allocate more space to ju if needed */
12310c27e3dSHong Zhang     if (iu[k+1] > umax) {
12410c27e3dSHong Zhang       /* estimate how much additional space we will need */
12510c27e3dSHong Zhang       /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */
12610c27e3dSHong Zhang       /* just double the memory each time */
12710c27e3dSHong Zhang       maxadd = umax;
12810c27e3dSHong Zhang       if (maxadd < nzk) maxadd = (mbs-k)*(nzk+1)/2;
12910c27e3dSHong Zhang       umax += maxadd;
13010c27e3dSHong Zhang 
13110c27e3dSHong Zhang       /* allocate a longer ju */
132*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscMalloc1(umax,&jutmp));
133*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscArraycpy(jutmp,ju,iu[k]));
134*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFree(ju));
13510c27e3dSHong Zhang       ju   = jutmp;
13610c27e3dSHong Zhang       reallocs++; /* count how many times we realloc */
13710c27e3dSHong Zhang     }
13810c27e3dSHong Zhang 
13910c27e3dSHong Zhang     /* save nonzero structure of k-th row in ju */
14010c27e3dSHong Zhang     i=k;
14110c27e3dSHong Zhang     while (nzk--) {
14210c27e3dSHong Zhang       i           = q[i];
14310c27e3dSHong Zhang       ju[juidx++] = i;
14410c27e3dSHong Zhang     }
14510c27e3dSHong Zhang   }
14610c27e3dSHong Zhang 
1476cf91177SBarry Smith #if defined(PETSC_USE_INFO)
14810c27e3dSHong Zhang   if (ai[mbs] != 0) {
14910c27e3dSHong Zhang     PetscReal af = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]);
150*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscInfo(A,"Reallocs %" PetscInt_FMT " Fill ratio:given %g needed %g\n",reallocs,(double)f,(double)af));
151*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscInfo(A,"Run with -pc_factor_fill %g or use \n",(double)af));
152*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscInfo(A,"PCFactorSetFill(pc,%g);\n",(double)af));
153*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscInfo(A,"for best performance.\n"));
15410c27e3dSHong Zhang   } else {
155*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscInfo(A,"Empty matrix\n"));
15610c27e3dSHong Zhang   }
15763ba0a88SBarry Smith #endif
15810c27e3dSHong Zhang 
159*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISRestoreIndices(perm,&rip));
160*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree2(jl,q));
16110c27e3dSHong Zhang 
16210c27e3dSHong Zhang   /* put together the new matrix */
163*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqSBAIJSetPreallocation(F,bs,MAT_SKIP_ALLOCATION,NULL));
16410c27e3dSHong Zhang 
165*5f80ce2aSJacob Faibussowitsch   /* CHKERRQ(PetscLogObjectParent((PetscObject)B,(PetscObject)iperm)); */
166719d5645SBarry Smith   b                = (Mat_SeqSBAIJ*)(F)->data;
16710c27e3dSHong Zhang   b->singlemalloc  = PETSC_FALSE;
168e6b907acSBarry Smith   b->free_a        = PETSC_TRUE;
169e6b907acSBarry Smith   b->free_ij       = PETSC_TRUE;
17026fbe8dcSKarl Rupp 
171*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1((iu[mbs]+1)*bs2,&b->a));
17210c27e3dSHong Zhang   b->j    = ju;
17310c27e3dSHong Zhang   b->i    = iu;
174f4259b30SLisandro Dalcin   b->diag = NULL;
175f4259b30SLisandro Dalcin   b->ilen = NULL;
176f4259b30SLisandro Dalcin   b->imax = NULL;
17710c27e3dSHong Zhang   b->row  = perm;
17826fbe8dcSKarl Rupp 
17910c27e3dSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
18026fbe8dcSKarl Rupp 
181*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectReference((PetscObject)perm));
18226fbe8dcSKarl Rupp 
18310c27e3dSHong Zhang   b->icol = perm;
184*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectReference((PetscObject)perm));
185*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(bs*mbs+bs,&b->solve_work));
18610c27e3dSHong Zhang   /* In b structure:  Free imax, ilen, old a, old j.
18710c27e3dSHong Zhang      Allocate idnew, solve_work, new a, new j */
188*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogObjectMemory((PetscObject)F,(iu[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar))));
18910c27e3dSHong Zhang   b->maxnz = b->nz = iu[mbs];
19010c27e3dSHong Zhang 
191719d5645SBarry Smith   (F)->info.factor_mallocs   = reallocs;
192719d5645SBarry Smith   (F)->info.fill_ratio_given = f;
19310c27e3dSHong Zhang   if (ai[mbs] != 0) {
194719d5645SBarry Smith     (F)->info.fill_ratio_needed = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]);
19510c27e3dSHong Zhang   } else {
196719d5645SBarry Smith     (F)->info.fill_ratio_needed = 0.0;
19710c27e3dSHong Zhang   }
198*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqSBAIJSetNumericFactorization_inplace(F,perm_identity));
19910c27e3dSHong Zhang   PetscFunctionReturn(0);
20010c27e3dSHong Zhang }
20110c27e3dSHong Zhang /*
20210c27e3dSHong Zhang     Symbolic U^T*D*U factorization for SBAIJ format.
20380d3d5a7SHong Zhang     See MatICCFactorSymbolic_SeqAIJ() for description of its data structure.
20410c27e3dSHong Zhang */
205c6db04a5SJed Brown #include <petscbt.h>
206c6db04a5SJed Brown #include <../src/mat/utils/freespace.h>
2070481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
20849b5e25fSSatish Balay {
20998a8e78dSHong Zhang   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data;
21098a8e78dSHong Zhang   Mat_SeqSBAIJ       *b;
211ace3abfcSBarry Smith   PetscBool          perm_identity,missing;
21298a8e78dSHong Zhang   PetscReal          fill = info->fill;
2137b056e98SHong Zhang   const PetscInt     *rip,*ai=a->i,*aj=a->j;
2149186b105SHong Zhang   PetscInt           i,mbs=a->mbs,bs=A->rmap->bs,reallocs=0,prow;
21598a8e78dSHong Zhang   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
216c6d0d4f0SHong Zhang   PetscInt           nlnk,*lnk,ncols,*cols,*uj,**ui_ptr,*uj_ptr,*udiag;
2170298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
21898a8e78dSHong Zhang   PetscBT            lnkbt;
219d595f711SHong Zhang 
220d595f711SHong Zhang   PetscFunctionBegin;
221*5f80ce2aSJacob Faibussowitsch   PetscCheck(A->rmap->n == A->cmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Must be square matrix, rows %" PetscInt_FMT " columns %" PetscInt_FMT,A->rmap->n,A->cmap->n);
222*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMissingDiagonal(A,&missing,&i));
223*5f80ce2aSJacob Faibussowitsch   PetscCheck(!missing,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %" PetscInt_FMT,i);
2246d07c2b1SHong Zhang   if (bs > 1) {
225*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCholeskyFactorSymbolic_SeqSBAIJ_inplace(fact,A,perm,info));
2266d07c2b1SHong Zhang     PetscFunctionReturn(0);
2276d07c2b1SHong Zhang   }
228d595f711SHong Zhang 
229d595f711SHong Zhang   /* check whether perm is the identity mapping */
230*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISIdentity(perm,&perm_identity));
231*5f80ce2aSJacob Faibussowitsch   PetscCheck(perm_identity,PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix reordering is not supported for sbaij matrix. Use aij format");
232d595f711SHong Zhang   a->permute = PETSC_FALSE;
233*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetIndices(perm,&rip));
234d595f711SHong Zhang 
235d595f711SHong Zhang   /* initialization */
236*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(mbs+1,&ui));
237*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(mbs+1,&udiag));
238d595f711SHong Zhang   ui[0] = 0;
239d595f711SHong Zhang 
240d595f711SHong Zhang   /* jl: linked list for storing indices of the pivot rows
241d595f711SHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
242*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc4(mbs,&ui_ptr,mbs,&il,mbs,&jl,mbs,&cols));
243d595f711SHong Zhang   for (i=0; i<mbs; i++) {
244d595f711SHong Zhang     jl[i] = mbs; il[i] = 0;
245d595f711SHong Zhang   }
246d595f711SHong Zhang 
247d595f711SHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
248d595f711SHong Zhang   nlnk = mbs + 1;
249*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt));
250d595f711SHong Zhang 
251d595f711SHong Zhang   /* initial FreeSpace size is fill*(ai[mbs]+1) */
252*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,ai[mbs]+1),&free_space));
253d595f711SHong Zhang   current_space = free_space;
254d595f711SHong Zhang 
255d595f711SHong Zhang   for (k=0; k<mbs; k++) {  /* for each active row k */
256d595f711SHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
257d595f711SHong Zhang     nzk   = 0;
258c6d0d4f0SHong Zhang     ncols = ai[k+1] - ai[k];
259*5f80ce2aSJacob Faibussowitsch     PetscCheck(ncols,PETSC_COMM_SELF,PETSC_ERR_MAT_CH_ZRPVT,"Empty row %" PetscInt_FMT " in matrix ",k);
260d595f711SHong Zhang     for (j=0; j<ncols; j++) {
261c6d0d4f0SHong Zhang       i       = *(aj + ai[k] + j);
262c6d0d4f0SHong Zhang       cols[j] = i;
263d595f711SHong Zhang     }
264*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscLLAdd(ncols,cols,mbs,nlnk,lnk,lnkbt));
265d595f711SHong Zhang     nzk += nlnk;
266d595f711SHong Zhang 
267d595f711SHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
268d595f711SHong Zhang     prow = jl[k]; /* 1st pivot row */
269d595f711SHong Zhang 
270d595f711SHong Zhang     while (prow < k) {
271d595f711SHong Zhang       nextprow = jl[prow];
272d595f711SHong Zhang       /* merge prow into k-th row */
273d595f711SHong Zhang       jmin   = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
274d595f711SHong Zhang       jmax   = ui[prow+1];
275d595f711SHong Zhang       ncols  = jmax-jmin;
276d595f711SHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
277*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt));
278d595f711SHong Zhang       nzk   += nlnk;
279d595f711SHong Zhang 
280d595f711SHong Zhang       /* update il and jl for prow */
281d595f711SHong Zhang       if (jmin < jmax) {
282d595f711SHong Zhang         il[prow] = jmin;
283d595f711SHong Zhang         j        = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
284d595f711SHong Zhang       }
285d595f711SHong Zhang       prow = nextprow;
286d595f711SHong Zhang     }
287d595f711SHong Zhang 
288d595f711SHong Zhang     /* if free space is not available, make more free space */
289d595f711SHong Zhang     if (current_space->local_remaining<nzk) {
290d595f711SHong Zhang       i    = mbs - k + 1; /* num of unfactored rows */
291f91af8c7SBarry Smith       i    = PetscIntMultTruncate(i,PetscMin(nzk, i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
292*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFreeSpaceGet(i,&current_space));
293d595f711SHong Zhang       reallocs++;
294d595f711SHong Zhang     }
295d595f711SHong Zhang 
296d595f711SHong Zhang     /* copy data into free space, then initialize lnk */
297*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt));
298d595f711SHong Zhang 
299d595f711SHong Zhang     /* add the k-th row into il and jl */
3007b056e98SHong Zhang     if (nzk > 1) {
301d595f711SHong Zhang       i     = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
302d595f711SHong Zhang       jl[k] = jl[i]; jl[i] = k;
303d595f711SHong Zhang       il[k] = ui[k] + 1;
304d595f711SHong Zhang     }
305d595f711SHong Zhang     ui_ptr[k] = current_space->array;
30626fbe8dcSKarl Rupp 
307d595f711SHong Zhang     current_space->array           += nzk;
308d595f711SHong Zhang     current_space->local_used      += nzk;
309d595f711SHong Zhang     current_space->local_remaining -= nzk;
310d595f711SHong Zhang 
311d595f711SHong Zhang     ui[k+1] = ui[k] + nzk;
312d595f711SHong Zhang   }
313d595f711SHong Zhang 
314*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISRestoreIndices(perm,&rip));
315*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree4(ui_ptr,il,jl,cols));
316d595f711SHong Zhang 
317d595f711SHong Zhang   /* destroy list of free space and other temporary array(s) */
318*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(ui[mbs]+1,&uj));
319*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFreeSpaceContiguous_Cholesky(&free_space,uj,mbs,ui,udiag)); /* store matrix factor */
320*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLLDestroy(lnk,lnkbt));
321d595f711SHong Zhang 
322d595f711SHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
323*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqSBAIJSetPreallocation(fact,bs,MAT_SKIP_ALLOCATION,NULL));
324d595f711SHong Zhang 
3257b056e98SHong Zhang   b               = (Mat_SeqSBAIJ*)fact->data;
326d595f711SHong Zhang   b->singlemalloc = PETSC_FALSE;
327d595f711SHong Zhang   b->free_a       = PETSC_TRUE;
328d595f711SHong Zhang   b->free_ij      = PETSC_TRUE;
32926fbe8dcSKarl Rupp 
330*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(ui[mbs]+1,&b->a));
33126fbe8dcSKarl Rupp 
332d595f711SHong Zhang   b->j         = uj;
333d595f711SHong Zhang   b->i         = ui;
334c6d0d4f0SHong Zhang   b->diag      = udiag;
335c6d0d4f0SHong Zhang   b->free_diag = PETSC_TRUE;
336f4259b30SLisandro Dalcin   b->ilen      = NULL;
337f4259b30SLisandro Dalcin   b->imax      = NULL;
338d595f711SHong Zhang   b->row       = perm;
339d595f711SHong Zhang   b->icol      = perm;
34026fbe8dcSKarl Rupp 
341*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectReference((PetscObject)perm));
342*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectReference((PetscObject)perm));
34326fbe8dcSKarl Rupp 
3447b056e98SHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
34526fbe8dcSKarl Rupp 
346*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(mbs+1,&b->solve_work));
347*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogObjectMemory((PetscObject)fact,ui[mbs]*(sizeof(PetscInt)+sizeof(MatScalar))));
34826fbe8dcSKarl Rupp 
349d595f711SHong Zhang   b->maxnz = b->nz = ui[mbs];
350d595f711SHong Zhang 
35195b5ac22SHong Zhang   fact->info.factor_mallocs   = reallocs;
35295b5ac22SHong Zhang   fact->info.fill_ratio_given = fill;
353d595f711SHong Zhang   if (ai[mbs] != 0) {
35495b5ac22SHong Zhang     fact->info.fill_ratio_needed = ((PetscReal)ui[mbs])/ai[mbs];
355d595f711SHong Zhang   } else {
35695b5ac22SHong Zhang     fact->info.fill_ratio_needed = 0.0;
357d595f711SHong Zhang   }
35895b5ac22SHong Zhang #if defined(PETSC_USE_INFO)
35995b5ac22SHong Zhang   if (ai[mbs] != 0) {
36095b5ac22SHong Zhang     PetscReal af = fact->info.fill_ratio_needed;
361*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscInfo(A,"Reallocs %" PetscInt_FMT " Fill ratio:given %g needed %g\n",reallocs,(double)fill,(double)af));
362*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscInfo(A,"Run with -pc_factor_fill %g or use \n",(double)af));
363*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscInfo(A,"PCFactorSetFill(pc,%g) for best performance.\n",(double)af));
36495b5ac22SHong Zhang   } else {
365*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscInfo(A,"Empty matrix\n"));
36695b5ac22SHong Zhang   }
36795b5ac22SHong Zhang #endif
368c6d0d4f0SHong Zhang   fact->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering;
369d595f711SHong Zhang   PetscFunctionReturn(0);
370d595f711SHong Zhang }
371d595f711SHong Zhang 
372d595f711SHong Zhang PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ_inplace(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
373d595f711SHong Zhang {
374d595f711SHong Zhang   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data;
375d595f711SHong Zhang   Mat_SeqSBAIJ       *b;
376ace3abfcSBarry Smith   PetscBool          perm_identity,missing;
377d595f711SHong Zhang   PetscReal          fill = info->fill;
378d595f711SHong Zhang   const PetscInt     *rip,*ai,*aj;
379d595f711SHong Zhang   PetscInt           i,mbs=a->mbs,bs=A->rmap->bs,reallocs=0,prow,d;
380d595f711SHong Zhang   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
381d595f711SHong Zhang   PetscInt           nlnk,*lnk,ncols,*cols,*uj,**ui_ptr,*uj_ptr;
3820298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
383d595f711SHong Zhang   PetscBT            lnkbt;
38449b5e25fSSatish Balay 
38549b5e25fSSatish Balay   PetscFunctionBegin;
386*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMissingDiagonal(A,&missing,&d));
387*5f80ce2aSJacob Faibussowitsch   PetscCheck(!missing,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %" PetscInt_FMT,d);
38858ebbce7SBarry Smith 
38910c27e3dSHong Zhang   /*
39010c27e3dSHong Zhang    This code originally uses Modified Sparse Row (MSR) storage
39110c27e3dSHong Zhang    (see page 85, "Iterative Methods ..." by Saad) for the output matrix B - bad choise!
39210c27e3dSHong Zhang    Then it is rewritten so the factor B takes seqsbaij format. However the associated
39310c27e3dSHong Zhang    MatCholeskyFactorNumeric_() have not been modified for the cases of bs>1 or !perm_identity,
39410c27e3dSHong Zhang    thus the original code in MSR format is still used for these cases.
39510c27e3dSHong Zhang    The code below should replace MatCholeskyFactorSymbolic_SeqSBAIJ_MSR() whenever
39610c27e3dSHong Zhang    MatCholeskyFactorNumeric_() is modified for using sbaij symbolic factor.
39710c27e3dSHong Zhang   */
398fff829cfSHong Zhang   if (bs > 1) {
399*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCholeskyFactorSymbolic_SeqSBAIJ_MSR(fact,A,perm,info));
400fff829cfSHong Zhang     PetscFunctionReturn(0);
401fff829cfSHong Zhang   }
40210c27e3dSHong Zhang 
40398a8e78dSHong Zhang   /* check whether perm is the identity mapping */
404*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISIdentity(perm,&perm_identity));
405*5f80ce2aSJacob Faibussowitsch   PetscCheck(perm_identity,PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix reordering is not supported for sbaij matrix. Use aij format");
406fff829cfSHong Zhang   a->permute = PETSC_FALSE;
407fff829cfSHong Zhang   ai = a->i; aj = a->j;
408*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetIndices(perm,&rip));
409fff829cfSHong Zhang 
410fff829cfSHong Zhang   /* initialization */
411*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(mbs+1,&ui));
41298a8e78dSHong Zhang   ui[0] = 0;
41398a8e78dSHong Zhang 
41498a8e78dSHong Zhang   /* jl: linked list for storing indices of the pivot rows
4157625dc9aSHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
416*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc4(mbs,&ui_ptr,mbs,&il,mbs,&jl,mbs,&cols));
4177625dc9aSHong Zhang   for (i=0; i<mbs; i++) {
4187625dc9aSHong Zhang     jl[i] = mbs; il[i] = 0;
419fff829cfSHong Zhang   }
420fff829cfSHong Zhang 
42198a8e78dSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
4227625dc9aSHong Zhang   nlnk = mbs + 1;
423*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt));
424fff829cfSHong Zhang 
4257625dc9aSHong Zhang   /* initial FreeSpace size is fill*(ai[mbs]+1) */
426*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,ai[mbs]+1),&free_space));
42798a8e78dSHong Zhang   current_space = free_space;
42898a8e78dSHong Zhang 
4297625dc9aSHong Zhang   for (k=0; k<mbs; k++) {  /* for each active row k */
43098a8e78dSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
43198a8e78dSHong Zhang     nzk   = 0;
43298a8e78dSHong Zhang     ncols = ai[rip[k]+1] - ai[rip[k]];
4337625dc9aSHong Zhang     for (j=0; j<ncols; j++) {
4347625dc9aSHong Zhang       i       = *(aj + ai[rip[k]] + j);
4357625dc9aSHong Zhang       cols[j] = rip[i];
4367625dc9aSHong Zhang     }
437*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscLLAdd(ncols,cols,mbs,nlnk,lnk,lnkbt));
43898a8e78dSHong Zhang     nzk += nlnk;
43998a8e78dSHong Zhang 
44098a8e78dSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
44198a8e78dSHong Zhang     prow = jl[k]; /* 1st pivot row */
442fff829cfSHong Zhang 
443fff829cfSHong Zhang     while (prow < k) {
444fff829cfSHong Zhang       nextprow = jl[prow];
44598a8e78dSHong Zhang       /* merge prow into k-th row */
4467625dc9aSHong Zhang       jmin   = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
44798a8e78dSHong Zhang       jmax   = ui[prow+1];
44898a8e78dSHong Zhang       ncols  = jmax-jmin;
4497625dc9aSHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
450*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt));
45198a8e78dSHong Zhang       nzk   += nlnk;
452fff829cfSHong Zhang 
45398a8e78dSHong Zhang       /* update il and jl for prow */
454fff829cfSHong Zhang       if (jmin < jmax) {
455fff829cfSHong Zhang         il[prow] = jmin;
45626fbe8dcSKarl Rupp 
4577625dc9aSHong Zhang         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
458fff829cfSHong Zhang       }
459fff829cfSHong Zhang       prow = nextprow;
460fff829cfSHong Zhang     }
461fff829cfSHong Zhang 
46298a8e78dSHong Zhang     /* if free space is not available, make more free space */
46398a8e78dSHong Zhang     if (current_space->local_remaining<nzk) {
4647625dc9aSHong Zhang       i    = mbs - k + 1; /* num of unfactored rows */
465f91af8c7SBarry Smith       i    = PetscMin(PetscIntMultTruncate(i,nzk), PetscIntMultTruncate(i,i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
466*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFreeSpaceGet(i,&current_space));
46798a8e78dSHong Zhang       reallocs++;
46898a8e78dSHong Zhang     }
46998a8e78dSHong Zhang 
47098a8e78dSHong Zhang     /* copy data into free space, then initialize lnk */
471*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt));
47298a8e78dSHong Zhang 
47398a8e78dSHong Zhang     /* add the k-th row into il and jl */
47498a8e78dSHong Zhang     if (nzk-1 > 0) {
4757625dc9aSHong Zhang       i     = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
476fff829cfSHong Zhang       jl[k] = jl[i]; jl[i] = k;
47798a8e78dSHong Zhang       il[k] = ui[k] + 1;
478fff829cfSHong Zhang     }
4797625dc9aSHong Zhang     ui_ptr[k] = current_space->array;
48026fbe8dcSKarl Rupp 
48198a8e78dSHong Zhang     current_space->array           += nzk;
48298a8e78dSHong Zhang     current_space->local_used      += nzk;
48398a8e78dSHong Zhang     current_space->local_remaining -= nzk;
484fff829cfSHong Zhang 
48598a8e78dSHong Zhang     ui[k+1] = ui[k] + nzk;
486fff829cfSHong Zhang   }
487fff829cfSHong Zhang 
488*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISRestoreIndices(perm,&rip));
489*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree4(ui_ptr,il,jl,cols));
490fff829cfSHong Zhang 
49198a8e78dSHong Zhang   /* destroy list of free space and other temporary array(s) */
492*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(ui[mbs]+1,&uj));
493*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFreeSpaceContiguous(&free_space,uj));
494*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLLDestroy(lnk,lnkbt));
495fff829cfSHong Zhang 
49698a8e78dSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
497*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqSBAIJSetPreallocation(fact,bs,MAT_SKIP_ALLOCATION,NULL));
49898a8e78dSHong Zhang 
49995b5ac22SHong Zhang   b               = (Mat_SeqSBAIJ*)fact->data;
500fff829cfSHong Zhang   b->singlemalloc = PETSC_FALSE;
501e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
502e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
50326fbe8dcSKarl Rupp 
504*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(ui[mbs]+1,&b->a));
50526fbe8dcSKarl Rupp 
50698a8e78dSHong Zhang   b->j    = uj;
50798a8e78dSHong Zhang   b->i    = ui;
508f4259b30SLisandro Dalcin   b->diag = NULL;
509f4259b30SLisandro Dalcin   b->ilen = NULL;
510f4259b30SLisandro Dalcin   b->imax = NULL;
511fff829cfSHong Zhang   b->row  = perm;
51226fbe8dcSKarl Rupp 
513fff829cfSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
51426fbe8dcSKarl Rupp 
515*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectReference((PetscObject)perm));
516fff829cfSHong Zhang   b->icol  = perm;
517*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectReference((PetscObject)perm));
518*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(mbs+1,&b->solve_work));
519*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogObjectMemory((PetscObject)fact,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar))));
5207625dc9aSHong Zhang   b->maxnz = b->nz = ui[mbs];
521fff829cfSHong Zhang 
52295b5ac22SHong Zhang   fact->info.factor_mallocs   = reallocs;
52395b5ac22SHong Zhang   fact->info.fill_ratio_given = fill;
5247625dc9aSHong Zhang   if (ai[mbs] != 0) {
52595b5ac22SHong Zhang     fact->info.fill_ratio_needed = ((PetscReal)ui[mbs])/ai[mbs];
526fff829cfSHong Zhang   } else {
52795b5ac22SHong Zhang     fact->info.fill_ratio_needed = 0.0;
528fff829cfSHong Zhang   }
52995b5ac22SHong Zhang #if defined(PETSC_USE_INFO)
53095b5ac22SHong Zhang   if (ai[mbs] != 0) {
53195b5ac22SHong Zhang     PetscReal af = fact->info.fill_ratio_needed;
532*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscInfo(A,"Reallocs %" PetscInt_FMT " Fill ratio:given %g needed %g\n",reallocs,(double)fill,(double)af));
533*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscInfo(A,"Run with -pc_factor_fill %g or use \n",(double)af));
534*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscInfo(A,"PCFactorSetFill(pc,%g) for best performance.\n",(double)af));
53595b5ac22SHong Zhang   } else {
536*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscInfo(A,"Empty matrix\n"));
53795b5ac22SHong Zhang   }
53895b5ac22SHong Zhang #endif
539*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqSBAIJSetNumericFactorization_inplace(fact,perm_identity));
540064503c5SHong Zhang   PetscFunctionReturn(0);
541064503c5SHong Zhang }
542d595f711SHong Zhang 
5430481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat C,Mat A,const MatFactorInfo *info)
54449b5e25fSSatish Balay {
5454c16a6a6SHong Zhang   Mat_SeqSBAIJ   *a   = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data;
5464c16a6a6SHong Zhang   IS             perm = b->row;
5475d0c19d7SBarry Smith   const PetscInt *ai,*aj,*perm_ptr,mbs=a->mbs,*bi=b->i,*bj=b->j;
5485d0c19d7SBarry Smith   PetscInt       i,j;
5495d0c19d7SBarry Smith   PetscInt       *a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
5503bc0b13bSBarry Smith   PetscInt       bs  =A->rmap->bs,bs2 = a->bs2;
5514c16a6a6SHong Zhang   MatScalar      *ba = b->a,*aa,*ap,*dk,*uik;
5524c16a6a6SHong Zhang   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
55328de702eSHong Zhang   MatScalar      *work;
55413f74950SBarry Smith   PetscInt       *pivots;
5555f8bbccaSHong Zhang   PetscBool      allowzeropivot,zeropivotdetected;
5564c16a6a6SHong Zhang 
5574c16a6a6SHong Zhang   PetscFunctionBegin;
5584c16a6a6SHong Zhang   /* initialization */
559*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCalloc1(bs2*mbs,&rtmp));
560*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc2(mbs,&il,mbs,&jl));
5615f8bbccaSHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
5625f8bbccaSHong Zhang 
5636df5ee2eSHong Zhang   il[0] = 0;
5646df5ee2eSHong Zhang   for (i=0; i<mbs; i++) jl[i] = mbs;
5656df5ee2eSHong Zhang 
566*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc3(bs2,&dk,bs2,&uik,bs,&work));
567*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(bs,&pivots));
5684c16a6a6SHong Zhang 
569*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetIndices(perm,&perm_ptr));
5704c16a6a6SHong Zhang 
5714c16a6a6SHong Zhang   /* check permutation */
5724c16a6a6SHong Zhang   if (!a->permute) {
5734c16a6a6SHong Zhang     ai = a->i; aj = a->j; aa = a->a;
5744c16a6a6SHong Zhang   } else {
5754c16a6a6SHong Zhang     ai   = a->inew; aj = a->jnew;
576*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(bs2*ai[mbs],&aa));
577*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscArraycpy(aa,a->a,bs2*ai[mbs]));
578*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(ai[mbs],&a2anew));
579*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscArraycpy(a2anew,a->a2anew,ai[mbs]));
5804c16a6a6SHong Zhang 
5814c16a6a6SHong Zhang     for (i=0; i<mbs; i++) {
5824c16a6a6SHong Zhang       jmin = ai[i]; jmax = ai[i+1];
5834c16a6a6SHong Zhang       for (j=jmin; j<jmax; j++) {
5844c16a6a6SHong Zhang         while (a2anew[j] != j) {
5854c16a6a6SHong Zhang           k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
5864c16a6a6SHong Zhang           for (k1=0; k1<bs2; k1++) {
5874c16a6a6SHong Zhang             dk[k1]       = aa[k*bs2+k1];
5884c16a6a6SHong Zhang             aa[k*bs2+k1] = aa[j*bs2+k1];
5894c16a6a6SHong Zhang             aa[j*bs2+k1] = dk[k1];
5904c16a6a6SHong Zhang           }
5914c16a6a6SHong Zhang         }
5924c16a6a6SHong Zhang         /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
5934c16a6a6SHong Zhang         if (i > aj[j]) {
5944c16a6a6SHong Zhang           ap = aa + j*bs2;                     /* ptr to the beginning of j-th block of aa */
5954c16a6a6SHong Zhang           for (k=0; k<bs2; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */
5964c16a6a6SHong Zhang           for (k=0; k<bs; k++) {               /* j-th block of aa <- dk^T */
5974c16a6a6SHong Zhang             for (k1=0; k1<bs; k1++) *ap++ = dk[k + bs*k1];
5984c16a6a6SHong Zhang           }
5994c16a6a6SHong Zhang         }
6004c16a6a6SHong Zhang       }
6014c16a6a6SHong Zhang     }
602*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(a2anew));
6034c16a6a6SHong Zhang   }
6044c16a6a6SHong Zhang 
6054c16a6a6SHong Zhang   /* for each row k */
6064c16a6a6SHong Zhang   for (k = 0; k<mbs; k++) {
6074c16a6a6SHong Zhang 
6084c16a6a6SHong Zhang     /*initialize k-th row with elements nonzero in row perm(k) of A */
6094c16a6a6SHong Zhang     jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];
610057f5ba7SHong Zhang 
6114c16a6a6SHong Zhang     ap = aa + jmin*bs2;
6124c16a6a6SHong Zhang     for (j = jmin; j < jmax; j++) {
6134c16a6a6SHong Zhang       vj       = perm_ptr[aj[j]];   /* block col. index */
6144c16a6a6SHong Zhang       rtmp_ptr = rtmp + vj*bs2;
6154c16a6a6SHong Zhang       for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++;
6164c16a6a6SHong Zhang     }
6174c16a6a6SHong Zhang 
6184c16a6a6SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
619*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscArraycpy(dk,rtmp+k*bs2,bs2));
6204c16a6a6SHong Zhang     i    = jl[k]; /* first row to be added to k_th row  */
6214c16a6a6SHong Zhang 
622057f5ba7SHong Zhang     while (i < k) {
6234c16a6a6SHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
6244c16a6a6SHong Zhang 
6254c16a6a6SHong Zhang       /* compute multiplier */
6264c16a6a6SHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
6274c16a6a6SHong Zhang 
6284c16a6a6SHong Zhang       /* uik = -inv(Di)*U_bar(i,k) */
6294c16a6a6SHong Zhang       diag = ba + i*bs2;
6304c16a6a6SHong Zhang       u    = ba + ili*bs2;
631*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscArrayzero(uik,bs2));
63296b95a6bSBarry Smith       PetscKernel_A_gets_A_minus_B_times_C(bs,uik,diag,u);
6334c16a6a6SHong Zhang 
6344c16a6a6SHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
63596b95a6bSBarry Smith       PetscKernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u);
636*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscLogFlops(4.0*bs*bs2));
6374c16a6a6SHong Zhang 
6384c16a6a6SHong Zhang       /* update -U(i,k) */
639*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscArraycpy(ba+ili*bs2,uik,bs2));
6404c16a6a6SHong Zhang 
6414c16a6a6SHong Zhang       /* add multiple of row i to k-th row ... */
6424c16a6a6SHong Zhang       jmin = ili + 1; jmax = bi[i+1];
6434c16a6a6SHong Zhang       if (jmin < jmax) {
6444c16a6a6SHong Zhang         for (j=jmin; j<jmax; j++) {
6454c16a6a6SHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j) */
6464c16a6a6SHong Zhang           rtmp_ptr = rtmp + bj[j]*bs2;
6474c16a6a6SHong Zhang           u        = ba + j*bs2;
64896b95a6bSBarry Smith           PetscKernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u);
6494c16a6a6SHong Zhang         }
650*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscLogFlops(2.0*bs*bs2*(jmax-jmin)));
6514c16a6a6SHong Zhang 
6524c16a6a6SHong Zhang         /* ... add i to row list for next nonzero entry */
6534c16a6a6SHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
6544c16a6a6SHong Zhang         j     = bj[jmin];
6554c16a6a6SHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
6564c16a6a6SHong Zhang       }
6574c16a6a6SHong Zhang       i = nexti;
6584c16a6a6SHong Zhang     }
6594c16a6a6SHong Zhang 
6604c16a6a6SHong Zhang     /* save nonzero entries in k-th row of U ... */
6614c16a6a6SHong Zhang 
6624c16a6a6SHong Zhang     /* invert diagonal block */
6634c16a6a6SHong Zhang     diag = ba+k*bs2;
664*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscArraycpy(diag,dk,bs2));
6655f8bbccaSHong Zhang 
666*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscKernel_A_gets_inverse_A(bs,diag,pivots,work,allowzeropivot,&zeropivotdetected));
6677b6c816cSBarry Smith     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
6684c16a6a6SHong Zhang 
6694c16a6a6SHong Zhang     jmin = bi[k]; jmax = bi[k+1];
6704c16a6a6SHong Zhang     if (jmin < jmax) {
6714c16a6a6SHong Zhang       for (j=jmin; j<jmax; j++) {
6724c16a6a6SHong Zhang         vj       = bj[j];      /* block col. index of U */
6734c16a6a6SHong Zhang         u        = ba + j*bs2;
6744c16a6a6SHong Zhang         rtmp_ptr = rtmp + vj*bs2;
6754c16a6a6SHong Zhang         for (k1=0; k1<bs2; k1++) {
6764c16a6a6SHong Zhang           *u++        = *rtmp_ptr;
6774c16a6a6SHong Zhang           *rtmp_ptr++ = 0.0;
6784c16a6a6SHong Zhang         }
6794c16a6a6SHong Zhang       }
6804c16a6a6SHong Zhang 
6814c16a6a6SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
6824c16a6a6SHong Zhang       il[k] = jmin;
6834c16a6a6SHong Zhang       i     = bj[jmin];
6844c16a6a6SHong Zhang       jl[k] = jl[i]; jl[i] = k;
6854c16a6a6SHong Zhang     }
6864c16a6a6SHong Zhang   }
6874c16a6a6SHong Zhang 
688*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(rtmp));
689*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree2(il,jl));
690*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree3(dk,uik,work));
691*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(pivots));
6924c16a6a6SHong Zhang   if (a->permute) {
693*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(aa));
6944c16a6a6SHong Zhang   }
6954c16a6a6SHong Zhang 
696*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISRestoreIndices(perm,&perm_ptr));
69726fbe8dcSKarl Rupp 
6984f79d315SHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_N_inplace;
6994f79d315SHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_N_inplace;
7004f79d315SHong Zhang   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_N_inplace;
7014f79d315SHong Zhang   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_N_inplace;
702db4efbfdSBarry Smith 
7034c16a6a6SHong Zhang   C->assembled    = PETSC_TRUE;
7044c16a6a6SHong Zhang   C->preallocated = PETSC_TRUE;
70526fbe8dcSKarl Rupp 
706*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops(1.3333*bs*bs2*b->mbs)); /* from inverting diagonal blocks */
7074c16a6a6SHong Zhang   PetscFunctionReturn(0);
7084c16a6a6SHong Zhang }
709d4edadd4SHong Zhang 
7100481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info)
711671cb588SHong Zhang {
712671cb588SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data;
71313f74950SBarry Smith   PetscInt       i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
71413f74950SBarry Smith   PetscInt       *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
7153bc0b13bSBarry Smith   PetscInt       bs  =A->rmap->bs,bs2 = a->bs2;
716671cb588SHong Zhang   MatScalar      *ba = b->a,*aa,*ap,*dk,*uik;
717671cb588SHong Zhang   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
71828de702eSHong Zhang   MatScalar      *work;
71913f74950SBarry Smith   PetscInt       *pivots;
7205f8bbccaSHong Zhang   PetscBool      allowzeropivot,zeropivotdetected;
721671cb588SHong Zhang 
722671cb588SHong Zhang   PetscFunctionBegin;
723*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCalloc1(bs2*mbs,&rtmp));
724*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc2(mbs,&il,mbs,&jl));
7256df5ee2eSHong Zhang   il[0] = 0;
7266df5ee2eSHong Zhang   for (i=0; i<mbs; i++) jl[i] = mbs;
7276df5ee2eSHong Zhang 
728*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc3(bs2,&dk,bs2,&uik,bs,&work));
729*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(bs,&pivots));
7305f8bbccaSHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
731671cb588SHong Zhang 
732671cb588SHong Zhang   ai = a->i; aj = a->j; aa = a->a;
733671cb588SHong Zhang 
734671cb588SHong Zhang   /* for each row k */
735671cb588SHong Zhang   for (k = 0; k<mbs; k++) {
736671cb588SHong Zhang 
737671cb588SHong Zhang     /*initialize k-th row with elements nonzero in row k of A */
738671cb588SHong Zhang     jmin = ai[k]; jmax = ai[k+1];
739671cb588SHong Zhang     ap   = aa + jmin*bs2;
740671cb588SHong Zhang     for (j = jmin; j < jmax; j++) {
741671cb588SHong Zhang       vj       = aj[j];   /* block col. index */
742671cb588SHong Zhang       rtmp_ptr = rtmp + vj*bs2;
743671cb588SHong Zhang       for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++;
744671cb588SHong Zhang     }
745671cb588SHong Zhang 
746671cb588SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
747*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscArraycpy(dk,rtmp+k*bs2,bs2));
748671cb588SHong Zhang     i    = jl[k]; /* first row to be added to k_th row  */
749671cb588SHong Zhang 
750057f5ba7SHong Zhang     while (i < k) {
751671cb588SHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
752671cb588SHong Zhang 
753671cb588SHong Zhang       /* compute multiplier */
754671cb588SHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
755671cb588SHong Zhang 
756671cb588SHong Zhang       /* uik = -inv(Di)*U_bar(i,k) */
757671cb588SHong Zhang       diag = ba + i*bs2;
758671cb588SHong Zhang       u    = ba + ili*bs2;
759*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscArrayzero(uik,bs2));
76096b95a6bSBarry Smith       PetscKernel_A_gets_A_minus_B_times_C(bs,uik,diag,u);
761671cb588SHong Zhang 
762671cb588SHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
76396b95a6bSBarry Smith       PetscKernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u);
764*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscLogFlops(2.0*bs*bs2));
765671cb588SHong Zhang 
766671cb588SHong Zhang       /* update -U(i,k) */
767*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscArraycpy(ba+ili*bs2,uik,bs2));
768671cb588SHong Zhang 
769671cb588SHong Zhang       /* add multiple of row i to k-th row ... */
770671cb588SHong Zhang       jmin = ili + 1; jmax = bi[i+1];
771671cb588SHong Zhang       if (jmin < jmax) {
772671cb588SHong Zhang         for (j=jmin; j<jmax; j++) {
773671cb588SHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j) */
774671cb588SHong Zhang           rtmp_ptr = rtmp + bj[j]*bs2;
775671cb588SHong Zhang           u        = ba + j*bs2;
77696b95a6bSBarry Smith           PetscKernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u);
777671cb588SHong Zhang         }
778*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscLogFlops(2.0*bs*bs2*(jmax-jmin)));
779671cb588SHong Zhang 
780671cb588SHong Zhang         /* ... add i to row list for next nonzero entry */
781671cb588SHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
782671cb588SHong Zhang         j     = bj[jmin];
783671cb588SHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
784671cb588SHong Zhang       }
785671cb588SHong Zhang       i = nexti;
786671cb588SHong Zhang     }
787671cb588SHong Zhang 
788671cb588SHong Zhang     /* save nonzero entries in k-th row of U ... */
789671cb588SHong Zhang 
790671cb588SHong Zhang     /* invert diagonal block */
791671cb588SHong Zhang     diag = ba+k*bs2;
792*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscArraycpy(diag,dk,bs2));
7935f8bbccaSHong Zhang 
794*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscKernel_A_gets_inverse_A(bs,diag,pivots,work,allowzeropivot,&zeropivotdetected));
7957b6c816cSBarry Smith     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
796671cb588SHong Zhang 
797671cb588SHong Zhang     jmin = bi[k]; jmax = bi[k+1];
798671cb588SHong Zhang     if (jmin < jmax) {
799671cb588SHong Zhang       for (j=jmin; j<jmax; j++) {
800671cb588SHong Zhang         vj       = bj[j];      /* block col. index of U */
801671cb588SHong Zhang         u        = ba + j*bs2;
802671cb588SHong Zhang         rtmp_ptr = rtmp + vj*bs2;
803671cb588SHong Zhang         for (k1=0; k1<bs2; k1++) {
804671cb588SHong Zhang           *u++        = *rtmp_ptr;
805671cb588SHong Zhang           *rtmp_ptr++ = 0.0;
806671cb588SHong Zhang         }
807671cb588SHong Zhang       }
808671cb588SHong Zhang 
809671cb588SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
810671cb588SHong Zhang       il[k] = jmin;
811671cb588SHong Zhang       i     = bj[jmin];
812671cb588SHong Zhang       jl[k] = jl[i]; jl[i] = k;
813671cb588SHong Zhang     }
814671cb588SHong Zhang   }
815671cb588SHong Zhang 
816*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(rtmp));
817*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree2(il,jl));
818*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree3(dk,uik,work));
819*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(pivots));
820671cb588SHong Zhang 
8214f79d315SHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_N_NaturalOrdering_inplace;
8224f79d315SHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_N_NaturalOrdering_inplace;
8234f79d315SHong Zhang   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_N_NaturalOrdering_inplace;
8244f79d315SHong Zhang   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_N_NaturalOrdering_inplace;
825671cb588SHong Zhang   C->assembled           = PETSC_TRUE;
826671cb588SHong Zhang   C->preallocated        = PETSC_TRUE;
82726fbe8dcSKarl Rupp 
828*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops(1.3333*bs*bs2*b->mbs)); /* from inverting diagonal blocks */
829671cb588SHong Zhang   PetscFunctionReturn(0);
830671cb588SHong Zhang }
831671cb588SHong Zhang 
83249b5e25fSSatish Balay /*
833fcf159c0SHong Zhang     Numeric U^T*D*U factorization for SBAIJ format. Modified from SNF of YSMP.
834cc0c071aSHong Zhang     Version for blocks 2 by 2.
83549b5e25fSSatish Balay */
8360481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat C,Mat A,const MatFactorInfo *info)
83749b5e25fSSatish Balay {
83891602c66SHong Zhang   Mat_SeqSBAIJ   *a   = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data;
839cc0c071aSHong Zhang   IS             perm = b->row;
8405d0c19d7SBarry Smith   const PetscInt *ai,*aj,*perm_ptr;
8415d0c19d7SBarry Smith   PetscInt       i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
8425d0c19d7SBarry Smith   PetscInt       *a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
843d8c74875SBarry Smith   MatScalar      *ba = b->a,*aa,*ap;
844d8c74875SBarry Smith   MatScalar      *u,*diag,*rtmp,*rtmp_ptr,dk[4],uik[4];
84514d2772eSHong Zhang   PetscReal      shift = info->shiftamount;
846a455e926SHong Zhang   PetscBool      allowzeropivot,zeropivotdetected;
84749b5e25fSSatish Balay 
84849b5e25fSSatish Balay   PetscFunctionBegin;
8490164db54SHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
8500164db54SHong Zhang 
85191602c66SHong Zhang   /* initialization */
85291602c66SHong Zhang   /* il and jl record the first nonzero element in each row of the accessing
85391602c66SHong Zhang      window U(0:k, k:mbs-1).
85491602c66SHong Zhang      jl:    list of rows to be added to uneliminated rows
85591602c66SHong Zhang             i>= k: jl(i) is the first row to be added to row i
85691602c66SHong Zhang             i<  k: jl(i) is the row following row i in some list of rows
85791602c66SHong Zhang             jl(i) = mbs indicates the end of a list
85891602c66SHong Zhang      il(i): points to the first nonzero element in columns k,...,mbs-1 of
85991602c66SHong Zhang             row i of U */
860*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCalloc1(4*mbs,&rtmp));
861*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc2(mbs,&il,mbs,&jl));
8626df5ee2eSHong Zhang   il[0] = 0;
8636df5ee2eSHong Zhang   for (i=0; i<mbs; i++) jl[i] = mbs;
8646df5ee2eSHong Zhang 
865*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetIndices(perm,&perm_ptr));
866cc0c071aSHong Zhang 
867cc0c071aSHong Zhang   /* check permutation */
868cc0c071aSHong Zhang   if (!a->permute) {
869cc0c071aSHong Zhang     ai = a->i; aj = a->j; aa = a->a;
870cc0c071aSHong Zhang   } else {
871cc0c071aSHong Zhang     ai   = a->inew; aj = a->jnew;
872*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(4*ai[mbs],&aa));
873*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscArraycpy(aa,a->a,4*ai[mbs]));
874*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(ai[mbs],&a2anew));
875*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscArraycpy(a2anew,a->a2anew,ai[mbs]));
876cc0c071aSHong Zhang 
877cc0c071aSHong Zhang     for (i=0; i<mbs; i++) {
878cc0c071aSHong Zhang       jmin = ai[i]; jmax = ai[i+1];
879cc0c071aSHong Zhang       for (j=jmin; j<jmax; j++) {
880cc0c071aSHong Zhang         while (a2anew[j] != j) {
881cc0c071aSHong Zhang           k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
882cc0c071aSHong Zhang           for (k1=0; k1<4; k1++) {
883cc0c071aSHong Zhang             dk[k1]     = aa[k*4+k1];
884cc0c071aSHong Zhang             aa[k*4+k1] = aa[j*4+k1];
885cc0c071aSHong Zhang             aa[j*4+k1] = dk[k1];
886cc0c071aSHong Zhang           }
887cc0c071aSHong Zhang         }
888cc0c071aSHong Zhang         /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
889cc0c071aSHong Zhang         if (i > aj[j]) {
890cc0c071aSHong Zhang           ap    = aa + j*4;  /* ptr to the beginning of the block */
891cc0c071aSHong Zhang           dk[1] = ap[1];     /* swap ap[1] and ap[2] */
892cc0c071aSHong Zhang           ap[1] = ap[2];
893cc0c071aSHong Zhang           ap[2] = dk[1];
894cc0c071aSHong Zhang         }
895cc0c071aSHong Zhang       }
896cc0c071aSHong Zhang     }
897*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(a2anew));
898cc0c071aSHong Zhang   }
8993845f261SHong Zhang 
90091602c66SHong Zhang   /* for each row k */
90191602c66SHong Zhang   for (k = 0; k<mbs; k++) {
90291602c66SHong Zhang 
90391602c66SHong Zhang     /*initialize k-th row with elements nonzero in row perm(k) of A */
904cc0c071aSHong Zhang     jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];
905cc0c071aSHong Zhang     ap   = aa + jmin*4;
90691602c66SHong Zhang     for (j = jmin; j < jmax; j++) {
907cc0c071aSHong Zhang       vj       = perm_ptr[aj[j]];   /* block col. index */
908cc0c071aSHong Zhang       rtmp_ptr = rtmp + vj*4;
909cc0c071aSHong Zhang       for (i=0; i<4; i++) *rtmp_ptr++ = *ap++;
91091602c66SHong Zhang     }
91191602c66SHong Zhang 
91291602c66SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
913*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscArraycpy(dk,rtmp+k*4,4));
91491602c66SHong Zhang     i    = jl[k]; /* first row to be added to k_th row  */
91591602c66SHong Zhang 
916057f5ba7SHong Zhang     while (i < k) {
91791602c66SHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
91891602c66SHong Zhang 
9193845f261SHong Zhang       /* compute multiplier */
92091602c66SHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
9213845f261SHong Zhang 
9223845f261SHong Zhang       /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */
923cc0c071aSHong Zhang       diag   = ba + i*4;
924cc0c071aSHong Zhang       u      = ba + ili*4;
925cc0c071aSHong Zhang       uik[0] = -(diag[0]*u[0] + diag[2]*u[1]);
926cc0c071aSHong Zhang       uik[1] = -(diag[1]*u[0] + diag[3]*u[1]);
927cc0c071aSHong Zhang       uik[2] = -(diag[0]*u[2] + diag[2]*u[3]);
928cc0c071aSHong Zhang       uik[3] = -(diag[1]*u[2] + diag[3]*u[3]);
9293845f261SHong Zhang 
9303845f261SHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */
931cc0c071aSHong Zhang       dk[0] += uik[0]*u[0] + uik[1]*u[1];
932cc0c071aSHong Zhang       dk[1] += uik[2]*u[0] + uik[3]*u[1];
933cc0c071aSHong Zhang       dk[2] += uik[0]*u[2] + uik[1]*u[3];
934cc0c071aSHong Zhang       dk[3] += uik[2]*u[2] + uik[3]*u[3];
9353845f261SHong Zhang 
936*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscLogFlops(16.0*2.0));
937187a9f4bSHong Zhang 
9383845f261SHong Zhang       /* update -U(i,k): ba[ili] = uik */
939*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscArraycpy(ba+ili*4,uik,4));
94091602c66SHong Zhang 
94191602c66SHong Zhang       /* add multiple of row i to k-th row ... */
94291602c66SHong Zhang       jmin = ili + 1; jmax = bi[i+1];
94391602c66SHong Zhang       if (jmin < jmax) {
9443845f261SHong Zhang         for (j=jmin; j<jmax; j++) {
9453845f261SHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */
946cc0c071aSHong Zhang           rtmp_ptr     = rtmp + bj[j]*4;
947cc0c071aSHong Zhang           u            = ba + j*4;
948cc0c071aSHong Zhang           rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1];
949cc0c071aSHong Zhang           rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1];
950cc0c071aSHong Zhang           rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3];
951cc0c071aSHong Zhang           rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3];
9523845f261SHong Zhang         }
953*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscLogFlops(16.0*(jmax-jmin)));
9543845f261SHong Zhang 
95591602c66SHong Zhang         /* ... add i to row list for next nonzero entry */
95691602c66SHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
95791602c66SHong Zhang         j     = bj[jmin];
95891602c66SHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
95991602c66SHong Zhang       }
960a1723e09SHong Zhang       i = nexti;
96191602c66SHong Zhang     }
962cc0c071aSHong Zhang 
96391602c66SHong Zhang     /* save nonzero entries in k-th row of U ... */
9643845f261SHong Zhang 
965cc0c071aSHong Zhang     /* invert diagonal block */
966cc0c071aSHong Zhang     diag = ba+k*4;
967*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscArraycpy(diag,dk,4));
968*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscKernel_A_gets_inverse_A_2(diag,shift,allowzeropivot,&zeropivotdetected));
9697b6c816cSBarry Smith     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
9703845f261SHong Zhang 
97191602c66SHong Zhang     jmin = bi[k]; jmax = bi[k+1];
97291602c66SHong Zhang     if (jmin < jmax) {
97391602c66SHong Zhang       for (j=jmin; j<jmax; j++) {
974cc0c071aSHong Zhang         vj       = bj[j];      /* block col. index of U */
975cc0c071aSHong Zhang         u        = ba + j*4;
976cc0c071aSHong Zhang         rtmp_ptr = rtmp + vj*4;
977cc0c071aSHong Zhang         for (k1=0; k1<4; k1++) {
978cc0c071aSHong Zhang           *u++        = *rtmp_ptr;
979cc0c071aSHong Zhang           *rtmp_ptr++ = 0.0;
980cc0c071aSHong Zhang         }
981cc0c071aSHong Zhang       }
9823845f261SHong Zhang 
98391602c66SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
98491602c66SHong Zhang       il[k] = jmin;
98591602c66SHong Zhang       i     = bj[jmin];
98691602c66SHong Zhang       jl[k] = jl[i]; jl[i] = k;
98791602c66SHong Zhang     }
98891602c66SHong Zhang   }
9893845f261SHong Zhang 
990*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(rtmp));
991*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree2(il,jl));
99291602c66SHong Zhang   if (a->permute) {
993*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(aa));
99491602c66SHong Zhang   }
995*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISRestoreIndices(perm,&perm_ptr));
99626fbe8dcSKarl Rupp 
9974f79d315SHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_2_inplace;
9984f79d315SHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_2_inplace;
99949b5e25fSSatish Balay   C->assembled           = PETSC_TRUE;
10005c0bcdfcSHong Zhang   C->preallocated        = PETSC_TRUE;
100126fbe8dcSKarl Rupp 
1002*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops(1.3333*8*b->mbs)); /* from inverting diagonal blocks */
100349b5e25fSSatish Balay   PetscFunctionReturn(0);
100449b5e25fSSatish Balay }
100591602c66SHong Zhang 
100649b5e25fSSatish Balay /*
100749b5e25fSSatish Balay       Version for when blocks are 2 by 2 Using natural ordering
100849b5e25fSSatish Balay */
10090481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info)
101049b5e25fSSatish Balay {
1011ab27746eSHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data;
101213f74950SBarry Smith   PetscInt       i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
101313f74950SBarry Smith   PetscInt       *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
1014d8c74875SBarry Smith   MatScalar      *ba = b->a,*aa,*ap,dk[8],uik[8];
1015ab27746eSHong Zhang   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
101614d2772eSHong Zhang   PetscReal      shift = info->shiftamount;
1017a455e926SHong Zhang   PetscBool      allowzeropivot,zeropivotdetected;
101849b5e25fSSatish Balay 
101949b5e25fSSatish Balay   PetscFunctionBegin;
10200164db54SHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
10210164db54SHong Zhang 
1022ab27746eSHong Zhang   /* initialization */
1023ab27746eSHong Zhang   /* il and jl record the first nonzero element in each row of the accessing
1024ab27746eSHong Zhang      window U(0:k, k:mbs-1).
1025ab27746eSHong Zhang      jl:    list of rows to be added to uneliminated rows
1026ab27746eSHong Zhang             i>= k: jl(i) is the first row to be added to row i
1027ab27746eSHong Zhang             i<  k: jl(i) is the row following row i in some list of rows
1028ab27746eSHong Zhang             jl(i) = mbs indicates the end of a list
1029ab27746eSHong Zhang      il(i): points to the first nonzero element in columns k,...,mbs-1 of
1030ab27746eSHong Zhang             row i of U */
1031*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCalloc1(4*mbs,&rtmp));
1032*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc2(mbs,&il,mbs,&jl));
10336df5ee2eSHong Zhang   il[0] = 0;
10346df5ee2eSHong Zhang   for (i=0; i<mbs; i++) jl[i] = mbs;
10356df5ee2eSHong Zhang 
1036ab27746eSHong Zhang   ai = a->i; aj = a->j; aa = a->a;
1037ab27746eSHong Zhang 
1038ab27746eSHong Zhang   /* for each row k */
1039ab27746eSHong Zhang   for (k = 0; k<mbs; k++) {
1040ab27746eSHong Zhang 
1041ab27746eSHong Zhang     /*initialize k-th row with elements nonzero in row k of A */
1042ab27746eSHong Zhang     jmin = ai[k]; jmax = ai[k+1];
1043ab27746eSHong Zhang     ap   = aa + jmin*4;
1044ab27746eSHong Zhang     for (j = jmin; j < jmax; j++) {
1045ab27746eSHong Zhang       vj       = aj[j];   /* block col. index */
1046ab27746eSHong Zhang       rtmp_ptr = rtmp + vj*4;
1047ab27746eSHong Zhang       for (i=0; i<4; i++) *rtmp_ptr++ = *ap++;
104849b5e25fSSatish Balay     }
1049ab27746eSHong Zhang 
1050ab27746eSHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
1051*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscArraycpy(dk,rtmp+k*4,4));
1052ab27746eSHong Zhang     i    = jl[k]; /* first row to be added to k_th row  */
1053ab27746eSHong Zhang 
1054057f5ba7SHong Zhang     while (i < k) {
1055ab27746eSHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
1056ab27746eSHong Zhang 
1057ab27746eSHong Zhang       /* compute multiplier */
1058ab27746eSHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
1059ab27746eSHong Zhang 
1060ab27746eSHong Zhang       /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */
1061ab27746eSHong Zhang       diag   = ba + i*4;
1062ab27746eSHong Zhang       u      = ba + ili*4;
1063ab27746eSHong Zhang       uik[0] = -(diag[0]*u[0] + diag[2]*u[1]);
1064ab27746eSHong Zhang       uik[1] = -(diag[1]*u[0] + diag[3]*u[1]);
1065ab27746eSHong Zhang       uik[2] = -(diag[0]*u[2] + diag[2]*u[3]);
1066ab27746eSHong Zhang       uik[3] = -(diag[1]*u[2] + diag[3]*u[3]);
1067ab27746eSHong Zhang 
1068ab27746eSHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */
1069ab27746eSHong Zhang       dk[0] += uik[0]*u[0] + uik[1]*u[1];
1070ab27746eSHong Zhang       dk[1] += uik[2]*u[0] + uik[3]*u[1];
1071ab27746eSHong Zhang       dk[2] += uik[0]*u[2] + uik[1]*u[3];
1072ab27746eSHong Zhang       dk[3] += uik[2]*u[2] + uik[3]*u[3];
1073ab27746eSHong Zhang 
1074*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscLogFlops(16.0*2.0));
1075187a9f4bSHong Zhang 
1076ab27746eSHong Zhang       /* update -U(i,k): ba[ili] = uik */
1077*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscArraycpy(ba+ili*4,uik,4));
1078ab27746eSHong Zhang 
1079ab27746eSHong Zhang       /* add multiple of row i to k-th row ... */
1080ab27746eSHong Zhang       jmin = ili + 1; jmax = bi[i+1];
1081ab27746eSHong Zhang       if (jmin < jmax) {
1082ab27746eSHong Zhang         for (j=jmin; j<jmax; j++) {
1083ab27746eSHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */
1084ab27746eSHong Zhang           rtmp_ptr     = rtmp + bj[j]*4;
1085ab27746eSHong Zhang           u            = ba + j*4;
1086ab27746eSHong Zhang           rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1];
1087ab27746eSHong Zhang           rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1];
1088ab27746eSHong Zhang           rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3];
1089ab27746eSHong Zhang           rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3];
109049b5e25fSSatish Balay         }
1091*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscLogFlops(16.0*(jmax-jmin)));
1092ab27746eSHong Zhang 
1093ab27746eSHong Zhang         /* ... add i to row list for next nonzero entry */
1094ab27746eSHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
1095ab27746eSHong Zhang         j     = bj[jmin];
1096ab27746eSHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
109749b5e25fSSatish Balay       }
1098ab27746eSHong Zhang       i = nexti;
109949b5e25fSSatish Balay     }
1100ab27746eSHong Zhang 
1101ab27746eSHong Zhang     /* save nonzero entries in k-th row of U ... */
1102ab27746eSHong Zhang 
110349b5e25fSSatish Balay     /* invert diagonal block */
1104ab27746eSHong Zhang     diag = ba+k*4;
1105*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscArraycpy(diag,dk,4));
1106*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscKernel_A_gets_inverse_A_2(diag,shift,allowzeropivot,&zeropivotdetected));
11077b6c816cSBarry Smith     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1108ab27746eSHong Zhang 
1109ab27746eSHong Zhang     jmin = bi[k]; jmax = bi[k+1];
1110ab27746eSHong Zhang     if (jmin < jmax) {
1111ab27746eSHong Zhang       for (j=jmin; j<jmax; j++) {
1112ab27746eSHong Zhang         vj       = bj[j];      /* block col. index of U */
1113ab27746eSHong Zhang         u        = ba + j*4;
1114ab27746eSHong Zhang         rtmp_ptr = rtmp + vj*4;
1115ab27746eSHong Zhang         for (k1=0; k1<4; k1++) {
1116ab27746eSHong Zhang           *u++        = *rtmp_ptr;
1117ab27746eSHong Zhang           *rtmp_ptr++ = 0.0;
1118ab27746eSHong Zhang         }
1119ab27746eSHong Zhang       }
1120ab27746eSHong Zhang 
1121ab27746eSHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
1122ab27746eSHong Zhang       il[k] = jmin;
1123ab27746eSHong Zhang       i     = bj[jmin];
1124ab27746eSHong Zhang       jl[k] = jl[i]; jl[i] = k;
1125ab27746eSHong Zhang     }
112649b5e25fSSatish Balay   }
112749b5e25fSSatish Balay 
1128*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(rtmp));
1129*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree2(il,jl));
1130ab27746eSHong Zhang 
11314f79d315SHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_2_NaturalOrdering_inplace;
11324f79d315SHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_2_NaturalOrdering_inplace;
11334f79d315SHong Zhang   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_2_NaturalOrdering_inplace;
11344f79d315SHong Zhang   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_2_NaturalOrdering_inplace;
113549b5e25fSSatish Balay   C->assembled           = PETSC_TRUE;
11365c0bcdfcSHong Zhang   C->preallocated        = PETSC_TRUE;
113726fbe8dcSKarl Rupp 
1138*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops(1.3333*8*b->mbs)); /* from inverting diagonal blocks */
113949b5e25fSSatish Balay   PetscFunctionReturn(0);
114049b5e25fSSatish Balay }
114149b5e25fSSatish Balay 
114249b5e25fSSatish Balay /*
114398a8e78dSHong Zhang     Numeric U^T*D*U factorization for SBAIJ format.
114491602c66SHong Zhang     Version for blocks are 1 by 1.
114549b5e25fSSatish Balay */
1146d595f711SHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_inplace(Mat C,Mat A,const MatFactorInfo *info)
114749b5e25fSSatish Balay {
114849b5e25fSSatish Balay   Mat_SeqSBAIJ   *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ*)C->data;
114949b5e25fSSatish Balay   IS             ip=b->row;
11505d0c19d7SBarry Smith   const PetscInt *ai,*aj,*rip;
11515d0c19d7SBarry Smith   PetscInt       *a2anew,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j,*bcol;
1152997a0949SHong Zhang   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
1153997a0949SHong Zhang   MatScalar      *rtmp,*ba=b->a,*bval,*aa,dk,uikdi;
11540e95ead3SHong Zhang   PetscReal      rs;
11550e95ead3SHong Zhang   FactorShiftCtx sctx;
115649b5e25fSSatish Balay 
115749b5e25fSSatish Balay   PetscFunctionBegin;
11580e95ead3SHong Zhang   /* MatPivotSetUp(): initialize shift context sctx */
1159*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMemzero(&sctx,sizeof(FactorShiftCtx)));
11603cea4cbeSHong Zhang 
1161*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetIndices(ip,&rip));
1162cb718733SHong Zhang   if (!a->permute) {
11632d007305SHong Zhang     ai = a->i; aj = a->j; aa = a->a;
11642d007305SHong Zhang   } else {
11652d007305SHong Zhang     ai     = a->inew; aj = a->jnew;
1166fff829cfSHong Zhang     nz     = ai[mbs];
1167*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(nz,&aa));
1168fff829cfSHong Zhang     a2anew = a->a2anew;
1169997a0949SHong Zhang     bval   = a->a;
1170fff829cfSHong Zhang     for (j=0; j<nz; j++) {
1171997a0949SHong Zhang       aa[a2anew[j]] = *(bval++);
11722d007305SHong Zhang     }
11732d007305SHong Zhang   }
11742d007305SHong Zhang 
117591602c66SHong Zhang   /* initialization */
117649b5e25fSSatish Balay   /* il and jl record the first nonzero element in each row of the accessing
117749b5e25fSSatish Balay      window U(0:k, k:mbs-1).
117849b5e25fSSatish Balay      jl:    list of rows to be added to uneliminated rows
117949b5e25fSSatish Balay             i>= k: jl(i) is the first row to be added to row i
118049b5e25fSSatish Balay             i<  k: jl(i) is the row following row i in some list of rows
118149b5e25fSSatish Balay             jl(i) = mbs indicates the end of a list
118249b5e25fSSatish Balay      il(i): points to the first nonzero element in columns k,...,mbs-1 of
118349b5e25fSSatish Balay             row i of U */
1184*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc3(mbs,&rtmp,mbs,&il,mbs,&jl));
1185997a0949SHong Zhang 
1186997a0949SHong Zhang   do {
118707b50cabSHong Zhang     sctx.newshift = PETSC_FALSE;
11886df5ee2eSHong Zhang     il[0] = 0;
118949b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
11906df5ee2eSHong Zhang       rtmp[i] = 0.0; jl[i] = mbs;
119149b5e25fSSatish Balay     }
1192997a0949SHong Zhang 
119349b5e25fSSatish Balay     for (k = 0; k<mbs; k++) {
1194997a0949SHong Zhang       /*initialize k-th row by the perm[k]-th row of A */
119549b5e25fSSatish Balay       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
11967701de02SHong Zhang       bval = ba + bi[k];
119749b5e25fSSatish Balay       for (j = jmin; j < jmax; j++) {
1198997a0949SHong Zhang         col       = rip[aj[j]];
1199997a0949SHong Zhang         rtmp[col] = aa[j];
12007701de02SHong Zhang         *bval++   = 0.0; /* for in-place factorization */
120149b5e25fSSatish Balay       }
12023cea4cbeSHong Zhang 
12033cea4cbeSHong Zhang       /* shift the diagonal of the matrix */
12043cea4cbeSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
120549b5e25fSSatish Balay 
120691602c66SHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
120749b5e25fSSatish Balay       dk = rtmp[k];
120849b5e25fSSatish Balay       i  = jl[k]; /* first row to be added to k_th row  */
120949b5e25fSSatish Balay 
1210057f5ba7SHong Zhang       while (i < k) {
121149b5e25fSSatish Balay         nexti = jl[i]; /* next row to be added to k_th row */
121249b5e25fSSatish Balay 
1213fff829cfSHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
121449b5e25fSSatish Balay         ili     = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
1215997a0949SHong Zhang         uikdi   = -ba[ili]*ba[bi[i]]; /* diagonal(k) */
121649b5e25fSSatish Balay         dk     += uikdi*ba[ili];
1217658e7b3eSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
121849b5e25fSSatish Balay 
1219997a0949SHong Zhang         /* add multiple of row i to k-th row */
122049b5e25fSSatish Balay         jmin = ili + 1; jmax = bi[i+1];
122149b5e25fSSatish Balay         if (jmin < jmax) {
122249b5e25fSSatish Balay           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
1223*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscLogFlops(2.0*(jmax-jmin)));
1224187a9f4bSHong Zhang 
1225fff829cfSHong Zhang           /* update il and jl for row i */
1226fff829cfSHong Zhang           il[i] = jmin;
1227fff829cfSHong Zhang           j     = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
122849b5e25fSSatish Balay         }
1229ab27746eSHong Zhang         i = nexti;
123049b5e25fSSatish Balay       }
123149b5e25fSSatish Balay 
12323cea4cbeSHong Zhang       /* shift the diagonals when zero pivot is detected */
12333cea4cbeSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
12343cea4cbeSHong Zhang       rs   = 0.0;
1235997a0949SHong Zhang       jmin = bi[k]+1;
1236997a0949SHong Zhang       nz   = bi[k+1] - jmin;
1237997a0949SHong Zhang       if (nz) {
1238997a0949SHong Zhang         bcol = bj + jmin;
1239997a0949SHong Zhang         while (nz--) {
124089f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
124189f845c8SHong Zhang           bcol++;
1242997a0949SHong Zhang         }
1243997a0949SHong Zhang       }
12443cea4cbeSHong Zhang 
12453cea4cbeSHong Zhang       sctx.rs = rs;
12463cea4cbeSHong Zhang       sctx.pv = dk;
1247*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatPivotCheck(C,A,info,&sctx,k));
124807b50cabSHong Zhang       if (sctx.newshift) break;    /* sctx.shift_amount is updated */
12490e95ead3SHong Zhang       dk = sctx.pv;
125049b5e25fSSatish Balay 
1251997a0949SHong Zhang       /* copy data into U(k,:) */
125298a8e78dSHong Zhang       ba[bi[k]] = 1.0/dk; /* U(k,k) */
1253fff829cfSHong Zhang       jmin      = bi[k]+1; jmax = bi[k+1];
125449b5e25fSSatish Balay       if (jmin < jmax) {
125549b5e25fSSatish Balay         for (j=jmin; j<jmax; j++) {
1256997a0949SHong Zhang           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
125749b5e25fSSatish Balay         }
125898a8e78dSHong Zhang         /* add the k-th row into il and jl */
125949b5e25fSSatish Balay         il[k] = jmin;
126098a8e78dSHong Zhang         i     = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
126149b5e25fSSatish Balay       }
126249b5e25fSSatish Balay     }
126307b50cabSHong Zhang   } while (sctx.newshift);
1264*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree3(rtmp,il,jl));
1265*5f80ce2aSJacob Faibussowitsch   if (a->permute) CHKERRQ(PetscFree(aa));
126649b5e25fSSatish Balay 
1267*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISRestoreIndices(ip,&rip));
126826fbe8dcSKarl Rupp 
12690a3c351aSHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_1_inplace;
12704f79d315SHong Zhang   C->ops->solves         = MatSolves_SeqSBAIJ_1_inplace;
12710a3c351aSHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_1_inplace;
12720a3c351aSHong Zhang   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_1_inplace;
12730a3c351aSHong Zhang   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_1_inplace;
127449b5e25fSSatish Balay   C->assembled           = PETSC_TRUE;
12755c0bcdfcSHong Zhang   C->preallocated        = PETSC_TRUE;
127626fbe8dcSKarl Rupp 
1277*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops(C->rmap->N));
12783cea4cbeSHong Zhang   if (sctx.nshift) {
1279f4db908eSBarry Smith     if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
1280*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscInfo(A,"number of shiftnz tries %" PetscInt_FMT ", shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount));
1281f4db908eSBarry Smith     } else if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
1282*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscInfo(A,"number of shiftpd tries %" PetscInt_FMT ", shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount));
1283997a0949SHong Zhang     }
1284997a0949SHong Zhang   }
128549b5e25fSSatish Balay   PetscFunctionReturn(0);
128649b5e25fSSatish Balay }
128749b5e25fSSatish Balay 
1288671cb588SHong Zhang /*
128980d3d5a7SHong Zhang   Version for when blocks are 1 by 1 Using natural ordering under new datastructure
129080d3d5a7SHong Zhang   Modified from MatCholeskyFactorNumeric_SeqAIJ()
1291671cb588SHong Zhang */
1292d595f711SHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info)
1293d595f711SHong Zhang {
1294d595f711SHong Zhang   Mat_SeqSBAIJ   *a=(Mat_SeqSBAIJ*)A->data;
12957b056e98SHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)B->data;
1296d595f711SHong Zhang   PetscInt       i,j,mbs=A->rmap->n,*bi=b->i,*bj=b->j,*bdiag=b->diag,*bjtmp;
1297545dd064SHong Zhang   PetscInt       *ai=a->i,*aj=a->j,*ajtmp;
1298d595f711SHong Zhang   PetscInt       k,jmin,jmax,*c2r,*il,col,nexti,ili,nz;
1299d595f711SHong Zhang   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
1300d90ac83dSHong Zhang   FactorShiftCtx sctx;
1301d595f711SHong Zhang   PetscReal      rs;
1302d595f711SHong Zhang   MatScalar      d,*v;
1303d595f711SHong Zhang 
1304d595f711SHong Zhang   PetscFunctionBegin;
1305*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc3(mbs,&rtmp,mbs,&il,mbs,&c2r));
1306545dd064SHong Zhang 
1307d595f711SHong Zhang   /* MatPivotSetUp(): initialize shift context sctx */
1308*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMemzero(&sctx,sizeof(FactorShiftCtx)));
1309d595f711SHong Zhang 
1310f4db908eSBarry Smith   if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) { /* set sctx.shift_top=max{rs} */
1311d595f711SHong Zhang     sctx.shift_top = info->zeropivot;
131226fbe8dcSKarl Rupp 
1313*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscArrayzero(rtmp,mbs));
131426fbe8dcSKarl Rupp 
1315d595f711SHong Zhang     for (i=0; i<mbs; i++) {
1316d595f711SHong Zhang       /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
1317d595f711SHong Zhang       d        = (aa)[a->diag[i]];
1318545dd064SHong Zhang       rtmp[i] += -PetscRealPart(d);  /* diagonal entry */
1319545dd064SHong Zhang       ajtmp    = aj + ai[i] + 1;     /* exclude diagonal */
1320545dd064SHong Zhang       v        = aa + ai[i] + 1;
1321545dd064SHong Zhang       nz       = ai[i+1] - ai[i] - 1;
1322545dd064SHong Zhang       for (j=0; j<nz; j++) {
1323545dd064SHong Zhang         rtmp[i]        += PetscAbsScalar(v[j]);
1324545dd064SHong Zhang         rtmp[ajtmp[j]] += PetscAbsScalar(v[j]);
1325545dd064SHong Zhang       }
13260838c725SBarry Smith       if (PetscRealPart(rtmp[i]) > sctx.shift_top) sctx.shift_top = PetscRealPart(rtmp[i]);
1327d595f711SHong Zhang     }
1328d595f711SHong Zhang     sctx.shift_top *= 1.1;
1329d595f711SHong Zhang     sctx.nshift_max = 5;
1330d595f711SHong Zhang     sctx.shift_lo   = 0.;
1331d595f711SHong Zhang     sctx.shift_hi   = 1.;
1332d595f711SHong Zhang   }
1333d595f711SHong Zhang 
1334d595f711SHong Zhang   /* allocate working arrays
1335d595f711SHong Zhang      c2r: linked list, keep track of pivot rows for a given column. c2r[col]: head of the list for a given col
1336d595f711SHong 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
1337d595f711SHong Zhang   */
1338d595f711SHong Zhang   do {
133907b50cabSHong Zhang     sctx.newshift = PETSC_FALSE;
1340d595f711SHong Zhang 
1341d595f711SHong Zhang     for (i=0; i<mbs; i++) c2r[i] = mbs;
13429e95ef84SSatish Balay     if (mbs) il[0] = 0;
1343d595f711SHong Zhang 
1344d595f711SHong Zhang     for (k = 0; k<mbs; k++) {
1345d595f711SHong Zhang       /* zero rtmp */
1346d595f711SHong Zhang       nz    = bi[k+1] - bi[k];
1347d595f711SHong Zhang       bjtmp = bj + bi[k];
13487b056e98SHong Zhang       for (j=0; j<nz; j++) rtmp[bjtmp[j]] = 0.0;
1349d595f711SHong Zhang 
1350d595f711SHong Zhang       /* load in initial unfactored row */
1351d595f711SHong Zhang       bval = ba + bi[k];
1352d595f711SHong Zhang       jmin = ai[k]; jmax = ai[k+1];
1353d595f711SHong Zhang       for (j = jmin; j < jmax; j++) {
1354d595f711SHong Zhang         col       = aj[j];
1355d595f711SHong Zhang         rtmp[col] = aa[j];
1356d595f711SHong Zhang         *bval++   = 0.0; /* for in-place factorization */
1357d595f711SHong Zhang       }
1358d595f711SHong Zhang       /* shift the diagonal of the matrix: ZeropivotApply() */
1359d595f711SHong Zhang       rtmp[k] += sctx.shift_amount;  /* shift the diagonal of the matrix */
1360d595f711SHong Zhang 
1361d595f711SHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
1362d595f711SHong Zhang       dk = rtmp[k];
1363d595f711SHong Zhang       i  = c2r[k]; /* first row to be added to k_th row  */
1364d595f711SHong Zhang 
1365d595f711SHong Zhang       while (i < k) {
1366d595f711SHong Zhang         nexti = c2r[i]; /* next row to be added to k_th row */
1367d595f711SHong Zhang 
1368d595f711SHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
1369d595f711SHong Zhang         ili     = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
1370d595f711SHong Zhang         uikdi   = -ba[ili]*ba[bdiag[i]]; /* diagonal(k) */
1371d595f711SHong Zhang         dk     += uikdi*ba[ili]; /* update diag[k] */
1372d595f711SHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
1373d595f711SHong Zhang 
1374d595f711SHong Zhang         /* add multiple of row i to k-th row */
1375d595f711SHong Zhang         jmin = ili + 1; jmax = bi[i+1];
1376d595f711SHong Zhang         if (jmin < jmax) {
1377d595f711SHong Zhang           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
1378d595f711SHong Zhang           /* update il and c2r for row i */
1379d595f711SHong Zhang           il[i] = jmin;
1380d595f711SHong Zhang           j     = bj[jmin]; c2r[i] = c2r[j]; c2r[j] = i;
1381d595f711SHong Zhang         }
1382d595f711SHong Zhang         i = nexti;
1383d595f711SHong Zhang       }
1384d595f711SHong Zhang 
1385d595f711SHong Zhang       /* copy data into U(k,:) */
1386d595f711SHong Zhang       rs   = 0.0;
1387d595f711SHong Zhang       jmin = bi[k]; jmax = bi[k+1]-1;
1388d595f711SHong Zhang       if (jmin < jmax) {
1389d595f711SHong Zhang         for (j=jmin; j<jmax; j++) {
1390d595f711SHong Zhang           col = bj[j]; ba[j] = rtmp[col]; rs += PetscAbsScalar(ba[j]);
1391d595f711SHong Zhang         }
1392d595f711SHong Zhang         /* add the k-th row into il and c2r */
1393d595f711SHong Zhang         il[k] = jmin;
1394d595f711SHong Zhang         i     = bj[jmin]; c2r[k] = c2r[i]; c2r[i] = k;
1395d595f711SHong Zhang       }
1396d595f711SHong Zhang 
1397d595f711SHong Zhang       sctx.rs = rs;
1398d595f711SHong Zhang       sctx.pv = dk;
1399*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatPivotCheck(B,A,info,&sctx,k));
140007b50cabSHong Zhang       if (sctx.newshift) break;
1401d595f711SHong Zhang       dk = sctx.pv;
1402d595f711SHong Zhang 
1403d595f711SHong Zhang       ba[bdiag[k]] = 1.0/dk; /* U(k,k) */
1404d595f711SHong Zhang     }
140507b50cabSHong Zhang   } while (sctx.newshift);
1406d595f711SHong Zhang 
1407*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree3(rtmp,il,c2r));
1408d595f711SHong Zhang 
1409d595f711SHong Zhang   B->ops->solve          = MatSolve_SeqSBAIJ_1_NaturalOrdering;
14109dcdb97aSHong Zhang   B->ops->solves         = MatSolves_SeqSBAIJ_1;
1411d595f711SHong Zhang   B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1412910cf402Sprj-   B->ops->matsolve       = MatMatSolve_SeqSBAIJ_1_NaturalOrdering;
141380d3d5a7SHong Zhang   B->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering;
141480d3d5a7SHong Zhang   B->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering;
1415d595f711SHong Zhang 
14167b056e98SHong Zhang   B->assembled    = PETSC_TRUE;
14177b056e98SHong Zhang   B->preallocated = PETSC_TRUE;
141826fbe8dcSKarl Rupp 
1419*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops(B->rmap->n));
1420d595f711SHong Zhang 
1421d595f711SHong Zhang   /* MatPivotView() */
1422d595f711SHong Zhang   if (sctx.nshift) {
1423f4db908eSBarry Smith     if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
1424*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscInfo(A,"number of shift_pd tries %" PetscInt_FMT ", 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));
1425f4db908eSBarry Smith     } else if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
1426*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscInfo(A,"number of shift_nz tries %" PetscInt_FMT ", shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount));
1427f4db908eSBarry Smith     } else if (info->shifttype == (PetscReal)MAT_SHIFT_INBLOCKS) {
1428*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscInfo(A,"number of shift_inblocks applied %" PetscInt_FMT ", each shift_amount %g\n",sctx.nshift,(double)info->shiftamount));
1429d595f711SHong Zhang     }
1430d595f711SHong Zhang   }
1431d595f711SHong Zhang   PetscFunctionReturn(0);
1432d595f711SHong Zhang }
1433d595f711SHong Zhang 
1434d595f711SHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering_inplace(Mat C,Mat A,const MatFactorInfo *info)
1435671cb588SHong Zhang {
1436671cb588SHong Zhang   Mat_SeqSBAIJ   *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ*)C->data;
143713f74950SBarry Smith   PetscInt       i,j,mbs = a->mbs;
143813f74950SBarry Smith   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
14393cea4cbeSHong Zhang   PetscInt       k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz;
1440653a6975SHong Zhang   MatScalar      *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval;
144114d2772eSHong Zhang   PetscReal      rs;
14420e95ead3SHong Zhang   FactorShiftCtx sctx;
1443653a6975SHong Zhang 
1444653a6975SHong Zhang   PetscFunctionBegin;
14450e95ead3SHong Zhang   /* MatPivotSetUp(): initialize shift context sctx */
1446*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMemzero(&sctx,sizeof(FactorShiftCtx)));
14470e95ead3SHong Zhang 
1448653a6975SHong Zhang   /* initialization */
1449653a6975SHong Zhang   /* il and jl record the first nonzero element in each row of the accessing
1450653a6975SHong Zhang      window U(0:k, k:mbs-1).
1451653a6975SHong Zhang      jl:    list of rows to be added to uneliminated rows
1452653a6975SHong Zhang             i>= k: jl(i) is the first row to be added to row i
1453653a6975SHong Zhang             i<  k: jl(i) is the row following row i in some list of rows
1454653a6975SHong Zhang             jl(i) = mbs indicates the end of a list
1455653a6975SHong Zhang      il(i): points to the first nonzero element in U(i,k:mbs-1)
1456653a6975SHong Zhang   */
1457*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(mbs,&rtmp));
1458*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc2(mbs,&il,mbs,&jl));
1459b00f7748SHong Zhang 
1460b00f7748SHong Zhang   do {
146107b50cabSHong Zhang     sctx.newshift = PETSC_FALSE;
14626df5ee2eSHong Zhang     il[0] = 0;
1463653a6975SHong Zhang     for (i=0; i<mbs; i++) {
14646df5ee2eSHong Zhang       rtmp[i] = 0.0; jl[i] = mbs;
1465653a6975SHong Zhang     }
1466653a6975SHong Zhang 
1467997a0949SHong Zhang     for (k = 0; k<mbs; k++) {
1468653a6975SHong Zhang       /*initialize k-th row with elements nonzero in row perm(k) of A */
1469653a6975SHong Zhang       nz   = ai[k+1] - ai[k];
1470653a6975SHong Zhang       acol = aj + ai[k];
1471653a6975SHong Zhang       aval = aa + ai[k];
1472653a6975SHong Zhang       bval = ba + bi[k];
1473653a6975SHong Zhang       while (nz--) {
1474653a6975SHong Zhang         rtmp[*acol++] = *aval++;
1475653a6975SHong Zhang         *bval++       = 0.0; /* for in-place factorization */
1476653a6975SHong Zhang       }
14773cea4cbeSHong Zhang 
14783cea4cbeSHong Zhang       /* shift the diagonal of the matrix */
14793cea4cbeSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
1480653a6975SHong Zhang 
1481653a6975SHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
1482653a6975SHong Zhang       dk = rtmp[k];
1483653a6975SHong Zhang       i  = jl[k]; /* first row to be added to k_th row  */
1484653a6975SHong Zhang 
1485653a6975SHong Zhang       while (i < k) {
1486653a6975SHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
1487653a6975SHong Zhang         /* compute multiplier, update D(k) and U(i,k) */
1488653a6975SHong Zhang         ili     = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
1489653a6975SHong Zhang         uikdi   = -ba[ili]*ba[bi[i]];
1490653a6975SHong Zhang         dk     += uikdi*ba[ili];
1491653a6975SHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
1492653a6975SHong Zhang 
1493653a6975SHong Zhang         /* add multiple of row i to k-th row ... */
1494653a6975SHong Zhang         jmin = ili + 1;
1495653a6975SHong Zhang         nz   = bi[i+1] - jmin;
1496653a6975SHong Zhang         if (nz > 0) {
1497653a6975SHong Zhang           bcol = bj + jmin;
1498653a6975SHong Zhang           bval = ba + jmin;
1499*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscLogFlops(2.0*nz));
1500653a6975SHong Zhang           while (nz--) rtmp[*bcol++] += uikdi*(*bval++);
1501187a9f4bSHong Zhang 
1502997a0949SHong Zhang           /* update il and jl for i-th row */
1503997a0949SHong Zhang           il[i] = jmin;
1504997a0949SHong Zhang           j     = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
1505653a6975SHong Zhang         }
1506653a6975SHong Zhang         i = nexti;
1507653a6975SHong Zhang       }
1508653a6975SHong Zhang 
15093cea4cbeSHong Zhang       /* shift the diagonals when zero pivot is detected */
15103cea4cbeSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
15113cea4cbeSHong Zhang       rs   = 0.0;
151221c26570Svictorle       jmin = bi[k]+1;
151321c26570Svictorle       nz   = bi[k+1] - jmin;
151421c26570Svictorle       if (nz) {
151521c26570Svictorle         bcol = bj + jmin;
151621c26570Svictorle         while (nz--) {
151789f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
151889f845c8SHong Zhang           bcol++;
151921c26570Svictorle         }
152021c26570Svictorle       }
15213cea4cbeSHong Zhang 
15223cea4cbeSHong Zhang       sctx.rs = rs;
15233cea4cbeSHong Zhang       sctx.pv = dk;
1524*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatPivotCheck(C,A,info,&sctx,k));
152507b50cabSHong Zhang       if (sctx.newshift) break;    /* sctx.shift_amount is updated */
15260e95ead3SHong Zhang       dk = sctx.pv;
1527653a6975SHong Zhang 
1528997a0949SHong Zhang       /* copy data into U(k,:) */
1529653a6975SHong Zhang       ba[bi[k]] = 1.0/dk;
1530653a6975SHong Zhang       jmin      = bi[k]+1;
1531653a6975SHong Zhang       nz        = bi[k+1] - jmin;
1532653a6975SHong Zhang       if (nz) {
1533653a6975SHong Zhang         bcol = bj + jmin;
1534653a6975SHong Zhang         bval = ba + jmin;
1535653a6975SHong Zhang         while (nz--) {
1536653a6975SHong Zhang           *bval++       = rtmp[*bcol];
1537653a6975SHong Zhang           rtmp[*bcol++] = 0.0;
1538653a6975SHong Zhang         }
1539997a0949SHong Zhang         /* add k-th row into il and jl */
1540653a6975SHong Zhang         il[k] = jmin;
1541997a0949SHong Zhang         i     = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
1542653a6975SHong Zhang       }
1543b00f7748SHong Zhang     } /* end of for (k = 0; k<mbs; k++) */
154407b50cabSHong Zhang   } while (sctx.newshift);
1545*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(rtmp));
1546*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree2(il,jl));
1547653a6975SHong Zhang 
15480a3c351aSHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
15494f79d315SHong Zhang   C->ops->solves         = MatSolves_SeqSBAIJ_1_inplace;
15500a3c351aSHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
15510a3c351aSHong Zhang   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
15520a3c351aSHong Zhang   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
1553db4efbfdSBarry Smith 
1554653a6975SHong Zhang   C->assembled    = PETSC_TRUE;
1555653a6975SHong Zhang   C->preallocated = PETSC_TRUE;
155626fbe8dcSKarl Rupp 
1557*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops(C->rmap->N));
15583cea4cbeSHong Zhang   if (sctx.nshift) {
1559f4db908eSBarry Smith     if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
1560*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscInfo(A,"number of shiftnz tries %" PetscInt_FMT ", shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount));
1561f4db908eSBarry Smith     } else if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
1562*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscInfo(A,"number of shiftpd tries %" PetscInt_FMT ", shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount));
1563b00f7748SHong Zhang     }
1564fb10cecfSBarry Smith   }
1565653a6975SHong Zhang   PetscFunctionReturn(0);
1566653a6975SHong Zhang }
1567653a6975SHong Zhang 
15680481f469SBarry Smith PetscErrorCode MatCholeskyFactor_SeqSBAIJ(Mat A,IS perm,const MatFactorInfo *info)
156949b5e25fSSatish Balay {
157049b5e25fSSatish Balay   Mat            C;
157149b5e25fSSatish Balay 
157249b5e25fSSatish Balay   PetscFunctionBegin;
1573*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetFactor(A,"petsc",MAT_FACTOR_CHOLESKY,&C));
1574*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCholeskyFactorSymbolic(C,A,perm,info));
1575*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCholeskyFactorNumeric(C,A,info));
157626fbe8dcSKarl Rupp 
1577db4efbfdSBarry Smith   A->ops->solve          = C->ops->solve;
1578db4efbfdSBarry Smith   A->ops->solvetranspose = C->ops->solvetranspose;
157926fbe8dcSKarl Rupp 
1580*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatHeaderMerge(A,&C));
158149b5e25fSSatish Balay   PetscFunctionReturn(0);
158249b5e25fSSatish Balay }
1583