xref: /petsc/src/mat/impls/sbaij/seq/sbaijfact.c (revision 96b95a6bf9ed541503cb02ce9899069044d5ddfa)
149b5e25fSSatish Balay 
2c6db04a5SJed Brown #include <../src/mat/impls/baij/seq/baij.h>
3c6db04a5SJed Brown #include <../src/mat/impls/sbaij/seq/sbaij.h>
4c6db04a5SJed Brown #include <../src/mat/blockinvert.h>
5c6db04a5SJed Brown #include <petscis.h>
649b5e25fSSatish Balay 
75f9f512dSHong Zhang /*
85f9f512dSHong Zhang   input:
9c037c3f7SHong Zhang    F -- numeric factor
105f9f512dSHong Zhang   output:
11c037c3f7SHong Zhang    nneg, nzero, npos: matrix inertia
125f9f512dSHong Zhang */
135f9f512dSHong Zhang 
145f9f512dSHong Zhang #undef __FUNCT__
155f9f512dSHong Zhang #define __FUNCT__ "MatGetInertia_SeqSBAIJ"
1613f74950SBarry Smith PetscErrorCode MatGetInertia_SeqSBAIJ(Mat F,PetscInt *nneig,PetscInt *nzero,PetscInt *npos)
175f9f512dSHong Zhang {
18638f5ce0SDinesh Kaushik   Mat_SeqSBAIJ *fact_ptr = (Mat_SeqSBAIJ*)F->data;
19dd6ea824SBarry Smith   MatScalar    *dd = fact_ptr->a;
20e70ad169SHong Zhang   PetscInt     mbs=fact_ptr->mbs,bs=F->rmap->bs,i,nneig_tmp,npos_tmp,*fi = fact_ptr->diag;
215f9f512dSHong Zhang 
225f9f512dSHong Zhang   PetscFunctionBegin;
23e32f2f54SBarry Smith   if (bs != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for bs: %D >1 yet",bs);
24eeeff2ecSHong Zhang   nneig_tmp = 0; npos_tmp = 0;
25eeeff2ecSHong Zhang   for (i=0; i<mbs; i++){
26eeeff2ecSHong Zhang     if (PetscRealPart(dd[*fi]) > 0.0){
27eeeff2ecSHong Zhang       npos_tmp++;
28eeeff2ecSHong Zhang     } else if (PetscRealPart(dd[*fi]) < 0.0){
29eeeff2ecSHong Zhang       nneig_tmp++;
305f9f512dSHong Zhang     }
31eeeff2ecSHong Zhang     fi++;
323e0d88b5SBarry Smith   }
33eeeff2ecSHong Zhang   if (nneig) *nneig = nneig_tmp;
34eeeff2ecSHong Zhang   if (npos)  *npos  = npos_tmp;
35eeeff2ecSHong Zhang   if (nzero) *nzero = mbs - nneig_tmp - npos_tmp;
36eeeff2ecSHong Zhang 
375f9f512dSHong Zhang   PetscFunctionReturn(0);
385f9f512dSHong Zhang }
395f9f512dSHong Zhang 
405f9f512dSHong Zhang /*
415f9f512dSHong Zhang   Symbolic U^T*D*U factorization for SBAIJ format. Modified from SSF of YSMP.
4210c27e3dSHong Zhang   Use Modified Sparse Row (MSR) storage for u and ju. See page 85, "Iterative Methods ..." by Saad.
435f9f512dSHong Zhang */
4410c27e3dSHong Zhang #undef __FUNCT__
4510c27e3dSHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqSBAIJ_MSR"
460481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ_MSR(Mat F,Mat A,IS perm,const MatFactorInfo *info)
4710c27e3dSHong Zhang {
4810c27e3dSHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b;
4910c27e3dSHong Zhang   PetscErrorCode ierr;
505d0c19d7SBarry Smith   const PetscInt *rip,*ai,*aj;
515d0c19d7SBarry Smith   PetscInt       i,mbs = a->mbs,*jutmp,bs = A->rmap->bs,bs2=a->bs2;
5210c27e3dSHong Zhang   PetscInt       m,reallocs = 0,prow;
5310c27e3dSHong Zhang   PetscInt       *jl,*q,jmin,jmax,juidx,nzk,qm,*iu,*ju,k,j,vj,umax,maxadd;
5410c27e3dSHong Zhang   PetscReal      f = info->fill;
55ace3abfcSBarry Smith   PetscBool      perm_identity;
5610c27e3dSHong Zhang 
5710c27e3dSHong Zhang   PetscFunctionBegin;
5810c27e3dSHong Zhang   /* check whether perm is the identity mapping */
5910c27e3dSHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
6010c27e3dSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
6110c27e3dSHong Zhang 
6210c27e3dSHong Zhang   if (perm_identity){ /* without permutation */
6310c27e3dSHong Zhang     a->permute = PETSC_FALSE;
6410c27e3dSHong Zhang     ai = a->i; aj = a->j;
6510c27e3dSHong Zhang   } else {            /* non-trivial permutation */
6610c27e3dSHong Zhang     a->permute = PETSC_TRUE;
6710c27e3dSHong Zhang     ierr = MatReorderingSeqSBAIJ(A,perm);CHKERRQ(ierr);
6810c27e3dSHong Zhang     ai = a->inew; aj = a->jnew;
6910c27e3dSHong Zhang   }
7010c27e3dSHong Zhang 
7110c27e3dSHong Zhang   /* initialization */
7210c27e3dSHong Zhang   ierr  = PetscMalloc((mbs+1)*sizeof(PetscInt),&iu);CHKERRQ(ierr);
7310c27e3dSHong Zhang   umax  = (PetscInt)(f*ai[mbs] + 1); umax += mbs + 1;
7410c27e3dSHong Zhang   ierr  = PetscMalloc(umax*sizeof(PetscInt),&ju);CHKERRQ(ierr);
7510c27e3dSHong Zhang   iu[0] = mbs+1;
7610c27e3dSHong Zhang   juidx = mbs + 1; /* index for ju */
77d8c74875SBarry Smith   /* jl linked list for pivot row -- linked list for col index */
78d8c74875SBarry Smith   ierr  = PetscMalloc2(mbs,PetscInt,&jl,mbs,PetscInt,&q);CHKERRQ(ierr);
7910c27e3dSHong Zhang   for (i=0; i<mbs; i++){
8010c27e3dSHong Zhang     jl[i] = mbs;
8110c27e3dSHong Zhang     q[i] = 0;
8210c27e3dSHong Zhang   }
8310c27e3dSHong Zhang 
8410c27e3dSHong Zhang   /* for each row k */
8510c27e3dSHong Zhang   for (k=0; k<mbs; k++){
8610c27e3dSHong Zhang     for (i=0; i<mbs; i++) q[i] = 0;  /* to be removed! */
8710c27e3dSHong Zhang     nzk  = 0; /* num. of nz blocks in k-th block row with diagonal block excluded */
8810c27e3dSHong Zhang     q[k] = mbs;
8910c27e3dSHong Zhang     /* initialize nonzero structure of k-th row to row rip[k] of A */
9010c27e3dSHong Zhang     jmin = ai[rip[k]] +1; /* exclude diag[k] */
9110c27e3dSHong Zhang     jmax = ai[rip[k]+1];
9210c27e3dSHong Zhang     for (j=jmin; j<jmax; j++){
9310c27e3dSHong Zhang       vj = rip[aj[j]]; /* col. value */
9410c27e3dSHong Zhang       if(vj > k){
9510c27e3dSHong Zhang         qm = k;
9610c27e3dSHong Zhang         do {
9710c27e3dSHong Zhang           m  = qm; qm = q[m];
9810c27e3dSHong Zhang         } while(qm < vj);
9910c27e3dSHong Zhang         if (qm == vj) {
100e32f2f54SBarry Smith           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Duplicate entry in A\n");
10110c27e3dSHong Zhang         }
10210c27e3dSHong Zhang         nzk++;
10310c27e3dSHong Zhang         q[m]  = vj;
10410c27e3dSHong Zhang         q[vj] = qm;
10510c27e3dSHong Zhang       } /* if(vj > k) */
10610c27e3dSHong Zhang     } /* for (j=jmin; j<jmax; j++) */
10710c27e3dSHong Zhang 
10810c27e3dSHong Zhang     /* modify nonzero structure of k-th row by computing fill-in
10910c27e3dSHong Zhang        for each row i to be merged in */
11010c27e3dSHong Zhang     prow = k;
11110c27e3dSHong Zhang     prow = jl[prow]; /* next pivot row (== mbs for symbolic factorization) */
11210c27e3dSHong Zhang 
11310c27e3dSHong Zhang     while (prow < k){
11410c27e3dSHong Zhang       /* merge row prow into k-th row */
11510c27e3dSHong Zhang       jmin = iu[prow] + 1; jmax = iu[prow+1];
11610c27e3dSHong Zhang       qm = k;
11710c27e3dSHong Zhang       for (j=jmin; j<jmax; j++){
11810c27e3dSHong Zhang         vj = ju[j];
11910c27e3dSHong Zhang         do {
12010c27e3dSHong Zhang           m = qm; qm = q[m];
12110c27e3dSHong Zhang         } while (qm < vj);
12210c27e3dSHong Zhang         if (qm != vj){
12310c27e3dSHong Zhang          nzk++; q[m] = vj; q[vj] = qm; qm = vj;
12410c27e3dSHong Zhang         }
12510c27e3dSHong Zhang       }
12610c27e3dSHong Zhang       prow = jl[prow]; /* next pivot row */
12710c27e3dSHong Zhang     }
12810c27e3dSHong Zhang 
12910c27e3dSHong Zhang     /* add k to row list for first nonzero element in k-th row */
13010c27e3dSHong Zhang     if (nzk > 0){
13110c27e3dSHong Zhang       i = q[k]; /* col value of first nonzero element in U(k, k+1:mbs-1) */
13210c27e3dSHong Zhang       jl[k] = jl[i]; jl[i] = k;
13310c27e3dSHong Zhang     }
13410c27e3dSHong Zhang     iu[k+1] = iu[k] + nzk;
13510c27e3dSHong Zhang 
13610c27e3dSHong Zhang     /* allocate more space to ju if needed */
13710c27e3dSHong Zhang     if (iu[k+1] > umax) {
13810c27e3dSHong Zhang       /* estimate how much additional space we will need */
13910c27e3dSHong Zhang       /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */
14010c27e3dSHong Zhang       /* just double the memory each time */
14110c27e3dSHong Zhang       maxadd = umax;
14210c27e3dSHong Zhang       if (maxadd < nzk) maxadd = (mbs-k)*(nzk+1)/2;
14310c27e3dSHong Zhang       umax += maxadd;
14410c27e3dSHong Zhang 
14510c27e3dSHong Zhang       /* allocate a longer ju */
14610c27e3dSHong Zhang       ierr = PetscMalloc(umax*sizeof(PetscInt),&jutmp);CHKERRQ(ierr);
14710c27e3dSHong Zhang       ierr = PetscMemcpy(jutmp,ju,iu[k]*sizeof(PetscInt));CHKERRQ(ierr);
14810c27e3dSHong Zhang       ierr = PetscFree(ju);CHKERRQ(ierr);
14910c27e3dSHong Zhang       ju   = jutmp;
15010c27e3dSHong Zhang       reallocs++; /* count how many times we realloc */
15110c27e3dSHong Zhang     }
15210c27e3dSHong Zhang 
15310c27e3dSHong Zhang     /* save nonzero structure of k-th row in ju */
15410c27e3dSHong Zhang     i=k;
15510c27e3dSHong Zhang     while (nzk --) {
15610c27e3dSHong Zhang       i           = q[i];
15710c27e3dSHong Zhang       ju[juidx++] = i;
15810c27e3dSHong Zhang     }
15910c27e3dSHong Zhang   }
16010c27e3dSHong Zhang 
1616cf91177SBarry Smith #if defined(PETSC_USE_INFO)
16210c27e3dSHong Zhang   if (ai[mbs] != 0) {
16310c27e3dSHong Zhang     PetscReal af = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]);
164ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr);
165ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
166ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr);
167ae15b995SBarry Smith     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
16810c27e3dSHong Zhang   } else {
169ae15b995SBarry Smith     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
17010c27e3dSHong Zhang   }
17163ba0a88SBarry Smith #endif
17210c27e3dSHong Zhang 
17310c27e3dSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
174d8c74875SBarry Smith   ierr = PetscFree2(jl,q);CHKERRQ(ierr);
17510c27e3dSHong Zhang 
17610c27e3dSHong Zhang   /* put together the new matrix */
177719d5645SBarry Smith   ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(F,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
17810c27e3dSHong Zhang 
179719d5645SBarry Smith   /* ierr = PetscLogObjectParent(B,iperm);CHKERRQ(ierr); */
180719d5645SBarry Smith   b = (Mat_SeqSBAIJ*)(F)->data;
18110c27e3dSHong Zhang   b->singlemalloc = PETSC_FALSE;
182e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
183e6b907acSBarry Smith   b->free_ij       = PETSC_TRUE;
18410c27e3dSHong Zhang   ierr = PetscMalloc((iu[mbs]+1)*sizeof(MatScalar)*bs2,&b->a);CHKERRQ(ierr);
18510c27e3dSHong Zhang   b->j    = ju;
18610c27e3dSHong Zhang   b->i    = iu;
18710c27e3dSHong Zhang   b->diag = 0;
18810c27e3dSHong Zhang   b->ilen = 0;
18910c27e3dSHong Zhang   b->imax = 0;
19010c27e3dSHong Zhang   b->row  = perm;
19110c27e3dSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
19210c27e3dSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
19310c27e3dSHong Zhang   b->icol = perm;
19410c27e3dSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
19510c27e3dSHong Zhang   ierr    = PetscMalloc((bs*mbs+bs)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
19610c27e3dSHong Zhang   /* In b structure:  Free imax, ilen, old a, old j.
19710c27e3dSHong Zhang      Allocate idnew, solve_work, new a, new j */
198719d5645SBarry Smith   ierr = PetscLogObjectMemory(F,(iu[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
19910c27e3dSHong Zhang   b->maxnz = b->nz = iu[mbs];
20010c27e3dSHong Zhang 
201719d5645SBarry Smith   (F)->info.factor_mallocs    = reallocs;
202719d5645SBarry Smith   (F)->info.fill_ratio_given  = f;
20310c27e3dSHong Zhang   if (ai[mbs] != 0) {
204719d5645SBarry Smith     (F)->info.fill_ratio_needed = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]);
20510c27e3dSHong Zhang   } else {
206719d5645SBarry Smith     (F)->info.fill_ratio_needed = 0.0;
20710c27e3dSHong Zhang   }
208d595f711SHong Zhang   ierr = MatSeqSBAIJSetNumericFactorization_inplace(F,perm_identity);CHKERRQ(ierr);
20910c27e3dSHong Zhang   PetscFunctionReturn(0);
21010c27e3dSHong Zhang }
21110c27e3dSHong Zhang /*
21210c27e3dSHong Zhang     Symbolic U^T*D*U factorization for SBAIJ format.
21380d3d5a7SHong Zhang     See MatICCFactorSymbolic_SeqAIJ() for description of its data structure.
21410c27e3dSHong Zhang */
215c6db04a5SJed Brown #include <petscbt.h>
216c6db04a5SJed Brown #include <../src/mat/utils/freespace.h>
2174a2ae208SSatish Balay #undef __FUNCT__
2184a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqSBAIJ"
2190481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
22049b5e25fSSatish Balay {
22198a8e78dSHong Zhang   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data;
22298a8e78dSHong Zhang   Mat_SeqSBAIJ       *b;
2236849ba73SBarry Smith   PetscErrorCode     ierr;
224ace3abfcSBarry Smith   PetscBool          perm_identity,missing;
22598a8e78dSHong Zhang   PetscReal          fill = info->fill;
2267b056e98SHong Zhang   const PetscInt     *rip,*ai=a->i,*aj=a->j;
2275d0c19d7SBarry Smith   PetscInt           i,mbs=a->mbs,bs=A->rmap->bs,reallocs=0,prow,d;
22898a8e78dSHong Zhang   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
229c6d0d4f0SHong Zhang   PetscInt           nlnk,*lnk,ncols,*cols,*uj,**ui_ptr,*uj_ptr,*udiag;
230a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
23198a8e78dSHong Zhang   PetscBT            lnkbt;
232d595f711SHong Zhang 
233d595f711SHong Zhang   PetscFunctionBegin;
2346d07c2b1SHong Zhang   if (bs > 1){
2356d07c2b1SHong Zhang     ierr = MatCholeskyFactorSymbolic_SeqSBAIJ_inplace(fact,A,perm,info);CHKERRQ(ierr);
2366d07c2b1SHong Zhang     PetscFunctionReturn(0);
2376d07c2b1SHong Zhang   }
238e32f2f54SBarry Smith   if (A->rmap->n != A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Must be square matrix, rows %D columns %D",A->rmap->n,A->cmap->n);
239d595f711SHong Zhang   ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr);
240e32f2f54SBarry Smith   if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
241d595f711SHong Zhang 
242d595f711SHong Zhang   /* check whether perm is the identity mapping */
243d595f711SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
244e32f2f54SBarry Smith   if (!perm_identity) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix reordering is not supported for sbaij matrix. Use aij format");
245d595f711SHong Zhang   a->permute = PETSC_FALSE;
246d595f711SHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
247d595f711SHong Zhang 
248d595f711SHong Zhang   /* initialization */
249d595f711SHong Zhang   ierr  = PetscMalloc((mbs+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
250c6d0d4f0SHong Zhang   ierr  = PetscMalloc((mbs+1)*sizeof(PetscInt),&udiag);CHKERRQ(ierr);
251d595f711SHong Zhang   ui[0] = 0;
252d595f711SHong Zhang 
253d595f711SHong Zhang   /* jl: linked list for storing indices of the pivot rows
254d595f711SHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
255d595f711SHong Zhang   ierr = PetscMalloc4(mbs,PetscInt*,&ui_ptr,mbs,PetscInt,&il,mbs,PetscInt,&jl,mbs,PetscInt,&cols);CHKERRQ(ierr);
256d595f711SHong Zhang   for (i=0; i<mbs; i++){
257d595f711SHong Zhang     jl[i] = mbs; il[i] = 0;
258d595f711SHong Zhang   }
259d595f711SHong Zhang 
260d595f711SHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
261d595f711SHong Zhang   nlnk = mbs + 1;
262d595f711SHong Zhang   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
263d595f711SHong Zhang 
264d595f711SHong Zhang   /* initial FreeSpace size is fill*(ai[mbs]+1) */
265d595f711SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[mbs]+1)),&free_space);CHKERRQ(ierr);
266d595f711SHong Zhang   current_space = free_space;
267d595f711SHong Zhang 
268d595f711SHong Zhang   for (k=0; k<mbs; k++){  /* for each active row k */
269d595f711SHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
270d595f711SHong Zhang     nzk   = 0;
271c6d0d4f0SHong Zhang     ncols = ai[k+1] - ai[k];
272e32f2f54SBarry Smith     if (!ncols) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_CH_ZRPVT,"Empty row %D in matrix ",k);
273d595f711SHong Zhang     for (j=0; j<ncols; j++){
274c6d0d4f0SHong Zhang       i = *(aj + ai[k] + j);
275c6d0d4f0SHong Zhang       cols[j] = i;
276d595f711SHong Zhang     }
277d595f711SHong Zhang     ierr = PetscLLAdd(ncols,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
278d595f711SHong Zhang     nzk += nlnk;
279d595f711SHong Zhang 
280d595f711SHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
281d595f711SHong Zhang     prow = jl[k]; /* 1st pivot row */
282d595f711SHong Zhang 
283d595f711SHong Zhang     while (prow < k){
284d595f711SHong Zhang       nextprow = jl[prow];
285d595f711SHong Zhang       /* merge prow into k-th row */
286d595f711SHong Zhang       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
287d595f711SHong Zhang       jmax = ui[prow+1];
288d595f711SHong Zhang       ncols = jmax-jmin;
289d595f711SHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
290d595f711SHong Zhang       ierr = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
291d595f711SHong Zhang       nzk += nlnk;
292d595f711SHong Zhang 
293d595f711SHong Zhang       /* update il and jl for prow */
294d595f711SHong Zhang       if (jmin < jmax){
295d595f711SHong Zhang         il[prow] = jmin;
296d595f711SHong Zhang         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
297d595f711SHong Zhang       }
298d595f711SHong Zhang       prow = nextprow;
299d595f711SHong Zhang     }
300d595f711SHong Zhang 
301d595f711SHong Zhang     /* if free space is not available, make more free space */
302d595f711SHong Zhang     if (current_space->local_remaining<nzk) {
303d595f711SHong Zhang       i  = mbs - k + 1; /* num of unfactored rows */
3047b056e98SHong Zhang       i *= PetscMin(nzk, i-1); /* i*nzk, i*(i-1): estimated and max additional space needed */
305d595f711SHong Zhang       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
306d595f711SHong Zhang       reallocs++;
307d595f711SHong Zhang     }
308d595f711SHong Zhang 
309d595f711SHong Zhang     /* copy data into free space, then initialize lnk */
310d595f711SHong Zhang     ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
311d595f711SHong Zhang 
312d595f711SHong Zhang     /* add the k-th row into il and jl */
3137b056e98SHong Zhang     if (nzk > 1){
314d595f711SHong Zhang       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
315d595f711SHong Zhang       jl[k] = jl[i]; jl[i] = k;
316d595f711SHong Zhang       il[k] = ui[k] + 1;
317d595f711SHong Zhang     }
318d595f711SHong Zhang     ui_ptr[k] = current_space->array;
319d595f711SHong Zhang     current_space->array           += nzk;
320d595f711SHong Zhang     current_space->local_used      += nzk;
321d595f711SHong Zhang     current_space->local_remaining -= nzk;
322d595f711SHong Zhang 
323d595f711SHong Zhang     ui[k+1] = ui[k] + nzk;
324d595f711SHong Zhang   }
325d595f711SHong Zhang 
326d595f711SHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
327d595f711SHong Zhang   ierr = PetscFree4(ui_ptr,il,jl,cols);CHKERRQ(ierr);
328d595f711SHong Zhang 
329d595f711SHong Zhang   /* destroy list of free space and other temporary array(s) */
330d595f711SHong Zhang   ierr = PetscMalloc((ui[mbs]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
3317b056e98SHong Zhang   ierr = PetscFreeSpaceContiguous_Cholesky(&free_space,uj,mbs,ui,udiag);CHKERRQ(ierr); /* store matrix factor */
332d595f711SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
333d595f711SHong Zhang 
334d595f711SHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
335d595f711SHong Zhang   ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(fact,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
336d595f711SHong Zhang 
3377b056e98SHong Zhang   b = (Mat_SeqSBAIJ*)fact->data;
338d595f711SHong Zhang   b->singlemalloc = PETSC_FALSE;
339d595f711SHong Zhang   b->free_a       = PETSC_TRUE;
340d595f711SHong Zhang   b->free_ij      = PETSC_TRUE;
341d595f711SHong Zhang   ierr = PetscMalloc((ui[mbs]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
342d595f711SHong Zhang   b->j    = uj;
343d595f711SHong Zhang   b->i    = ui;
344c6d0d4f0SHong Zhang   b->diag = udiag;
345c6d0d4f0SHong Zhang   b->free_diag = PETSC_TRUE;
346d595f711SHong Zhang   b->ilen = 0;
347d595f711SHong Zhang   b->imax = 0;
348d595f711SHong Zhang   b->row  = perm;
349d595f711SHong Zhang   b->icol = perm;
350d595f711SHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
3517b056e98SHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
3527b056e98SHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
353d595f711SHong Zhang   ierr    = PetscMalloc((mbs+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
3547b056e98SHong Zhang   ierr    = PetscLogObjectMemory(fact,ui[mbs]*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
355d595f711SHong Zhang   b->maxnz = b->nz = ui[mbs];
356d595f711SHong Zhang 
35795b5ac22SHong Zhang   fact->info.factor_mallocs    = reallocs;
35895b5ac22SHong Zhang   fact->info.fill_ratio_given  = fill;
359d595f711SHong Zhang   if (ai[mbs] != 0) {
36095b5ac22SHong Zhang     fact->info.fill_ratio_needed = ((PetscReal)ui[mbs])/ai[mbs];
361d595f711SHong Zhang   } else {
36295b5ac22SHong Zhang     fact->info.fill_ratio_needed = 0.0;
363d595f711SHong Zhang   }
36495b5ac22SHong Zhang #if defined(PETSC_USE_INFO)
36595b5ac22SHong Zhang   if (ai[mbs] != 0) {
36695b5ac22SHong Zhang     PetscReal af = fact->info.fill_ratio_needed;
36795b5ac22SHong Zhang     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
36895b5ac22SHong Zhang     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
36995b5ac22SHong Zhang     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
37095b5ac22SHong Zhang   } else {
37195b5ac22SHong Zhang     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
37295b5ac22SHong Zhang   }
37395b5ac22SHong Zhang #endif
374c6d0d4f0SHong Zhang   fact->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering;
375d595f711SHong Zhang   PetscFunctionReturn(0);
376d595f711SHong Zhang }
377d595f711SHong Zhang 
378d595f711SHong Zhang #undef __FUNCT__
379d595f711SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqSBAIJ_inplace"
380d595f711SHong Zhang PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ_inplace(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
381d595f711SHong Zhang {
382d595f711SHong Zhang   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data;
383d595f711SHong Zhang   Mat_SeqSBAIJ       *b;
384d595f711SHong Zhang   PetscErrorCode     ierr;
385ace3abfcSBarry Smith   PetscBool          perm_identity,missing;
386d595f711SHong Zhang   PetscReal          fill = info->fill;
387d595f711SHong Zhang   const PetscInt     *rip,*ai,*aj;
388d595f711SHong Zhang   PetscInt           i,mbs=a->mbs,bs=A->rmap->bs,reallocs=0,prow,d;
389d595f711SHong Zhang   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
390d595f711SHong Zhang   PetscInt           nlnk,*lnk,ncols,*cols,*uj,**ui_ptr,*uj_ptr;
391d595f711SHong Zhang   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
392d595f711SHong Zhang   PetscBT            lnkbt;
39349b5e25fSSatish Balay 
39449b5e25fSSatish Balay   PetscFunctionBegin;
39558ebbce7SBarry Smith   ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr);
396e32f2f54SBarry Smith   if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
39758ebbce7SBarry Smith 
39810c27e3dSHong Zhang   /*
39910c27e3dSHong Zhang    This code originally uses Modified Sparse Row (MSR) storage
40010c27e3dSHong Zhang    (see page 85, "Iterative Methods ..." by Saad) for the output matrix B - bad choise!
40110c27e3dSHong Zhang    Then it is rewritten so the factor B takes seqsbaij format. However the associated
40210c27e3dSHong Zhang    MatCholeskyFactorNumeric_() have not been modified for the cases of bs>1 or !perm_identity,
40310c27e3dSHong Zhang    thus the original code in MSR format is still used for these cases.
40410c27e3dSHong Zhang    The code below should replace MatCholeskyFactorSymbolic_SeqSBAIJ_MSR() whenever
40510c27e3dSHong Zhang    MatCholeskyFactorNumeric_() is modified for using sbaij symbolic factor.
40610c27e3dSHong Zhang   */
407fff829cfSHong Zhang   if (bs > 1){
408719d5645SBarry Smith     ierr = MatCholeskyFactorSymbolic_SeqSBAIJ_MSR(fact,A,perm,info);CHKERRQ(ierr);
409fff829cfSHong Zhang     PetscFunctionReturn(0);
410fff829cfSHong Zhang   }
41110c27e3dSHong Zhang 
41298a8e78dSHong Zhang   /* check whether perm is the identity mapping */
41398a8e78dSHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
41498a8e78dSHong Zhang 
415fff829cfSHong Zhang   if (perm_identity){
416fff829cfSHong Zhang     a->permute = PETSC_FALSE;
417fff829cfSHong Zhang     ai = a->i; aj = a->j;
418fff829cfSHong Zhang   } else {
419e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix reordering is not supported for sbaij matrix. Use aij format");
420d33dd558SHong Zhang     /* There are bugs for reordeing. Needs further work.
421d33dd558SHong Zhang        MatReordering for sbaij cannot be efficient. User should use aij formt! */
422fff829cfSHong Zhang     a->permute = PETSC_TRUE;
423fff829cfSHong Zhang     ierr = MatReorderingSeqSBAIJ(A,perm);CHKERRQ(ierr);
424fff829cfSHong Zhang     ai = a->inew; aj = a->jnew;
425fff829cfSHong Zhang   }
426fff829cfSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
427fff829cfSHong Zhang 
428fff829cfSHong Zhang   /* initialization */
4297625dc9aSHong Zhang   ierr  = PetscMalloc((mbs+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
43098a8e78dSHong Zhang   ui[0] = 0;
43198a8e78dSHong Zhang 
43298a8e78dSHong Zhang   /* jl: linked list for storing indices of the pivot rows
4337625dc9aSHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
434d8c74875SBarry Smith   ierr = PetscMalloc4(mbs,PetscInt*,&ui_ptr,mbs,PetscInt,&il,mbs,PetscInt,&jl,mbs,PetscInt,&cols);CHKERRQ(ierr);
4357625dc9aSHong Zhang   for (i=0; i<mbs; i++){
4367625dc9aSHong Zhang     jl[i] = mbs; il[i] = 0;
437fff829cfSHong Zhang   }
438fff829cfSHong Zhang 
43998a8e78dSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
4407625dc9aSHong Zhang   nlnk = mbs + 1;
4417625dc9aSHong Zhang   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
442fff829cfSHong Zhang 
4437625dc9aSHong Zhang   /* initial FreeSpace size is fill*(ai[mbs]+1) */
444a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[mbs]+1)),&free_space);CHKERRQ(ierr);
44598a8e78dSHong Zhang   current_space = free_space;
44698a8e78dSHong Zhang 
4477625dc9aSHong Zhang   for (k=0; k<mbs; k++){  /* for each active row k */
44898a8e78dSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
44998a8e78dSHong Zhang     nzk   = 0;
45098a8e78dSHong Zhang     ncols = ai[rip[k]+1] - ai[rip[k]];
4517625dc9aSHong Zhang     for (j=0; j<ncols; j++){
4527625dc9aSHong Zhang       i = *(aj + ai[rip[k]] + j);
4537625dc9aSHong Zhang       cols[j] = rip[i];
4547625dc9aSHong Zhang     }
4557625dc9aSHong Zhang     ierr = PetscLLAdd(ncols,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
45698a8e78dSHong Zhang     nzk += nlnk;
45798a8e78dSHong Zhang 
45898a8e78dSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
45998a8e78dSHong Zhang     prow = jl[k]; /* 1st pivot row */
460fff829cfSHong Zhang 
461fff829cfSHong Zhang     while (prow < k){
462fff829cfSHong Zhang       nextprow = jl[prow];
46398a8e78dSHong Zhang       /* merge prow into k-th row */
4647625dc9aSHong Zhang       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
46598a8e78dSHong Zhang       jmax = ui[prow+1];
46698a8e78dSHong Zhang       ncols = jmax-jmin;
4677625dc9aSHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
4685a8e39fbSHong Zhang       ierr = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
46998a8e78dSHong Zhang       nzk += nlnk;
470fff829cfSHong Zhang 
47198a8e78dSHong Zhang       /* update il and jl for prow */
472fff829cfSHong Zhang       if (jmin < jmax){
473fff829cfSHong Zhang         il[prow] = jmin;
4747625dc9aSHong Zhang         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
475fff829cfSHong Zhang       }
476fff829cfSHong Zhang       prow = nextprow;
477fff829cfSHong Zhang     }
478fff829cfSHong Zhang 
47998a8e78dSHong Zhang     /* if free space is not available, make more free space */
48098a8e78dSHong Zhang     if (current_space->local_remaining<nzk) {
4817625dc9aSHong Zhang       i = mbs - k + 1; /* num of unfactored rows */
48298a8e78dSHong Zhang       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
483a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
48498a8e78dSHong Zhang       reallocs++;
48598a8e78dSHong Zhang     }
48698a8e78dSHong Zhang 
48798a8e78dSHong Zhang     /* copy data into free space, then initialize lnk */
4887625dc9aSHong Zhang     ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
48998a8e78dSHong Zhang 
49098a8e78dSHong Zhang     /* add the k-th row into il and jl */
49198a8e78dSHong Zhang     if (nzk-1 > 0){
4927625dc9aSHong Zhang       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
493fff829cfSHong Zhang       jl[k] = jl[i]; jl[i] = k;
49498a8e78dSHong Zhang       il[k] = ui[k] + 1;
495fff829cfSHong Zhang     }
4967625dc9aSHong Zhang     ui_ptr[k] = current_space->array;
49798a8e78dSHong Zhang     current_space->array           += nzk;
49898a8e78dSHong Zhang     current_space->local_used      += nzk;
49998a8e78dSHong Zhang     current_space->local_remaining -= nzk;
500fff829cfSHong Zhang 
50198a8e78dSHong Zhang     ui[k+1] = ui[k] + nzk;
502fff829cfSHong Zhang   }
503fff829cfSHong Zhang 
504fff829cfSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
505d8c74875SBarry Smith   ierr = PetscFree4(ui_ptr,il,jl,cols);CHKERRQ(ierr);
506fff829cfSHong Zhang 
50798a8e78dSHong Zhang   /* destroy list of free space and other temporary array(s) */
5087625dc9aSHong Zhang   ierr = PetscMalloc((ui[mbs]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
509a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
51098a8e78dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
511fff829cfSHong Zhang 
51298a8e78dSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
513719d5645SBarry Smith   ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(fact,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
51498a8e78dSHong Zhang 
51595b5ac22SHong Zhang   b = (Mat_SeqSBAIJ*)fact->data;
516fff829cfSHong Zhang   b->singlemalloc = PETSC_FALSE;
517e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
518e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
5197625dc9aSHong Zhang   ierr = PetscMalloc((ui[mbs]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
52098a8e78dSHong Zhang   b->j    = uj;
52198a8e78dSHong Zhang   b->i    = ui;
522fff829cfSHong Zhang   b->diag = 0;
523fff829cfSHong Zhang   b->ilen = 0;
524fff829cfSHong Zhang   b->imax = 0;
525fff829cfSHong Zhang   b->row  = perm;
526fff829cfSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
527fff829cfSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
528fff829cfSHong Zhang   b->icol = perm;
529fff829cfSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
5307625dc9aSHong Zhang   ierr    = PetscMalloc((mbs+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
531719d5645SBarry Smith   ierr    = PetscLogObjectMemory(fact,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
5327625dc9aSHong Zhang   b->maxnz = b->nz = ui[mbs];
533fff829cfSHong Zhang 
53495b5ac22SHong Zhang   fact->info.factor_mallocs    = reallocs;
53595b5ac22SHong Zhang   fact->info.fill_ratio_given  = fill;
5367625dc9aSHong Zhang   if (ai[mbs] != 0) {
53795b5ac22SHong Zhang     fact->info.fill_ratio_needed = ((PetscReal)ui[mbs])/ai[mbs];
538fff829cfSHong Zhang   } else {
53995b5ac22SHong Zhang     fact->info.fill_ratio_needed = 0.0;
540fff829cfSHong Zhang   }
54195b5ac22SHong Zhang #if defined(PETSC_USE_INFO)
54295b5ac22SHong Zhang   if (ai[mbs] != 0) {
54395b5ac22SHong Zhang     PetscReal af = fact->info.fill_ratio_needed;
54495b5ac22SHong Zhang     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
54595b5ac22SHong Zhang     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
54695b5ac22SHong Zhang     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
54795b5ac22SHong Zhang   } else {
54895b5ac22SHong Zhang     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
54995b5ac22SHong Zhang   }
55095b5ac22SHong Zhang #endif
551d595f711SHong Zhang   ierr = MatSeqSBAIJSetNumericFactorization_inplace(fact,perm_identity);CHKERRQ(ierr);
552064503c5SHong Zhang   PetscFunctionReturn(0);
553064503c5SHong Zhang }
554d595f711SHong Zhang 
5554a2ae208SSatish Balay #undef __FUNCT__
5564a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N"
5570481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat C,Mat A,const MatFactorInfo *info)
55849b5e25fSSatish Balay {
5594c16a6a6SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
5604c16a6a6SHong Zhang   IS             perm = b->row;
5616849ba73SBarry Smith   PetscErrorCode ierr;
5625d0c19d7SBarry Smith   const PetscInt *ai,*aj,*perm_ptr,mbs=a->mbs,*bi=b->i,*bj=b->j;
5635d0c19d7SBarry Smith   PetscInt       i,j;
5645d0c19d7SBarry Smith   PetscInt       *a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
565d0f46423SBarry Smith   PetscInt       bs=A->rmap->bs,bs2 = a->bs2,bslog = 0;
5664c16a6a6SHong Zhang   MatScalar      *ba = b->a,*aa,*ap,*dk,*uik;
5674c16a6a6SHong Zhang   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
56828de702eSHong Zhang   MatScalar      *work;
56913f74950SBarry Smith   PetscInt       *pivots;
5704c16a6a6SHong Zhang 
5714c16a6a6SHong Zhang   PetscFunctionBegin;
5724c16a6a6SHong Zhang   /* initialization */
57382502324SSatish Balay   ierr = PetscMalloc(bs2*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
5744c16a6a6SHong Zhang   ierr = PetscMemzero(rtmp,bs2*mbs*sizeof(MatScalar));CHKERRQ(ierr);
575d8c74875SBarry Smith   ierr = PetscMalloc2(mbs,PetscInt,&il,mbs,PetscInt,&jl);CHKERRQ(ierr);
5764c16a6a6SHong Zhang   for (i=0; i<mbs; i++) {
5774c16a6a6SHong Zhang     jl[i] = mbs; il[0] = 0;
5784c16a6a6SHong Zhang   }
579d8c74875SBarry Smith   ierr = PetscMalloc3(bs2,MatScalar,&dk,bs2,MatScalar,&uik,bs,MatScalar,&work);CHKERRQ(ierr);
58013f74950SBarry Smith   ierr = PetscMalloc(bs*sizeof(PetscInt),&pivots);CHKERRQ(ierr);
5814c16a6a6SHong Zhang 
5824c16a6a6SHong Zhang   ierr  = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr);
5834c16a6a6SHong Zhang 
5844c16a6a6SHong Zhang   /* check permutation */
5854c16a6a6SHong Zhang   if (!a->permute){
5864c16a6a6SHong Zhang     ai = a->i; aj = a->j; aa = a->a;
5874c16a6a6SHong Zhang   } else {
5884c16a6a6SHong Zhang     ai   = a->inew; aj = a->jnew;
58982502324SSatish Balay     ierr = PetscMalloc(bs2*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr);
5904c16a6a6SHong Zhang     ierr = PetscMemcpy(aa,a->a,bs2*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr);
59113f74950SBarry Smith     ierr = PetscMalloc(ai[mbs]*sizeof(PetscInt),&a2anew);CHKERRQ(ierr);
59213f74950SBarry Smith     ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(PetscInt));CHKERRQ(ierr);
5934c16a6a6SHong Zhang 
594187a9f4bSHong Zhang     /* flops in while loop */
595187a9f4bSHong Zhang     bslog = 2*bs*bs2;
596187a9f4bSHong Zhang 
5974c16a6a6SHong Zhang     for (i=0; i<mbs; i++){
5984c16a6a6SHong Zhang       jmin = ai[i]; jmax = ai[i+1];
5994c16a6a6SHong Zhang       for (j=jmin; j<jmax; j++){
6004c16a6a6SHong Zhang         while (a2anew[j] != j){
6014c16a6a6SHong Zhang           k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
6024c16a6a6SHong Zhang           for (k1=0; k1<bs2; k1++){
6034c16a6a6SHong Zhang             dk[k1]       = aa[k*bs2+k1];
6044c16a6a6SHong Zhang             aa[k*bs2+k1] = aa[j*bs2+k1];
6054c16a6a6SHong Zhang             aa[j*bs2+k1] = dk[k1];
6064c16a6a6SHong Zhang           }
6074c16a6a6SHong Zhang         }
6084c16a6a6SHong Zhang         /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
6094c16a6a6SHong Zhang         if (i > aj[j]){
6104c16a6a6SHong Zhang           /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */
6114c16a6a6SHong Zhang           ap = aa + j*bs2;                     /* ptr to the beginning of j-th block of aa */
6124c16a6a6SHong Zhang           for (k=0; k<bs2; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */
6134c16a6a6SHong Zhang           for (k=0; k<bs; k++){               /* j-th block of aa <- dk^T */
6144c16a6a6SHong Zhang             for (k1=0; k1<bs; k1++) *ap++ = dk[k + bs*k1];
6154c16a6a6SHong Zhang           }
6164c16a6a6SHong Zhang         }
6174c16a6a6SHong Zhang       }
6184c16a6a6SHong Zhang     }
619323b833fSBarry Smith     ierr = PetscFree(a2anew);CHKERRQ(ierr);
6204c16a6a6SHong Zhang   }
6214c16a6a6SHong Zhang 
6224c16a6a6SHong Zhang   /* for each row k */
6234c16a6a6SHong Zhang   for (k = 0; k<mbs; k++){
6244c16a6a6SHong Zhang 
6254c16a6a6SHong Zhang     /*initialize k-th row with elements nonzero in row perm(k) of A */
6264c16a6a6SHong Zhang     jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];
627057f5ba7SHong Zhang 
6284c16a6a6SHong Zhang     ap = aa + jmin*bs2;
6294c16a6a6SHong Zhang     for (j = jmin; j < jmax; j++){
6304c16a6a6SHong Zhang       vj = perm_ptr[aj[j]];         /* block col. index */
6314c16a6a6SHong Zhang       rtmp_ptr = rtmp + vj*bs2;
6324c16a6a6SHong Zhang       for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++;
6334c16a6a6SHong Zhang     }
6344c16a6a6SHong Zhang 
6354c16a6a6SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
6364c16a6a6SHong Zhang     ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
6374c16a6a6SHong Zhang     i = jl[k]; /* first row to be added to k_th row  */
6384c16a6a6SHong Zhang 
639057f5ba7SHong Zhang     while (i < k){
6404c16a6a6SHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
6414c16a6a6SHong Zhang 
6424c16a6a6SHong Zhang       /* compute multiplier */
6434c16a6a6SHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
6444c16a6a6SHong Zhang 
6454c16a6a6SHong Zhang       /* uik = -inv(Di)*U_bar(i,k) */
6464c16a6a6SHong Zhang       diag = ba + i*bs2;
6474c16a6a6SHong Zhang       u    = ba + ili*bs2;
6484c16a6a6SHong Zhang       ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
649*96b95a6bSBarry Smith       PetscKernel_A_gets_A_minus_B_times_C(bs,uik,diag,u);
6504c16a6a6SHong Zhang 
6514c16a6a6SHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
652*96b95a6bSBarry Smith       PetscKernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u);
653dc0b31edSSatish Balay       ierr = PetscLogFlops(bslog*2.0);CHKERRQ(ierr);
6544c16a6a6SHong Zhang 
6554c16a6a6SHong Zhang       /* update -U(i,k) */
6564c16a6a6SHong Zhang       ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
6574c16a6a6SHong Zhang 
6584c16a6a6SHong Zhang       /* add multiple of row i to k-th row ... */
6594c16a6a6SHong Zhang       jmin = ili + 1; jmax = bi[i+1];
6604c16a6a6SHong Zhang       if (jmin < jmax){
6614c16a6a6SHong Zhang         for (j=jmin; j<jmax; j++) {
6624c16a6a6SHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j) */
6634c16a6a6SHong Zhang           rtmp_ptr = rtmp + bj[j]*bs2;
6644c16a6a6SHong Zhang           u = ba + j*bs2;
665*96b95a6bSBarry Smith           PetscKernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u);
6664c16a6a6SHong Zhang         }
667187a9f4bSHong Zhang         ierr = PetscLogFlops(bslog*(jmax-jmin));CHKERRQ(ierr);
6684c16a6a6SHong Zhang 
6694c16a6a6SHong Zhang         /* ... add i to row list for next nonzero entry */
6704c16a6a6SHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
6714c16a6a6SHong Zhang         j     = bj[jmin];
6724c16a6a6SHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
6734c16a6a6SHong Zhang       }
6744c16a6a6SHong Zhang       i = nexti;
6754c16a6a6SHong Zhang     }
6764c16a6a6SHong Zhang 
6774c16a6a6SHong Zhang     /* save nonzero entries in k-th row of U ... */
6784c16a6a6SHong Zhang 
6794c16a6a6SHong Zhang     /* invert diagonal block */
6804c16a6a6SHong Zhang     diag = ba+k*bs2;
6814c16a6a6SHong Zhang     ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr);
682*96b95a6bSBarry Smith     ierr = PetscKernel_A_gets_inverse_A(bs,diag,pivots,work);CHKERRQ(ierr);
6834c16a6a6SHong Zhang 
6844c16a6a6SHong Zhang     jmin = bi[k]; jmax = bi[k+1];
6854c16a6a6SHong Zhang     if (jmin < jmax) {
6864c16a6a6SHong Zhang       for (j=jmin; j<jmax; j++){
6874c16a6a6SHong Zhang          vj = bj[j];           /* block col. index of U */
6884c16a6a6SHong Zhang          u   = ba + j*bs2;
6894c16a6a6SHong Zhang          rtmp_ptr = rtmp + vj*bs2;
6904c16a6a6SHong Zhang          for (k1=0; k1<bs2; k1++){
6914c16a6a6SHong Zhang            *u++        = *rtmp_ptr;
6924c16a6a6SHong Zhang            *rtmp_ptr++ = 0.0;
6934c16a6a6SHong Zhang          }
6944c16a6a6SHong Zhang       }
6954c16a6a6SHong Zhang 
6964c16a6a6SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
6974c16a6a6SHong Zhang       il[k] = jmin;
6984c16a6a6SHong Zhang       i     = bj[jmin];
6994c16a6a6SHong Zhang       jl[k] = jl[i]; jl[i] = k;
7004c16a6a6SHong Zhang     }
7014c16a6a6SHong Zhang   }
7024c16a6a6SHong Zhang 
7034c16a6a6SHong Zhang   ierr = PetscFree(rtmp);CHKERRQ(ierr);
704d8c74875SBarry Smith   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
705d8c74875SBarry Smith   ierr = PetscFree3(dk,uik,work);CHKERRQ(ierr);
7064c16a6a6SHong Zhang   ierr = PetscFree(pivots);CHKERRQ(ierr);
7074c16a6a6SHong Zhang   if (a->permute){
7084c16a6a6SHong Zhang     ierr = PetscFree(aa);CHKERRQ(ierr);
7094c16a6a6SHong Zhang   }
7104c16a6a6SHong Zhang 
7114c16a6a6SHong Zhang   ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr);
7124f79d315SHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_N_inplace;
7134f79d315SHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_N_inplace;
7144f79d315SHong Zhang   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_N_inplace;
7154f79d315SHong Zhang   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_N_inplace;
716db4efbfdSBarry Smith 
7174c16a6a6SHong Zhang   C->assembled    = PETSC_TRUE;
7184c16a6a6SHong Zhang   C->preallocated = PETSC_TRUE;
719efee365bSSatish Balay   ierr = PetscLogFlops(1.3333*bs*bs2*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
7204c16a6a6SHong Zhang   PetscFunctionReturn(0);
7214c16a6a6SHong Zhang }
722d4edadd4SHong Zhang 
7234a2ae208SSatish Balay #undef __FUNCT__
7244a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering"
7250481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info)
726671cb588SHong Zhang {
727671cb588SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
728dfbe8321SBarry Smith   PetscErrorCode ierr;
72913f74950SBarry Smith   PetscInt       i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
73013f74950SBarry Smith   PetscInt       *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
731d0f46423SBarry Smith   PetscInt       bs=A->rmap->bs,bs2 = a->bs2,bslog;
732671cb588SHong Zhang   MatScalar      *ba = b->a,*aa,*ap,*dk,*uik;
733671cb588SHong Zhang   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
73428de702eSHong Zhang   MatScalar      *work;
73513f74950SBarry Smith   PetscInt       *pivots;
736671cb588SHong Zhang 
737671cb588SHong Zhang   PetscFunctionBegin;
73882502324SSatish Balay   ierr = PetscMalloc(bs2*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
739671cb588SHong Zhang   ierr = PetscMemzero(rtmp,bs2*mbs*sizeof(MatScalar));CHKERRQ(ierr);
740d8c74875SBarry Smith   ierr = PetscMalloc2(mbs,PetscInt,&il,mbs,PetscInt,&jl);CHKERRQ(ierr);
741671cb588SHong Zhang   for (i=0; i<mbs; i++) {
742671cb588SHong Zhang     jl[i] = mbs; il[0] = 0;
743671cb588SHong Zhang   }
744d8c74875SBarry Smith   ierr = PetscMalloc3(bs2,MatScalar,&dk,bs2,MatScalar,&uik,bs,MatScalar,&work);CHKERRQ(ierr);
74513f74950SBarry Smith   ierr = PetscMalloc(bs*sizeof(PetscInt),&pivots);CHKERRQ(ierr);
746671cb588SHong Zhang 
747671cb588SHong Zhang   ai = a->i; aj = a->j; aa = a->a;
748671cb588SHong Zhang 
749187a9f4bSHong Zhang   /* flops in while loop */
750187a9f4bSHong Zhang   bslog = 2*bs*bs2;
751187a9f4bSHong Zhang 
752671cb588SHong Zhang   /* for each row k */
753671cb588SHong Zhang   for (k = 0; k<mbs; k++){
754671cb588SHong Zhang 
755671cb588SHong Zhang     /*initialize k-th row with elements nonzero in row k of A */
756671cb588SHong Zhang     jmin = ai[k]; jmax = ai[k+1];
757671cb588SHong Zhang     ap = aa + jmin*bs2;
758671cb588SHong Zhang     for (j = jmin; j < jmax; j++){
759671cb588SHong Zhang       vj = aj[j];         /* block col. index */
760671cb588SHong Zhang       rtmp_ptr = rtmp + vj*bs2;
761671cb588SHong Zhang       for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++;
762671cb588SHong Zhang     }
763671cb588SHong Zhang 
764671cb588SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
765671cb588SHong Zhang     ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
766671cb588SHong Zhang     i = jl[k]; /* first row to be added to k_th row  */
767671cb588SHong Zhang 
768057f5ba7SHong Zhang     while (i < k){
769671cb588SHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
770671cb588SHong Zhang 
771671cb588SHong Zhang       /* compute multiplier */
772671cb588SHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
773671cb588SHong Zhang 
774671cb588SHong Zhang       /* uik = -inv(Di)*U_bar(i,k) */
775671cb588SHong Zhang       diag = ba + i*bs2;
776671cb588SHong Zhang       u    = ba + ili*bs2;
777671cb588SHong Zhang       ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
778*96b95a6bSBarry Smith       PetscKernel_A_gets_A_minus_B_times_C(bs,uik,diag,u);
779671cb588SHong Zhang 
780671cb588SHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
781*96b95a6bSBarry Smith       PetscKernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u);
782dc0b31edSSatish Balay       ierr = PetscLogFlops(bslog*2.0);CHKERRQ(ierr);
783671cb588SHong Zhang 
784671cb588SHong Zhang       /* update -U(i,k) */
785671cb588SHong Zhang       ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
786671cb588SHong Zhang 
787671cb588SHong Zhang       /* add multiple of row i to k-th row ... */
788671cb588SHong Zhang       jmin = ili + 1; jmax = bi[i+1];
789671cb588SHong Zhang       if (jmin < jmax){
790671cb588SHong Zhang         for (j=jmin; j<jmax; j++) {
791671cb588SHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j) */
792671cb588SHong Zhang           rtmp_ptr = rtmp + bj[j]*bs2;
793671cb588SHong Zhang           u = ba + j*bs2;
794*96b95a6bSBarry Smith           PetscKernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u);
795671cb588SHong Zhang         }
796187a9f4bSHong Zhang         ierr = PetscLogFlops(bslog*(jmax-jmin));CHKERRQ(ierr);
797671cb588SHong Zhang 
798671cb588SHong Zhang         /* ... add i to row list for next nonzero entry */
799671cb588SHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
800671cb588SHong Zhang         j     = bj[jmin];
801671cb588SHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
802671cb588SHong Zhang       }
803671cb588SHong Zhang       i = nexti;
804671cb588SHong Zhang     }
805671cb588SHong Zhang 
806671cb588SHong Zhang     /* save nonzero entries in k-th row of U ... */
807671cb588SHong Zhang 
808671cb588SHong Zhang     /* invert diagonal block */
809671cb588SHong Zhang     diag = ba+k*bs2;
810671cb588SHong Zhang     ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr);
811*96b95a6bSBarry Smith     ierr = PetscKernel_A_gets_inverse_A(bs,diag,pivots,work);CHKERRQ(ierr);
812671cb588SHong Zhang 
813671cb588SHong Zhang     jmin = bi[k]; jmax = bi[k+1];
814671cb588SHong Zhang     if (jmin < jmax) {
815671cb588SHong Zhang       for (j=jmin; j<jmax; j++){
816671cb588SHong Zhang          vj = bj[j];           /* block col. index of U */
817671cb588SHong Zhang          u   = ba + j*bs2;
818671cb588SHong Zhang          rtmp_ptr = rtmp + vj*bs2;
819671cb588SHong Zhang          for (k1=0; k1<bs2; k1++){
820671cb588SHong Zhang            *u++        = *rtmp_ptr;
821671cb588SHong Zhang            *rtmp_ptr++ = 0.0;
822671cb588SHong Zhang          }
823671cb588SHong Zhang       }
824671cb588SHong Zhang 
825671cb588SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
826671cb588SHong Zhang       il[k] = jmin;
827671cb588SHong Zhang       i     = bj[jmin];
828671cb588SHong Zhang       jl[k] = jl[i]; jl[i] = k;
829671cb588SHong Zhang     }
830671cb588SHong Zhang   }
831671cb588SHong Zhang 
832671cb588SHong Zhang   ierr = PetscFree(rtmp);CHKERRQ(ierr);
833d8c74875SBarry Smith   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
834d8c74875SBarry Smith   ierr = PetscFree3(dk,uik,work);CHKERRQ(ierr);
835671cb588SHong Zhang   ierr = PetscFree(pivots);CHKERRQ(ierr);
836671cb588SHong Zhang 
8374f79d315SHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_N_NaturalOrdering_inplace;
8384f79d315SHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_N_NaturalOrdering_inplace;
8394f79d315SHong Zhang   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_N_NaturalOrdering_inplace;
8404f79d315SHong Zhang   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_N_NaturalOrdering_inplace;
841671cb588SHong Zhang   C->assembled = PETSC_TRUE;
842671cb588SHong Zhang   C->preallocated = PETSC_TRUE;
843efee365bSSatish Balay   ierr = PetscLogFlops(1.3333*bs*bs2*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
844671cb588SHong Zhang   PetscFunctionReturn(0);
845671cb588SHong Zhang }
846671cb588SHong Zhang 
84749b5e25fSSatish Balay /*
848fcf159c0SHong Zhang     Numeric U^T*D*U factorization for SBAIJ format. Modified from SNF of YSMP.
849cc0c071aSHong Zhang     Version for blocks 2 by 2.
85049b5e25fSSatish Balay */
8514a2ae208SSatish Balay #undef __FUNCT__
8524a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2"
8530481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat C,Mat A,const MatFactorInfo *info)
85449b5e25fSSatish Balay {
85591602c66SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
856cc0c071aSHong Zhang   IS             perm = b->row;
8576849ba73SBarry Smith   PetscErrorCode ierr;
8585d0c19d7SBarry Smith   const PetscInt *ai,*aj,*perm_ptr;
8595d0c19d7SBarry Smith   PetscInt       i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
8605d0c19d7SBarry Smith   PetscInt       *a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
861d8c74875SBarry Smith   MatScalar      *ba = b->a,*aa,*ap;
862d8c74875SBarry Smith   MatScalar      *u,*diag,*rtmp,*rtmp_ptr,dk[4],uik[4];
86314d2772eSHong Zhang   PetscReal      shift = info->shiftamount;
86449b5e25fSSatish Balay 
86549b5e25fSSatish Balay   PetscFunctionBegin;
86691602c66SHong Zhang   /* initialization */
86791602c66SHong Zhang   /* il and jl record the first nonzero element in each row of the accessing
86891602c66SHong Zhang      window U(0:k, k:mbs-1).
86991602c66SHong Zhang      jl:    list of rows to be added to uneliminated rows
87091602c66SHong Zhang             i>= k: jl(i) is the first row to be added to row i
87191602c66SHong Zhang             i<  k: jl(i) is the row following row i in some list of rows
87291602c66SHong Zhang             jl(i) = mbs indicates the end of a list
87391602c66SHong Zhang      il(i): points to the first nonzero element in columns k,...,mbs-1 of
87491602c66SHong Zhang             row i of U */
87582502324SSatish Balay   ierr = PetscMalloc(4*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
876cc0c071aSHong Zhang   ierr = PetscMemzero(rtmp,4*mbs*sizeof(MatScalar));CHKERRQ(ierr);
877d8c74875SBarry Smith   ierr = PetscMalloc2(mbs,PetscInt,&il,mbs,PetscInt,&jl);CHKERRQ(ierr);
87891602c66SHong Zhang   for (i=0; i<mbs; i++) {
8793845f261SHong Zhang     jl[i] = mbs; il[0] = 0;
88091602c66SHong Zhang   }
881cc0c071aSHong Zhang   ierr = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr);
882cc0c071aSHong Zhang 
883cc0c071aSHong Zhang   /* check permutation */
884cc0c071aSHong Zhang   if (!a->permute){
885cc0c071aSHong Zhang     ai = a->i; aj = a->j; aa = a->a;
886cc0c071aSHong Zhang   } else {
887cc0c071aSHong Zhang     ai   = a->inew; aj = a->jnew;
88882502324SSatish Balay     ierr = PetscMalloc(4*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr);
889cc0c071aSHong Zhang     ierr = PetscMemcpy(aa,a->a,4*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr);
89013f74950SBarry Smith     ierr = PetscMalloc(ai[mbs]*sizeof(PetscInt),&a2anew);CHKERRQ(ierr);
89113f74950SBarry Smith     ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(PetscInt));CHKERRQ(ierr);
892cc0c071aSHong Zhang 
893cc0c071aSHong Zhang     for (i=0; i<mbs; i++){
894cc0c071aSHong Zhang       jmin = ai[i]; jmax = ai[i+1];
895cc0c071aSHong Zhang       for (j=jmin; j<jmax; j++){
896cc0c071aSHong Zhang         while (a2anew[j] != j){
897cc0c071aSHong Zhang           k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
898cc0c071aSHong Zhang           for (k1=0; k1<4; k1++){
899cc0c071aSHong Zhang             dk[k1]       = aa[k*4+k1];
900cc0c071aSHong Zhang             aa[k*4+k1] = aa[j*4+k1];
901cc0c071aSHong Zhang             aa[j*4+k1] = dk[k1];
902cc0c071aSHong Zhang           }
903cc0c071aSHong Zhang         }
904cc0c071aSHong Zhang         /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
905cc0c071aSHong Zhang         if (i > aj[j]){
906a1723e09SHong Zhang           /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */
907cc0c071aSHong Zhang           ap = aa + j*4;     /* ptr to the beginning of the block */
908cc0c071aSHong Zhang           dk[1] = ap[1];     /* swap ap[1] and ap[2] */
909cc0c071aSHong Zhang           ap[1] = ap[2];
910cc0c071aSHong Zhang           ap[2] = dk[1];
911cc0c071aSHong Zhang         }
912cc0c071aSHong Zhang       }
913cc0c071aSHong Zhang     }
914ac355199SBarry Smith     ierr = PetscFree(a2anew);CHKERRQ(ierr);
915cc0c071aSHong Zhang   }
9163845f261SHong Zhang 
91791602c66SHong Zhang   /* for each row k */
91891602c66SHong Zhang   for (k = 0; k<mbs; k++){
91991602c66SHong Zhang 
92091602c66SHong Zhang     /*initialize k-th row with elements nonzero in row perm(k) of A */
921cc0c071aSHong Zhang     jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];
922cc0c071aSHong Zhang     ap = aa + jmin*4;
92391602c66SHong Zhang     for (j = jmin; j < jmax; j++){
924cc0c071aSHong Zhang       vj = perm_ptr[aj[j]];         /* block col. index */
925cc0c071aSHong Zhang       rtmp_ptr = rtmp + vj*4;
926cc0c071aSHong Zhang       for (i=0; i<4; i++) *rtmp_ptr++ = *ap++;
92791602c66SHong Zhang     }
92891602c66SHong Zhang 
92991602c66SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
930cc0c071aSHong Zhang     ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr);
93191602c66SHong Zhang     i = jl[k]; /* first row to be added to k_th row  */
93291602c66SHong Zhang 
933057f5ba7SHong Zhang     while (i < k){
93491602c66SHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
93591602c66SHong Zhang 
9363845f261SHong Zhang       /* compute multiplier */
93791602c66SHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
9383845f261SHong Zhang 
9393845f261SHong Zhang       /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */
940cc0c071aSHong Zhang       diag = ba + i*4;
941cc0c071aSHong Zhang       u    = ba + ili*4;
942cc0c071aSHong Zhang       uik[0] = -(diag[0]*u[0] + diag[2]*u[1]);
943cc0c071aSHong Zhang       uik[1] = -(diag[1]*u[0] + diag[3]*u[1]);
944cc0c071aSHong Zhang       uik[2] = -(diag[0]*u[2] + diag[2]*u[3]);
945cc0c071aSHong Zhang       uik[3] = -(diag[1]*u[2] + diag[3]*u[3]);
9463845f261SHong Zhang 
9473845f261SHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */
948cc0c071aSHong Zhang       dk[0] += uik[0]*u[0] + uik[1]*u[1];
949cc0c071aSHong Zhang       dk[1] += uik[2]*u[0] + uik[3]*u[1];
950cc0c071aSHong Zhang       dk[2] += uik[0]*u[2] + uik[1]*u[3];
951cc0c071aSHong Zhang       dk[3] += uik[2]*u[2] + uik[3]*u[3];
9523845f261SHong Zhang 
953dc0b31edSSatish Balay       ierr = PetscLogFlops(16.0*2.0);CHKERRQ(ierr);
954187a9f4bSHong Zhang 
9553845f261SHong Zhang       /* update -U(i,k): ba[ili] = uik */
956cc0c071aSHong Zhang       ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr);
95791602c66SHong Zhang 
95891602c66SHong Zhang       /* add multiple of row i to k-th row ... */
95991602c66SHong Zhang       jmin = ili + 1; jmax = bi[i+1];
96091602c66SHong Zhang       if (jmin < jmax){
9613845f261SHong Zhang         for (j=jmin; j<jmax; j++) {
9623845f261SHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */
963cc0c071aSHong Zhang           rtmp_ptr = rtmp + bj[j]*4;
964cc0c071aSHong Zhang           u = ba + j*4;
965cc0c071aSHong Zhang           rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1];
966cc0c071aSHong Zhang           rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1];
967cc0c071aSHong Zhang           rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3];
968cc0c071aSHong Zhang           rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3];
9693845f261SHong Zhang         }
970dc0b31edSSatish Balay         ierr = PetscLogFlops(16.0*(jmax-jmin));CHKERRQ(ierr);
9713845f261SHong Zhang 
97291602c66SHong Zhang         /* ... add i to row list for next nonzero entry */
97391602c66SHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
97491602c66SHong Zhang         j     = bj[jmin];
97591602c66SHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
97691602c66SHong Zhang       }
977a1723e09SHong Zhang       i = nexti;
97891602c66SHong Zhang     }
979cc0c071aSHong Zhang 
98091602c66SHong Zhang     /* save nonzero entries in k-th row of U ... */
9813845f261SHong Zhang 
982cc0c071aSHong Zhang     /* invert diagonal block */
983cc0c071aSHong Zhang     diag = ba+k*4;
984cc0c071aSHong Zhang     ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr);
985*96b95a6bSBarry Smith     ierr = PetscKernel_A_gets_inverse_A_2(diag,shift);CHKERRQ(ierr);
9863845f261SHong Zhang 
98791602c66SHong Zhang     jmin = bi[k]; jmax = bi[k+1];
98891602c66SHong Zhang     if (jmin < jmax) {
98991602c66SHong Zhang       for (j=jmin; j<jmax; j++){
990cc0c071aSHong Zhang          vj = bj[j];           /* block col. index of U */
991cc0c071aSHong Zhang          u   = ba + j*4;
992cc0c071aSHong Zhang          rtmp_ptr = rtmp + vj*4;
993cc0c071aSHong Zhang          for (k1=0; k1<4; k1++){
994cc0c071aSHong Zhang            *u++        = *rtmp_ptr;
995cc0c071aSHong Zhang            *rtmp_ptr++ = 0.0;
996cc0c071aSHong Zhang          }
997cc0c071aSHong Zhang       }
9983845f261SHong Zhang 
99991602c66SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
100091602c66SHong Zhang       il[k] = jmin;
100191602c66SHong Zhang       i     = bj[jmin];
100291602c66SHong Zhang       jl[k] = jl[i]; jl[i] = k;
100391602c66SHong Zhang     }
100491602c66SHong Zhang   }
10053845f261SHong Zhang 
100649b5e25fSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
1007d8c74875SBarry Smith   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
100891602c66SHong Zhang   if (a->permute) {
100991602c66SHong Zhang     ierr = PetscFree(aa);CHKERRQ(ierr);
101091602c66SHong Zhang   }
1011cc0c071aSHong Zhang   ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr);
10124f79d315SHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_2_inplace;
10134f79d315SHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_2_inplace;
101449b5e25fSSatish Balay   C->assembled = PETSC_TRUE;
10155c0bcdfcSHong Zhang   C->preallocated = PETSC_TRUE;
1016efee365bSSatish Balay   ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
101749b5e25fSSatish Balay   PetscFunctionReturn(0);
101849b5e25fSSatish Balay }
101991602c66SHong Zhang 
102049b5e25fSSatish Balay /*
102149b5e25fSSatish Balay       Version for when blocks are 2 by 2 Using natural ordering
102249b5e25fSSatish Balay */
10234a2ae208SSatish Balay #undef __FUNCT__
10244a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering"
10250481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info)
102649b5e25fSSatish Balay {
1027ab27746eSHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
1028dfbe8321SBarry Smith   PetscErrorCode ierr;
102913f74950SBarry Smith   PetscInt       i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
103013f74950SBarry Smith   PetscInt       *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
1031d8c74875SBarry Smith   MatScalar      *ba = b->a,*aa,*ap,dk[8],uik[8];
1032ab27746eSHong Zhang   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
103314d2772eSHong Zhang   PetscReal      shift = info->shiftamount;
103449b5e25fSSatish Balay 
103549b5e25fSSatish Balay   PetscFunctionBegin;
1036ab27746eSHong Zhang   /* initialization */
1037ab27746eSHong Zhang   /* il and jl record the first nonzero element in each row of the accessing
1038ab27746eSHong Zhang      window U(0:k, k:mbs-1).
1039ab27746eSHong Zhang      jl:    list of rows to be added to uneliminated rows
1040ab27746eSHong Zhang             i>= k: jl(i) is the first row to be added to row i
1041ab27746eSHong Zhang             i<  k: jl(i) is the row following row i in some list of rows
1042ab27746eSHong Zhang             jl(i) = mbs indicates the end of a list
1043ab27746eSHong Zhang      il(i): points to the first nonzero element in columns k,...,mbs-1 of
1044ab27746eSHong Zhang             row i of U */
104582502324SSatish Balay   ierr = PetscMalloc(4*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
1046ab27746eSHong Zhang   ierr = PetscMemzero(rtmp,4*mbs*sizeof(MatScalar));CHKERRQ(ierr);
1047d8c74875SBarry Smith   ierr = PetscMalloc2(mbs,PetscInt,&il,mbs,PetscInt,&jl);CHKERRQ(ierr);
1048ab27746eSHong Zhang   for (i=0; i<mbs; i++) {
1049ab27746eSHong Zhang     jl[i] = mbs; il[0] = 0;
105049b5e25fSSatish Balay   }
1051ab27746eSHong Zhang   ai = a->i; aj = a->j; aa = a->a;
1052ab27746eSHong Zhang 
1053ab27746eSHong Zhang   /* for each row k */
1054ab27746eSHong Zhang   for (k = 0; k<mbs; k++){
1055ab27746eSHong Zhang 
1056ab27746eSHong Zhang     /*initialize k-th row with elements nonzero in row k of A */
1057ab27746eSHong Zhang     jmin = ai[k]; jmax = ai[k+1];
1058ab27746eSHong Zhang     ap = aa + jmin*4;
1059ab27746eSHong Zhang     for (j = jmin; j < jmax; j++){
1060ab27746eSHong Zhang       vj = aj[j];         /* block col. index */
1061ab27746eSHong Zhang       rtmp_ptr = rtmp + vj*4;
1062ab27746eSHong Zhang       for (i=0; i<4; i++) *rtmp_ptr++ = *ap++;
106349b5e25fSSatish Balay     }
1064ab27746eSHong Zhang 
1065ab27746eSHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
1066ab27746eSHong Zhang     ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr);
1067ab27746eSHong Zhang     i = jl[k]; /* first row to be added to k_th row  */
1068ab27746eSHong Zhang 
1069057f5ba7SHong Zhang     while (i < k){
1070ab27746eSHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
1071ab27746eSHong Zhang 
1072ab27746eSHong Zhang       /* compute multiplier */
1073ab27746eSHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
1074ab27746eSHong Zhang 
1075ab27746eSHong Zhang       /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */
1076ab27746eSHong Zhang       diag = ba + i*4;
1077ab27746eSHong Zhang       u    = ba + ili*4;
1078ab27746eSHong Zhang       uik[0] = -(diag[0]*u[0] + diag[2]*u[1]);
1079ab27746eSHong Zhang       uik[1] = -(diag[1]*u[0] + diag[3]*u[1]);
1080ab27746eSHong Zhang       uik[2] = -(diag[0]*u[2] + diag[2]*u[3]);
1081ab27746eSHong Zhang       uik[3] = -(diag[1]*u[2] + diag[3]*u[3]);
1082ab27746eSHong Zhang 
1083ab27746eSHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */
1084ab27746eSHong Zhang       dk[0] += uik[0]*u[0] + uik[1]*u[1];
1085ab27746eSHong Zhang       dk[1] += uik[2]*u[0] + uik[3]*u[1];
1086ab27746eSHong Zhang       dk[2] += uik[0]*u[2] + uik[1]*u[3];
1087ab27746eSHong Zhang       dk[3] += uik[2]*u[2] + uik[3]*u[3];
1088ab27746eSHong Zhang 
1089dc0b31edSSatish Balay       ierr = PetscLogFlops(16.0*2.0);CHKERRQ(ierr);
1090187a9f4bSHong Zhang 
1091ab27746eSHong Zhang       /* update -U(i,k): ba[ili] = uik */
1092ab27746eSHong Zhang       ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr);
1093ab27746eSHong Zhang 
1094ab27746eSHong Zhang       /* add multiple of row i to k-th row ... */
1095ab27746eSHong Zhang       jmin = ili + 1; jmax = bi[i+1];
1096ab27746eSHong Zhang       if (jmin < jmax){
1097ab27746eSHong Zhang         for (j=jmin; j<jmax; j++) {
1098ab27746eSHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */
1099ab27746eSHong Zhang           rtmp_ptr = rtmp + bj[j]*4;
1100ab27746eSHong Zhang           u = ba + j*4;
1101ab27746eSHong Zhang           rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1];
1102ab27746eSHong Zhang           rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1];
1103ab27746eSHong Zhang           rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3];
1104ab27746eSHong Zhang           rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3];
110549b5e25fSSatish Balay         }
1106dc0b31edSSatish Balay         ierr = PetscLogFlops(16.0*(jmax-jmin));CHKERRQ(ierr);
1107ab27746eSHong Zhang 
1108ab27746eSHong Zhang         /* ... add i to row list for next nonzero entry */
1109ab27746eSHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
1110ab27746eSHong Zhang         j     = bj[jmin];
1111ab27746eSHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
111249b5e25fSSatish Balay       }
1113ab27746eSHong Zhang       i = nexti;
111449b5e25fSSatish Balay     }
1115ab27746eSHong Zhang 
1116ab27746eSHong Zhang     /* save nonzero entries in k-th row of U ... */
1117ab27746eSHong Zhang 
111849b5e25fSSatish Balay     /* invert diagonal block */
1119ab27746eSHong Zhang     diag = ba+k*4;
1120ab27746eSHong Zhang     ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr);
1121*96b95a6bSBarry Smith     ierr = PetscKernel_A_gets_inverse_A_2(diag,shift);CHKERRQ(ierr);
1122ab27746eSHong Zhang 
1123ab27746eSHong Zhang     jmin = bi[k]; jmax = bi[k+1];
1124ab27746eSHong Zhang     if (jmin < jmax) {
1125ab27746eSHong Zhang       for (j=jmin; j<jmax; j++){
1126ab27746eSHong Zhang          vj = bj[j];           /* block col. index of U */
1127ab27746eSHong Zhang          u   = ba + j*4;
1128ab27746eSHong Zhang          rtmp_ptr = rtmp + vj*4;
1129ab27746eSHong Zhang          for (k1=0; k1<4; k1++){
1130ab27746eSHong Zhang            *u++        = *rtmp_ptr;
1131ab27746eSHong Zhang            *rtmp_ptr++ = 0.0;
1132ab27746eSHong Zhang          }
1133ab27746eSHong Zhang       }
1134ab27746eSHong Zhang 
1135ab27746eSHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
1136ab27746eSHong Zhang       il[k] = jmin;
1137ab27746eSHong Zhang       i     = bj[jmin];
1138ab27746eSHong Zhang       jl[k] = jl[i]; jl[i] = k;
1139ab27746eSHong Zhang     }
114049b5e25fSSatish Balay   }
114149b5e25fSSatish Balay 
114249b5e25fSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
1143d8c74875SBarry Smith   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
1144ab27746eSHong Zhang 
11454f79d315SHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_2_NaturalOrdering_inplace;
11464f79d315SHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_2_NaturalOrdering_inplace;
11474f79d315SHong Zhang   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_2_NaturalOrdering_inplace;
11484f79d315SHong Zhang   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_2_NaturalOrdering_inplace;
114949b5e25fSSatish Balay   C->assembled = PETSC_TRUE;
11505c0bcdfcSHong Zhang   C->preallocated = PETSC_TRUE;
1151efee365bSSatish Balay   ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
115249b5e25fSSatish Balay   PetscFunctionReturn(0);
115349b5e25fSSatish Balay }
115449b5e25fSSatish Balay 
115549b5e25fSSatish Balay /*
115698a8e78dSHong Zhang     Numeric U^T*D*U factorization for SBAIJ format.
115791602c66SHong Zhang     Version for blocks are 1 by 1.
115849b5e25fSSatish Balay */
11594a2ae208SSatish Balay #undef __FUNCT__
1160d595f711SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1_inplace"
1161d595f711SHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_inplace(Mat C,Mat A,const MatFactorInfo *info)
116249b5e25fSSatish Balay {
116349b5e25fSSatish Balay   Mat_SeqSBAIJ   *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ *)C->data;
116449b5e25fSSatish Balay   IS             ip=b->row;
11656849ba73SBarry Smith   PetscErrorCode ierr;
11665d0c19d7SBarry Smith   const PetscInt *ai,*aj,*rip;
11675d0c19d7SBarry Smith   PetscInt       *a2anew,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j,*bcol;
1168997a0949SHong Zhang   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
1169997a0949SHong Zhang   MatScalar      *rtmp,*ba=b->a,*bval,*aa,dk,uikdi;
11700e95ead3SHong Zhang   PetscReal      rs;
11710e95ead3SHong Zhang   FactorShiftCtx sctx;
117249b5e25fSSatish Balay 
117349b5e25fSSatish Balay   PetscFunctionBegin;
11740e95ead3SHong Zhang   /* MatPivotSetUp(): initialize shift context sctx */
11750e95ead3SHong Zhang   ierr = PetscMemzero(&sctx,sizeof(FactorShiftCtx));CHKERRQ(ierr);
11763cea4cbeSHong Zhang 
117749b5e25fSSatish Balay   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
1178cb718733SHong Zhang   if (!a->permute){
11792d007305SHong Zhang     ai = a->i; aj = a->j; aa = a->a;
11802d007305SHong Zhang   } else {
11812d007305SHong Zhang     ai = a->inew; aj = a->jnew;
1182fff829cfSHong Zhang     nz = ai[mbs];
1183fff829cfSHong Zhang     ierr = PetscMalloc(nz*sizeof(MatScalar),&aa);CHKERRQ(ierr);
1184fff829cfSHong Zhang     a2anew = a->a2anew;
1185997a0949SHong Zhang     bval   = a->a;
1186fff829cfSHong Zhang     for (j=0; j<nz; j++){
1187997a0949SHong Zhang       aa[a2anew[j]] = *(bval++);
11882d007305SHong Zhang     }
11892d007305SHong Zhang   }
11902d007305SHong Zhang 
119191602c66SHong Zhang   /* initialization */
119249b5e25fSSatish Balay   /* il and jl record the first nonzero element in each row of the accessing
119349b5e25fSSatish Balay      window U(0:k, k:mbs-1).
119449b5e25fSSatish Balay      jl:    list of rows to be added to uneliminated rows
119549b5e25fSSatish Balay             i>= k: jl(i) is the first row to be added to row i
119649b5e25fSSatish Balay             i<  k: jl(i) is the row following row i in some list of rows
119749b5e25fSSatish Balay             jl(i) = mbs indicates the end of a list
119849b5e25fSSatish Balay      il(i): points to the first nonzero element in columns k,...,mbs-1 of
119949b5e25fSSatish Balay             row i of U */
1200d8c74875SBarry Smith   ierr = PetscMalloc3(mbs,MatScalar,&rtmp,mbs,PetscInt,&il,mbs,PetscInt,&jl);CHKERRQ(ierr);
1201997a0949SHong Zhang 
1202997a0949SHong Zhang   do {
120307b50cabSHong Zhang     sctx.newshift = PETSC_FALSE;
120449b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
120549b5e25fSSatish Balay       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
120649b5e25fSSatish Balay     }
1207997a0949SHong Zhang 
120849b5e25fSSatish Balay     for (k = 0; k<mbs; k++){
1209997a0949SHong Zhang       /*initialize k-th row by the perm[k]-th row of A */
121049b5e25fSSatish Balay       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
12117701de02SHong Zhang       bval = ba + bi[k];
121249b5e25fSSatish Balay       for (j = jmin; j < jmax; j++){
1213997a0949SHong Zhang         col = rip[aj[j]];
1214997a0949SHong Zhang         rtmp[col] = aa[j];
12157701de02SHong Zhang         *bval++  = 0.0; /* for in-place factorization */
121649b5e25fSSatish Balay       }
12173cea4cbeSHong Zhang 
12183cea4cbeSHong Zhang       /* shift the diagonal of the matrix */
12193cea4cbeSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
122049b5e25fSSatish Balay 
122191602c66SHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
122249b5e25fSSatish Balay       dk = rtmp[k];
122349b5e25fSSatish Balay       i = jl[k]; /* first row to be added to k_th row  */
122449b5e25fSSatish Balay 
1225057f5ba7SHong Zhang       while (i < k){
122649b5e25fSSatish Balay         nexti = jl[i]; /* next row to be added to k_th row */
122749b5e25fSSatish Balay 
1228fff829cfSHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
122949b5e25fSSatish Balay         ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
1230997a0949SHong Zhang         uikdi = - ba[ili]*ba[bi[i]];  /* diagonal(k) */
123149b5e25fSSatish Balay         dk += uikdi*ba[ili];
1232658e7b3eSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
123349b5e25fSSatish Balay 
1234997a0949SHong Zhang         /* add multiple of row i to k-th row */
123549b5e25fSSatish Balay         jmin = ili + 1; jmax = bi[i+1];
123649b5e25fSSatish Balay         if (jmin < jmax){
123749b5e25fSSatish Balay           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
1238dc0b31edSSatish Balay           ierr = PetscLogFlops(2.0*(jmax-jmin));CHKERRQ(ierr);
1239187a9f4bSHong Zhang 
1240fff829cfSHong Zhang           /* update il and jl for row i */
1241fff829cfSHong Zhang           il[i] = jmin;
1242fff829cfSHong Zhang           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
124349b5e25fSSatish Balay         }
1244ab27746eSHong Zhang         i = nexti;
124549b5e25fSSatish Balay       }
124649b5e25fSSatish Balay 
12473cea4cbeSHong Zhang       /* shift the diagonals when zero pivot is detected */
12483cea4cbeSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
12493cea4cbeSHong Zhang       rs   = 0.0;
1250997a0949SHong Zhang       jmin = bi[k]+1;
1251997a0949SHong Zhang       nz   = bi[k+1] - jmin;
1252997a0949SHong Zhang       if (nz){
1253997a0949SHong Zhang         bcol = bj + jmin;
1254997a0949SHong Zhang         while (nz--){
125589f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
125689f845c8SHong Zhang           bcol++;
1257997a0949SHong Zhang         }
1258997a0949SHong Zhang       }
12593cea4cbeSHong Zhang 
12603cea4cbeSHong Zhang       sctx.rs = rs;
12613cea4cbeSHong Zhang       sctx.pv = dk;
1262d0660788SBarry Smith       ierr = MatPivotCheck(A,info,&sctx,k);CHKERRQ(ierr);
126307b50cabSHong Zhang       if (sctx.newshift) break;    /* sctx.shift_amount is updated */
12640e95ead3SHong Zhang       dk = sctx.pv;
126549b5e25fSSatish Balay 
1266997a0949SHong Zhang       /* copy data into U(k,:) */
126798a8e78dSHong Zhang       ba[bi[k]] = 1.0/dk; /* U(k,k) */
1268fff829cfSHong Zhang       jmin = bi[k]+1; jmax = bi[k+1];
126949b5e25fSSatish Balay       if (jmin < jmax) {
127049b5e25fSSatish Balay         for (j=jmin; j<jmax; j++){
1271997a0949SHong Zhang           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
127249b5e25fSSatish Balay         }
127398a8e78dSHong Zhang         /* add the k-th row into il and jl */
127449b5e25fSSatish Balay         il[k] = jmin;
127598a8e78dSHong Zhang         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
127649b5e25fSSatish Balay       }
127749b5e25fSSatish Balay     }
127807b50cabSHong Zhang   } while (sctx.newshift);
1279d8c74875SBarry Smith   ierr = PetscFree3(rtmp,il,jl);CHKERRQ(ierr);
128098a8e78dSHong Zhang   if (a->permute){ierr = PetscFree(aa);CHKERRQ(ierr);}
128149b5e25fSSatish Balay 
128249b5e25fSSatish Balay   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
12830a3c351aSHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_1_inplace;
12844f79d315SHong Zhang   C->ops->solves         = MatSolves_SeqSBAIJ_1_inplace;
12850a3c351aSHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_1_inplace;
12860a3c351aSHong Zhang   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_1_inplace;
12870a3c351aSHong Zhang   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_1_inplace;
128849b5e25fSSatish Balay   C->assembled    = PETSC_TRUE;
12895c0bcdfcSHong Zhang   C->preallocated = PETSC_TRUE;
1290d0f46423SBarry Smith   ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr);
12913cea4cbeSHong Zhang   if (sctx.nshift){
1292f4db908eSBarry Smith     if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
12931e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
1294f4db908eSBarry Smith     } else if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
12951e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
1296997a0949SHong Zhang     }
1297997a0949SHong Zhang   }
129849b5e25fSSatish Balay   PetscFunctionReturn(0);
129949b5e25fSSatish Balay }
130049b5e25fSSatish Balay 
1301671cb588SHong Zhang /*
130280d3d5a7SHong Zhang   Version for when blocks are 1 by 1 Using natural ordering under new datastructure
130380d3d5a7SHong Zhang   Modified from MatCholeskyFactorNumeric_SeqAIJ()
1304671cb588SHong Zhang */
13054a2ae208SSatish Balay #undef __FUNCT__
13064a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering"
1307d595f711SHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info)
1308d595f711SHong Zhang {
1309d595f711SHong Zhang   Mat_SeqSBAIJ   *a=(Mat_SeqSBAIJ*)A->data;
13107b056e98SHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)B->data;
1311d595f711SHong Zhang   PetscErrorCode ierr;
1312d595f711SHong Zhang   PetscInt       i,j,mbs=A->rmap->n,*bi=b->i,*bj=b->j,*bdiag=b->diag,*bjtmp;
1313545dd064SHong Zhang   PetscInt       *ai=a->i,*aj=a->j,*ajtmp;
1314d595f711SHong Zhang   PetscInt       k,jmin,jmax,*c2r,*il,col,nexti,ili,nz;
1315d595f711SHong Zhang   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
1316d90ac83dSHong Zhang   FactorShiftCtx sctx;
1317d595f711SHong Zhang   PetscReal      rs;
1318d595f711SHong Zhang   MatScalar      d,*v;
1319d595f711SHong Zhang 
1320d595f711SHong Zhang   PetscFunctionBegin;
1321545dd064SHong Zhang   ierr = PetscMalloc3(mbs,MatScalar,&rtmp,mbs,PetscInt,&il,mbs,PetscInt,&c2r);CHKERRQ(ierr);
1322545dd064SHong Zhang 
1323d595f711SHong Zhang   /* MatPivotSetUp(): initialize shift context sctx */
1324d90ac83dSHong Zhang   ierr = PetscMemzero(&sctx,sizeof(FactorShiftCtx));CHKERRQ(ierr);
1325d595f711SHong Zhang 
1326f4db908eSBarry Smith   if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) { /* set sctx.shift_top=max{rs} */
1327d595f711SHong Zhang     sctx.shift_top = info->zeropivot;
1328545dd064SHong Zhang     ierr = PetscMemzero(rtmp,mbs*sizeof(MatScalar));CHKERRQ(ierr);
1329d595f711SHong Zhang     for (i=0; i<mbs; i++) {
1330d595f711SHong Zhang       /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
1331d595f711SHong Zhang       d  = (aa)[a->diag[i]];
1332545dd064SHong Zhang       rtmp[i] += - PetscRealPart(d); /* diagonal entry */
1333545dd064SHong Zhang       ajtmp = aj + ai[i] + 1;        /* exclude diagonal */
1334545dd064SHong Zhang       v     = aa + ai[i] + 1;
1335545dd064SHong Zhang       nz    = ai[i+1] - ai[i] - 1 ;
1336545dd064SHong Zhang       for (j=0; j<nz; j++){
1337545dd064SHong Zhang 	rtmp[i] += PetscAbsScalar(v[j]);
1338545dd064SHong Zhang         rtmp[ajtmp[j]] += PetscAbsScalar(v[j]);
1339545dd064SHong Zhang       }
13400838c725SBarry Smith       if (PetscRealPart(rtmp[i]) > sctx.shift_top) sctx.shift_top = PetscRealPart(rtmp[i]);
1341d595f711SHong Zhang     }
1342d595f711SHong Zhang     sctx.shift_top   *= 1.1;
1343d595f711SHong Zhang     sctx.nshift_max   = 5;
1344d595f711SHong Zhang     sctx.shift_lo     = 0.;
1345d595f711SHong Zhang     sctx.shift_hi     = 1.;
1346d595f711SHong Zhang   }
1347d595f711SHong Zhang 
1348d595f711SHong Zhang   /* allocate working arrays
1349d595f711SHong Zhang      c2r: linked list, keep track of pivot rows for a given column. c2r[col]: head of the list for a given col
1350d595f711SHong 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
1351d595f711SHong Zhang   */
1352d595f711SHong Zhang   do {
135307b50cabSHong Zhang     sctx.newshift = PETSC_FALSE;
1354d595f711SHong Zhang 
1355d595f711SHong Zhang     for (i=0; i<mbs; i++)  c2r[i] = mbs;
1356d595f711SHong Zhang     il[0] = 0;
1357d595f711SHong Zhang 
1358d595f711SHong Zhang     for (k = 0; k<mbs; k++){
1359d595f711SHong Zhang       /* zero rtmp */
1360d595f711SHong Zhang       nz = bi[k+1] - bi[k];
1361d595f711SHong Zhang       bjtmp = bj + bi[k];
13627b056e98SHong Zhang       for (j=0; j<nz; j++) rtmp[bjtmp[j]] = 0.0;
1363d595f711SHong Zhang 
1364d595f711SHong Zhang       /* load in initial unfactored row */
1365d595f711SHong Zhang       bval = ba + bi[k];
1366d595f711SHong Zhang       jmin = ai[k]; jmax = ai[k+1];
1367d595f711SHong Zhang       for (j = jmin; j < jmax; j++){
1368d595f711SHong Zhang         col = aj[j];
1369d595f711SHong Zhang         rtmp[col] = aa[j];
1370d595f711SHong Zhang         *bval++   = 0.0; /* for in-place factorization */
1371d595f711SHong Zhang       }
1372d595f711SHong Zhang       /* shift the diagonal of the matrix: ZeropivotApply() */
1373d595f711SHong Zhang       rtmp[k] += sctx.shift_amount;  /* shift the diagonal of the matrix */
1374d595f711SHong Zhang 
1375d595f711SHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
1376d595f711SHong Zhang       dk = rtmp[k];
1377d595f711SHong Zhang       i  = c2r[k]; /* first row to be added to k_th row  */
1378d595f711SHong Zhang 
1379d595f711SHong Zhang       while (i < k){
1380d595f711SHong Zhang         nexti = c2r[i]; /* next row to be added to k_th row */
1381d595f711SHong Zhang 
1382d595f711SHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
1383d595f711SHong Zhang         ili   = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
1384d595f711SHong Zhang         uikdi = - ba[ili]*ba[bdiag[i]];  /* diagonal(k) */
1385d595f711SHong Zhang         dk   += uikdi*ba[ili]; /* update diag[k] */
1386d595f711SHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
1387d595f711SHong Zhang 
1388d595f711SHong Zhang         /* add multiple of row i to k-th row */
1389d595f711SHong Zhang         jmin = ili + 1; jmax = bi[i+1];
1390d595f711SHong Zhang         if (jmin < jmax){
1391d595f711SHong Zhang           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
1392d595f711SHong Zhang           /* update il and c2r for row i */
1393d595f711SHong Zhang           il[i] = jmin;
1394d595f711SHong Zhang           j = bj[jmin]; c2r[i] = c2r[j]; c2r[j] = i;
1395d595f711SHong Zhang         }
1396d595f711SHong Zhang         i = nexti;
1397d595f711SHong Zhang       }
1398d595f711SHong Zhang 
1399d595f711SHong Zhang       /* copy data into U(k,:) */
1400d595f711SHong Zhang       rs   = 0.0;
1401d595f711SHong Zhang       jmin = bi[k]; jmax = bi[k+1]-1;
1402d595f711SHong Zhang       if (jmin < jmax) {
1403d595f711SHong Zhang         for (j=jmin; j<jmax; j++){
1404d595f711SHong Zhang           col = bj[j]; ba[j] = rtmp[col]; rs += PetscAbsScalar(ba[j]);
1405d595f711SHong Zhang         }
1406d595f711SHong Zhang         /* add the k-th row into il and c2r */
1407d595f711SHong Zhang         il[k] = jmin;
1408d595f711SHong Zhang         i = bj[jmin]; c2r[k] = c2r[i]; c2r[i] = k;
1409d595f711SHong Zhang       }
1410d595f711SHong Zhang 
1411d595f711SHong Zhang       sctx.rs  = rs;
1412d595f711SHong Zhang       sctx.pv  = dk;
1413d0660788SBarry Smith       ierr = MatPivotCheck(A,info,&sctx,k);CHKERRQ(ierr);
141407b50cabSHong Zhang       if(sctx.newshift) break;
1415d595f711SHong Zhang       dk = sctx.pv;
1416d595f711SHong Zhang 
1417d595f711SHong Zhang       ba[bdiag[k]] = 1.0/dk; /* U(k,k) */
1418d595f711SHong Zhang     }
141907b50cabSHong Zhang   } while (sctx.newshift);
1420d595f711SHong Zhang 
1421d595f711SHong Zhang   ierr = PetscFree3(rtmp,il,c2r);CHKERRQ(ierr);
1422d595f711SHong Zhang 
1423d595f711SHong Zhang   B->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
14249dcdb97aSHong Zhang   B->ops->solves          = MatSolves_SeqSBAIJ_1;
1425d595f711SHong Zhang   B->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
142680d3d5a7SHong Zhang   B->ops->forwardsolve    = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering;
142780d3d5a7SHong Zhang   B->ops->backwardsolve   = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering;
1428d595f711SHong Zhang 
14297b056e98SHong Zhang   B->assembled    = PETSC_TRUE;
14307b056e98SHong Zhang   B->preallocated = PETSC_TRUE;
14317b056e98SHong Zhang   ierr = PetscLogFlops(B->rmap->n);CHKERRQ(ierr);
1432d595f711SHong Zhang 
1433d595f711SHong Zhang   /* MatPivotView() */
1434d595f711SHong Zhang   if (sctx.nshift){
1435f4db908eSBarry Smith     if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
1436d595f711SHong Zhang       ierr = PetscInfo4(A,"number of shift_pd tries %D, shift_amount %G, diagonal shifted up by %e fraction top_value %e\n",sctx.nshift,sctx.shift_amount,sctx.shift_fraction,sctx.shift_top);CHKERRQ(ierr);
1437f4db908eSBarry Smith     } else if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
1438d595f711SHong Zhang       ierr = PetscInfo2(A,"number of shift_nz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
1439f4db908eSBarry Smith     } else if (info->shifttype == (PetscReal)MAT_SHIFT_INBLOCKS){
144014d2772eSHong Zhang       ierr = PetscInfo2(A,"number of shift_inblocks applied %D, each shift_amount %G\n",sctx.nshift,info->shiftamount);CHKERRQ(ierr);
1441d595f711SHong Zhang     }
1442d595f711SHong Zhang   }
1443d595f711SHong Zhang   PetscFunctionReturn(0);
1444d595f711SHong Zhang }
1445d595f711SHong Zhang 
1446d595f711SHong Zhang #undef __FUNCT__
1447d595f711SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering_inplace"
1448d595f711SHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering_inplace(Mat C,Mat A,const MatFactorInfo *info)
1449671cb588SHong Zhang {
1450671cb588SHong Zhang   Mat_SeqSBAIJ   *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ *)C->data;
1451dfbe8321SBarry Smith   PetscErrorCode ierr;
145213f74950SBarry Smith   PetscInt       i,j,mbs = a->mbs;
145313f74950SBarry Smith   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
14543cea4cbeSHong Zhang   PetscInt       k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz;
1455653a6975SHong Zhang   MatScalar      *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval;
145614d2772eSHong Zhang   PetscReal      rs;
14570e95ead3SHong Zhang   FactorShiftCtx sctx;
1458653a6975SHong Zhang 
1459653a6975SHong Zhang   PetscFunctionBegin;
14600e95ead3SHong Zhang   /* MatPivotSetUp(): initialize shift context sctx */
14610e95ead3SHong Zhang   ierr = PetscMemzero(&sctx,sizeof(FactorShiftCtx));CHKERRQ(ierr);
14620e95ead3SHong Zhang 
1463653a6975SHong Zhang   /* initialization */
1464653a6975SHong Zhang   /* il and jl record the first nonzero element in each row of the accessing
1465653a6975SHong Zhang      window U(0:k, k:mbs-1).
1466653a6975SHong Zhang      jl:    list of rows to be added to uneliminated rows
1467653a6975SHong Zhang             i>= k: jl(i) is the first row to be added to row i
1468653a6975SHong Zhang             i<  k: jl(i) is the row following row i in some list of rows
1469653a6975SHong Zhang             jl(i) = mbs indicates the end of a list
1470653a6975SHong Zhang      il(i): points to the first nonzero element in U(i,k:mbs-1)
1471653a6975SHong Zhang   */
1472653a6975SHong Zhang   ierr = PetscMalloc(mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
1473d8c74875SBarry Smith   ierr = PetscMalloc2(mbs,PetscInt,&il,mbs,PetscInt,&jl);CHKERRQ(ierr);
1474b00f7748SHong Zhang 
1475b00f7748SHong Zhang   do {
147607b50cabSHong Zhang     sctx.newshift = PETSC_FALSE;
1477653a6975SHong Zhang     for (i=0; i<mbs; i++) {
1478653a6975SHong Zhang       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
1479653a6975SHong Zhang     }
1480653a6975SHong Zhang 
1481997a0949SHong Zhang     for (k = 0; k<mbs; k++){
1482653a6975SHong Zhang       /*initialize k-th row with elements nonzero in row perm(k) of A */
1483653a6975SHong Zhang       nz   = ai[k+1] - ai[k];
1484653a6975SHong Zhang       acol = aj + ai[k];
1485653a6975SHong Zhang       aval = aa + ai[k];
1486653a6975SHong Zhang       bval = ba + bi[k];
1487653a6975SHong Zhang       while (nz -- ){
1488653a6975SHong Zhang         rtmp[*acol++] = *aval++;
1489653a6975SHong Zhang         *bval++       = 0.0; /* for in-place factorization */
1490653a6975SHong Zhang       }
14913cea4cbeSHong Zhang 
14923cea4cbeSHong Zhang       /* shift the diagonal of the matrix */
14933cea4cbeSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
1494653a6975SHong Zhang 
1495653a6975SHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
1496653a6975SHong Zhang       dk = rtmp[k];
1497653a6975SHong Zhang       i  = jl[k]; /* first row to be added to k_th row  */
1498653a6975SHong Zhang 
1499653a6975SHong Zhang       while (i < k){
1500653a6975SHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
1501653a6975SHong Zhang         /* compute multiplier, update D(k) and U(i,k) */
1502653a6975SHong Zhang         ili   = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
1503653a6975SHong Zhang         uikdi = - ba[ili]*ba[bi[i]];
1504653a6975SHong Zhang         dk   += uikdi*ba[ili];
1505653a6975SHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
1506653a6975SHong Zhang 
1507653a6975SHong Zhang         /* add multiple of row i to k-th row ... */
1508653a6975SHong Zhang         jmin = ili + 1;
1509653a6975SHong Zhang         nz   = bi[i+1] - jmin;
1510653a6975SHong Zhang         if (nz > 0){
1511653a6975SHong Zhang           bcol = bj + jmin;
1512653a6975SHong Zhang           bval = ba + jmin;
1513dc0b31edSSatish Balay           ierr = PetscLogFlops(2.0*nz);CHKERRQ(ierr);
1514653a6975SHong Zhang           while (nz --) rtmp[*bcol++] += uikdi*(*bval++);
1515187a9f4bSHong Zhang 
1516997a0949SHong Zhang           /* update il and jl for i-th row */
1517997a0949SHong Zhang           il[i] = jmin;
1518997a0949SHong Zhang           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
1519653a6975SHong Zhang         }
1520653a6975SHong Zhang         i = nexti;
1521653a6975SHong Zhang       }
1522653a6975SHong Zhang 
15233cea4cbeSHong Zhang       /* shift the diagonals when zero pivot is detected */
15243cea4cbeSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
15253cea4cbeSHong Zhang       rs   = 0.0;
152621c26570Svictorle       jmin = bi[k]+1;
152721c26570Svictorle       nz   = bi[k+1] - jmin;
152821c26570Svictorle       if (nz){
152921c26570Svictorle         bcol = bj + jmin;
153021c26570Svictorle         while (nz--){
153189f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
153289f845c8SHong Zhang           bcol++;
153321c26570Svictorle         }
153421c26570Svictorle       }
15353cea4cbeSHong Zhang 
15363cea4cbeSHong Zhang       sctx.rs = rs;
15373cea4cbeSHong Zhang       sctx.pv = dk;
1538d0660788SBarry Smith       ierr = MatPivotCheck(A,info,&sctx,k);CHKERRQ(ierr);
153907b50cabSHong Zhang       if (sctx.newshift) break;    /* sctx.shift_amount is updated */
15400e95ead3SHong Zhang       dk = sctx.pv;
1541653a6975SHong Zhang 
1542997a0949SHong Zhang       /* copy data into U(k,:) */
1543653a6975SHong Zhang       ba[bi[k]] = 1.0/dk;
1544653a6975SHong Zhang       jmin      = bi[k]+1;
1545653a6975SHong Zhang       nz        = bi[k+1] - jmin;
1546653a6975SHong Zhang       if (nz){
1547653a6975SHong Zhang         bcol = bj + jmin;
1548653a6975SHong Zhang         bval = ba + jmin;
1549653a6975SHong Zhang         while (nz--){
1550653a6975SHong Zhang           *bval++       = rtmp[*bcol];
1551653a6975SHong Zhang           rtmp[*bcol++] = 0.0;
1552653a6975SHong Zhang         }
1553997a0949SHong Zhang         /* add k-th row into il and jl */
1554653a6975SHong Zhang         il[k] = jmin;
1555997a0949SHong Zhang         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
1556653a6975SHong Zhang       }
1557b00f7748SHong Zhang     } /* end of for (k = 0; k<mbs; k++) */
155807b50cabSHong Zhang   } while (sctx.newshift);
1559653a6975SHong Zhang   ierr = PetscFree(rtmp);CHKERRQ(ierr);
1560d8c74875SBarry Smith   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
1561653a6975SHong Zhang 
15620a3c351aSHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
15634f79d315SHong Zhang   C->ops->solves         = MatSolves_SeqSBAIJ_1_inplace;
15640a3c351aSHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
15650a3c351aSHong Zhang   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
15660a3c351aSHong Zhang   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
1567db4efbfdSBarry Smith 
1568653a6975SHong Zhang   C->assembled    = PETSC_TRUE;
1569653a6975SHong Zhang   C->preallocated = PETSC_TRUE;
1570d0f46423SBarry Smith   ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr);
15713cea4cbeSHong Zhang   if (sctx.nshift){
1572f4db908eSBarry Smith     if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
15731e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
1574f4db908eSBarry Smith     } else if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
15751e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
1576b00f7748SHong Zhang     }
1577fb10cecfSBarry Smith   }
1578653a6975SHong Zhang   PetscFunctionReturn(0);
1579653a6975SHong Zhang }
1580653a6975SHong Zhang 
1581653a6975SHong Zhang #undef __FUNCT__
15824a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactor_SeqSBAIJ"
15830481f469SBarry Smith PetscErrorCode MatCholeskyFactor_SeqSBAIJ(Mat A,IS perm,const MatFactorInfo *info)
158449b5e25fSSatish Balay {
1585dfbe8321SBarry Smith   PetscErrorCode ierr;
158649b5e25fSSatish Balay   Mat            C;
158749b5e25fSSatish Balay 
158849b5e25fSSatish Balay   PetscFunctionBegin;
1589719d5645SBarry Smith   ierr = MatGetFactor(A,"petsc",MAT_FACTOR_CHOLESKY,&C);CHKERRQ(ierr);
1590719d5645SBarry Smith   ierr = MatCholeskyFactorSymbolic(C,A,perm,info);CHKERRQ(ierr);
1591719d5645SBarry Smith   ierr = MatCholeskyFactorNumeric(C,A,info);CHKERRQ(ierr);
1592db4efbfdSBarry Smith   A->ops->solve            = C->ops->solve;
1593db4efbfdSBarry Smith   A->ops->solvetranspose   = C->ops->solvetranspose;
1594eb6b5d47SBarry Smith   ierr = MatHeaderMerge(A,C);CHKERRQ(ierr);
159549b5e25fSSatish Balay   PetscFunctionReturn(0);
159649b5e25fSSatish Balay }
159749b5e25fSSatish Balay 
159849b5e25fSSatish Balay 
1599