xref: /petsc/src/mat/impls/sbaij/seq/sbaijfact.c (revision 7b6c816cf7b5355c1764dbf81d50649d87d3fa60)
149b5e25fSSatish Balay 
2c6db04a5SJed Brown #include <../src/mat/impls/baij/seq/baij.h>
3c6db04a5SJed Brown #include <../src/mat/impls/sbaij/seq/sbaij.h>
4af0996ceSBarry Smith #include <petsc/private/kernels/blockinvert.h>
5c6db04a5SJed Brown #include <petscis.h>
649b5e25fSSatish Balay 
75f9f512dSHong Zhang #undef __FUNCT__
85f9f512dSHong Zhang #define __FUNCT__ "MatGetInertia_SeqSBAIJ"
94ff4e45dSHong Zhang PetscErrorCode MatGetInertia_SeqSBAIJ(Mat F,PetscInt *nneg,PetscInt *nzero,PetscInt *npos)
105f9f512dSHong Zhang {
114ff4e45dSHong Zhang   Mat_SeqSBAIJ *fact=(Mat_SeqSBAIJ*)F->data;
124ff4e45dSHong Zhang   MatScalar    *dd=fact->a;
134ff4e45dSHong Zhang   PetscInt     mbs=fact->mbs,bs=F->rmap->bs,i,nneg_tmp,npos_tmp,*fi=fact->diag;
145f9f512dSHong Zhang 
155f9f512dSHong Zhang   PetscFunctionBegin;
16e32f2f54SBarry Smith   if (bs != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for bs: %D >1 yet",bs);
17*7b6c816cSBarry Smith   if (F->factorerrortype==MAT_FACTOR_NUMERIC_ZEROPIVOT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatFactor fails with numeric zeropivot");
184ff4e45dSHong Zhang 
194ff4e45dSHong Zhang   nneg_tmp = 0; npos_tmp = 0;
20eeeff2ecSHong Zhang   for (i=0; i<mbs; i++) {
2126fbe8dcSKarl Rupp     if (PetscRealPart(dd[*fi]) > 0.0) npos_tmp++;
224ff4e45dSHong Zhang     else if (PetscRealPart(dd[*fi]) < 0.0) nneg_tmp++;
23eeeff2ecSHong Zhang     fi++;
243e0d88b5SBarry Smith   }
254ff4e45dSHong Zhang   if (nneg)  *nneg  = nneg_tmp;
26eeeff2ecSHong Zhang   if (npos)  *npos  = npos_tmp;
274ff4e45dSHong Zhang   if (nzero) *nzero = mbs - nneg_tmp - npos_tmp;
285f9f512dSHong Zhang   PetscFunctionReturn(0);
295f9f512dSHong Zhang }
305f9f512dSHong Zhang 
315f9f512dSHong Zhang /*
325f9f512dSHong Zhang   Symbolic U^T*D*U factorization for SBAIJ format. Modified from SSF of YSMP.
3310c27e3dSHong Zhang   Use Modified Sparse Row (MSR) storage for u and ju. See page 85, "Iterative Methods ..." by Saad.
345f9f512dSHong Zhang */
3510c27e3dSHong Zhang #undef __FUNCT__
3610c27e3dSHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqSBAIJ_MSR"
370481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ_MSR(Mat F,Mat A,IS perm,const MatFactorInfo *info)
3810c27e3dSHong Zhang {
3910c27e3dSHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b;
4010c27e3dSHong Zhang   PetscErrorCode ierr;
415d0c19d7SBarry Smith   const PetscInt *rip,*ai,*aj;
425d0c19d7SBarry Smith   PetscInt       i,mbs = a->mbs,*jutmp,bs = A->rmap->bs,bs2=a->bs2;
4310c27e3dSHong Zhang   PetscInt       m,reallocs = 0,prow;
4410c27e3dSHong Zhang   PetscInt       *jl,*q,jmin,jmax,juidx,nzk,qm,*iu,*ju,k,j,vj,umax,maxadd;
4510c27e3dSHong Zhang   PetscReal      f = info->fill;
46ace3abfcSBarry Smith   PetscBool      perm_identity;
4710c27e3dSHong Zhang 
4810c27e3dSHong Zhang   PetscFunctionBegin;
4910c27e3dSHong Zhang   /* check whether perm is the identity mapping */
5010c27e3dSHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
5110c27e3dSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
5210c27e3dSHong Zhang 
5310c27e3dSHong Zhang   if (perm_identity) { /* without permutation */
5410c27e3dSHong Zhang     a->permute = PETSC_FALSE;
5526fbe8dcSKarl Rupp 
5610c27e3dSHong Zhang     ai = a->i; aj = a->j;
5710c27e3dSHong Zhang   } else {            /* non-trivial permutation */
5810c27e3dSHong Zhang     a->permute = PETSC_TRUE;
5926fbe8dcSKarl Rupp 
6010c27e3dSHong Zhang     ierr = MatReorderingSeqSBAIJ(A,perm);CHKERRQ(ierr);
6126fbe8dcSKarl Rupp 
6210c27e3dSHong Zhang     ai = a->inew; aj = a->jnew;
6310c27e3dSHong Zhang   }
6410c27e3dSHong Zhang 
6510c27e3dSHong Zhang   /* initialization */
66854ce69bSBarry Smith   ierr  = PetscMalloc1(mbs+1,&iu);CHKERRQ(ierr);
6710c27e3dSHong Zhang   umax  = (PetscInt)(f*ai[mbs] + 1); umax += mbs + 1;
68785e854fSJed Brown   ierr  = PetscMalloc1(umax,&ju);CHKERRQ(ierr);
6910c27e3dSHong Zhang   iu[0] = mbs+1;
7010c27e3dSHong Zhang   juidx = mbs + 1; /* index for ju */
71d8c74875SBarry Smith   /* jl linked list for pivot row -- linked list for col index */
72dcca6d9dSJed Brown   ierr = PetscMalloc2(mbs,&jl,mbs,&q);CHKERRQ(ierr);
7310c27e3dSHong Zhang   for (i=0; i<mbs; i++) {
7410c27e3dSHong Zhang     jl[i] = mbs;
7510c27e3dSHong Zhang     q[i]  = 0;
7610c27e3dSHong Zhang   }
7710c27e3dSHong Zhang 
7810c27e3dSHong Zhang   /* for each row k */
7910c27e3dSHong Zhang   for (k=0; k<mbs; k++) {
8010c27e3dSHong Zhang     for (i=0; i<mbs; i++) q[i] = 0;  /* to be removed! */
8110c27e3dSHong Zhang     nzk  = 0; /* num. of nz blocks in k-th block row with diagonal block excluded */
8210c27e3dSHong Zhang     q[k] = mbs;
8310c27e3dSHong Zhang     /* initialize nonzero structure of k-th row to row rip[k] of A */
8410c27e3dSHong Zhang     jmin = ai[rip[k]] +1; /* exclude diag[k] */
8510c27e3dSHong Zhang     jmax = ai[rip[k]+1];
8610c27e3dSHong Zhang     for (j=jmin; j<jmax; j++) {
8710c27e3dSHong Zhang       vj = rip[aj[j]]; /* col. value */
8810c27e3dSHong Zhang       if (vj > k) {
8910c27e3dSHong Zhang         qm = k;
9010c27e3dSHong Zhang         do {
9110c27e3dSHong Zhang           m = qm; qm = q[m];
9210c27e3dSHong Zhang         } while (qm < vj);
93535b19f3SBarry Smith         if (qm == vj) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Duplicate entry in A\n");
9410c27e3dSHong Zhang         nzk++;
9510c27e3dSHong Zhang         q[m]  = vj;
9610c27e3dSHong Zhang         q[vj] = qm;
9710c27e3dSHong Zhang       } /* if (vj > k) */
9810c27e3dSHong Zhang     } /* for (j=jmin; j<jmax; j++) */
9910c27e3dSHong Zhang 
10010c27e3dSHong Zhang     /* modify nonzero structure of k-th row by computing fill-in
10110c27e3dSHong Zhang        for each row i to be merged in */
10210c27e3dSHong Zhang     prow = k;
10310c27e3dSHong Zhang     prow = jl[prow]; /* next pivot row (== mbs for symbolic factorization) */
10410c27e3dSHong Zhang 
10510c27e3dSHong Zhang     while (prow < k) {
10610c27e3dSHong Zhang       /* merge row prow into k-th row */
10710c27e3dSHong Zhang       jmin = iu[prow] + 1; jmax = iu[prow+1];
10810c27e3dSHong Zhang       qm   = k;
10910c27e3dSHong Zhang       for (j=jmin; j<jmax; j++) {
11010c27e3dSHong Zhang         vj = ju[j];
11110c27e3dSHong Zhang         do {
11210c27e3dSHong Zhang           m = qm; qm = q[m];
11310c27e3dSHong Zhang         } while (qm < vj);
11410c27e3dSHong Zhang         if (qm != vj) {
11510c27e3dSHong Zhang           nzk++; q[m] = vj; q[vj] = qm; qm = vj;
11610c27e3dSHong Zhang         }
11710c27e3dSHong Zhang       }
11810c27e3dSHong Zhang       prow = jl[prow]; /* next pivot row */
11910c27e3dSHong Zhang     }
12010c27e3dSHong Zhang 
12110c27e3dSHong Zhang     /* add k to row list for first nonzero element in k-th row */
12210c27e3dSHong Zhang     if (nzk > 0) {
12310c27e3dSHong Zhang       i     = q[k]; /* col value of first nonzero element in U(k, k+1:mbs-1) */
12410c27e3dSHong Zhang       jl[k] = jl[i]; jl[i] = k;
12510c27e3dSHong Zhang     }
12610c27e3dSHong Zhang     iu[k+1] = iu[k] + nzk;
12710c27e3dSHong Zhang 
12810c27e3dSHong Zhang     /* allocate more space to ju if needed */
12910c27e3dSHong Zhang     if (iu[k+1] > umax) {
13010c27e3dSHong Zhang       /* estimate how much additional space we will need */
13110c27e3dSHong Zhang       /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */
13210c27e3dSHong Zhang       /* just double the memory each time */
13310c27e3dSHong Zhang       maxadd = umax;
13410c27e3dSHong Zhang       if (maxadd < nzk) maxadd = (mbs-k)*(nzk+1)/2;
13510c27e3dSHong Zhang       umax += maxadd;
13610c27e3dSHong Zhang 
13710c27e3dSHong Zhang       /* allocate a longer ju */
138785e854fSJed Brown       ierr = PetscMalloc1(umax,&jutmp);CHKERRQ(ierr);
13910c27e3dSHong Zhang       ierr = PetscMemcpy(jutmp,ju,iu[k]*sizeof(PetscInt));CHKERRQ(ierr);
14010c27e3dSHong Zhang       ierr = PetscFree(ju);CHKERRQ(ierr);
14110c27e3dSHong Zhang       ju   = jutmp;
14210c27e3dSHong Zhang       reallocs++; /* count how many times we realloc */
14310c27e3dSHong Zhang     }
14410c27e3dSHong Zhang 
14510c27e3dSHong Zhang     /* save nonzero structure of k-th row in ju */
14610c27e3dSHong Zhang     i=k;
14710c27e3dSHong Zhang     while (nzk--) {
14810c27e3dSHong Zhang       i           = q[i];
14910c27e3dSHong Zhang       ju[juidx++] = i;
15010c27e3dSHong Zhang     }
15110c27e3dSHong Zhang   }
15210c27e3dSHong Zhang 
1536cf91177SBarry Smith #if defined(PETSC_USE_INFO)
15410c27e3dSHong Zhang   if (ai[mbs] != 0) {
15510c27e3dSHong Zhang     PetscReal af = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]);
15657622a8eSBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)f,(double)af);CHKERRQ(ierr);
15757622a8eSBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);CHKERRQ(ierr);
15857622a8eSBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%g);\n",(double)af);CHKERRQ(ierr);
159ae15b995SBarry Smith     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
16010c27e3dSHong Zhang   } else {
161ae15b995SBarry Smith     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
16210c27e3dSHong Zhang   }
16363ba0a88SBarry Smith #endif
16410c27e3dSHong Zhang 
16510c27e3dSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
166d8c74875SBarry Smith   ierr = PetscFree2(jl,q);CHKERRQ(ierr);
16710c27e3dSHong Zhang 
16810c27e3dSHong Zhang   /* put together the new matrix */
169367daffbSBarry Smith   ierr = MatSeqSBAIJSetPreallocation(F,bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr);
17010c27e3dSHong Zhang 
1713bb1ff40SBarry Smith   /* ierr = PetscLogObjectParent((PetscObject)B,(PetscObject)iperm);CHKERRQ(ierr); */
172719d5645SBarry Smith   b                = (Mat_SeqSBAIJ*)(F)->data;
17310c27e3dSHong Zhang   b->singlemalloc  = PETSC_FALSE;
174e6b907acSBarry Smith   b->free_a        = PETSC_TRUE;
175e6b907acSBarry Smith   b->free_ij       = PETSC_TRUE;
17626fbe8dcSKarl Rupp 
177854ce69bSBarry Smith   ierr    = PetscMalloc1((iu[mbs]+1)*bs2,&b->a);CHKERRQ(ierr);
17810c27e3dSHong Zhang   b->j    = ju;
17910c27e3dSHong Zhang   b->i    = iu;
18010c27e3dSHong Zhang   b->diag = 0;
18110c27e3dSHong Zhang   b->ilen = 0;
18210c27e3dSHong Zhang   b->imax = 0;
18310c27e3dSHong Zhang   b->row  = perm;
18426fbe8dcSKarl Rupp 
18510c27e3dSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
18626fbe8dcSKarl Rupp 
18710c27e3dSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
18826fbe8dcSKarl Rupp 
18910c27e3dSHong Zhang   b->icol = perm;
19010c27e3dSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
191854ce69bSBarry Smith   ierr    = PetscMalloc1(bs*mbs+bs,&b->solve_work);CHKERRQ(ierr);
19210c27e3dSHong Zhang   /* In b structure:  Free imax, ilen, old a, old j.
19310c27e3dSHong Zhang      Allocate idnew, solve_work, new a, new j */
1943bb1ff40SBarry Smith   ierr     = PetscLogObjectMemory((PetscObject)F,(iu[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
19510c27e3dSHong Zhang   b->maxnz = b->nz = iu[mbs];
19610c27e3dSHong Zhang 
197719d5645SBarry Smith   (F)->info.factor_mallocs   = reallocs;
198719d5645SBarry Smith   (F)->info.fill_ratio_given = f;
19910c27e3dSHong Zhang   if (ai[mbs] != 0) {
200719d5645SBarry Smith     (F)->info.fill_ratio_needed = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]);
20110c27e3dSHong Zhang   } else {
202719d5645SBarry Smith     (F)->info.fill_ratio_needed = 0.0;
20310c27e3dSHong Zhang   }
204d595f711SHong Zhang   ierr = MatSeqSBAIJSetNumericFactorization_inplace(F,perm_identity);CHKERRQ(ierr);
20510c27e3dSHong Zhang   PetscFunctionReturn(0);
20610c27e3dSHong Zhang }
20710c27e3dSHong Zhang /*
20810c27e3dSHong Zhang     Symbolic U^T*D*U factorization for SBAIJ format.
20980d3d5a7SHong Zhang     See MatICCFactorSymbolic_SeqAIJ() for description of its data structure.
21010c27e3dSHong Zhang */
211c6db04a5SJed Brown #include <petscbt.h>
212c6db04a5SJed Brown #include <../src/mat/utils/freespace.h>
2134a2ae208SSatish Balay #undef __FUNCT__
2144a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqSBAIJ"
2150481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
21649b5e25fSSatish Balay {
21798a8e78dSHong Zhang   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data;
21898a8e78dSHong Zhang   Mat_SeqSBAIJ       *b;
2196849ba73SBarry Smith   PetscErrorCode     ierr;
220ace3abfcSBarry Smith   PetscBool          perm_identity,missing;
22198a8e78dSHong Zhang   PetscReal          fill = info->fill;
2227b056e98SHong Zhang   const PetscInt     *rip,*ai=a->i,*aj=a->j;
2239186b105SHong Zhang   PetscInt           i,mbs=a->mbs,bs=A->rmap->bs,reallocs=0,prow;
22498a8e78dSHong Zhang   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
225c6d0d4f0SHong Zhang   PetscInt           nlnk,*lnk,ncols,*cols,*uj,**ui_ptr,*uj_ptr,*udiag;
2260298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
22798a8e78dSHong Zhang   PetscBT            lnkbt;
228d595f711SHong Zhang 
229d595f711SHong Zhang   PetscFunctionBegin;
2309186b105SHong Zhang   if (A->rmap->n != A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Must be square matrix, rows %D columns %D",A->rmap->n,A->cmap->n);
2319186b105SHong Zhang   ierr = MatMissingDiagonal(A,&missing,&i);CHKERRQ(ierr);
2329186b105SHong Zhang   if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i);
2336d07c2b1SHong Zhang   if (bs > 1) {
2346d07c2b1SHong Zhang     ierr = MatCholeskyFactorSymbolic_SeqSBAIJ_inplace(fact,A,perm,info);CHKERRQ(ierr);
2356d07c2b1SHong Zhang     PetscFunctionReturn(0);
2366d07c2b1SHong Zhang   }
237d595f711SHong Zhang 
238d595f711SHong Zhang   /* check whether perm is the identity mapping */
239d595f711SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
240e32f2f54SBarry Smith   if (!perm_identity) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix reordering is not supported for sbaij matrix. Use aij format");
241d595f711SHong Zhang   a->permute = PETSC_FALSE;
242d595f711SHong Zhang   ierr       = ISGetIndices(perm,&rip);CHKERRQ(ierr);
243d595f711SHong Zhang 
244d595f711SHong Zhang   /* initialization */
245854ce69bSBarry Smith   ierr  = PetscMalloc1(mbs+1,&ui);CHKERRQ(ierr);
246854ce69bSBarry Smith   ierr  = PetscMalloc1(mbs+1,&udiag);CHKERRQ(ierr);
247d595f711SHong Zhang   ui[0] = 0;
248d595f711SHong Zhang 
249d595f711SHong Zhang   /* jl: linked list for storing indices of the pivot rows
250d595f711SHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
251dcca6d9dSJed Brown   ierr = PetscMalloc4(mbs,&ui_ptr,mbs,&il,mbs,&jl,mbs,&cols);CHKERRQ(ierr);
252d595f711SHong Zhang   for (i=0; i<mbs; i++) {
253d595f711SHong Zhang     jl[i] = mbs; il[i] = 0;
254d595f711SHong Zhang   }
255d595f711SHong Zhang 
256d595f711SHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
257d595f711SHong Zhang   nlnk = mbs + 1;
258d595f711SHong Zhang   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
259d595f711SHong Zhang 
260d595f711SHong Zhang   /* initial FreeSpace size is fill*(ai[mbs]+1) */
261f91af8c7SBarry Smith   ierr          = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,ai[mbs]+1),&free_space);CHKERRQ(ierr);
262d595f711SHong Zhang   current_space = free_space;
263d595f711SHong Zhang 
264d595f711SHong Zhang   for (k=0; k<mbs; k++) {  /* for each active row k */
265d595f711SHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
266d595f711SHong Zhang     nzk   = 0;
267c6d0d4f0SHong Zhang     ncols = ai[k+1] - ai[k];
268e32f2f54SBarry Smith     if (!ncols) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_CH_ZRPVT,"Empty row %D in matrix ",k);
269d595f711SHong Zhang     for (j=0; j<ncols; j++) {
270c6d0d4f0SHong Zhang       i       = *(aj + ai[k] + j);
271c6d0d4f0SHong Zhang       cols[j] = i;
272d595f711SHong Zhang     }
273d595f711SHong Zhang     ierr = PetscLLAdd(ncols,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
274d595f711SHong Zhang     nzk += nlnk;
275d595f711SHong Zhang 
276d595f711SHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
277d595f711SHong Zhang     prow = jl[k]; /* 1st pivot row */
278d595f711SHong Zhang 
279d595f711SHong Zhang     while (prow < k) {
280d595f711SHong Zhang       nextprow = jl[prow];
281d595f711SHong Zhang       /* merge prow into k-th row */
282d595f711SHong Zhang       jmin   = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
283d595f711SHong Zhang       jmax   = ui[prow+1];
284d595f711SHong Zhang       ncols  = jmax-jmin;
285d595f711SHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
286d595f711SHong Zhang       ierr   = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
287d595f711SHong Zhang       nzk   += nlnk;
288d595f711SHong Zhang 
289d595f711SHong Zhang       /* update il and jl for prow */
290d595f711SHong Zhang       if (jmin < jmax) {
291d595f711SHong Zhang         il[prow] = jmin;
292d595f711SHong Zhang         j        = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
293d595f711SHong Zhang       }
294d595f711SHong Zhang       prow = nextprow;
295d595f711SHong Zhang     }
296d595f711SHong Zhang 
297d595f711SHong Zhang     /* if free space is not available, make more free space */
298d595f711SHong Zhang     if (current_space->local_remaining<nzk) {
299d595f711SHong Zhang       i    = mbs - k + 1; /* num of unfactored rows */
300f91af8c7SBarry Smith       i    = PetscIntMultTruncate(i,PetscMin(nzk, i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
301d595f711SHong Zhang       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
302d595f711SHong Zhang       reallocs++;
303d595f711SHong Zhang     }
304d595f711SHong Zhang 
305d595f711SHong Zhang     /* copy data into free space, then initialize lnk */
306d595f711SHong Zhang     ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
307d595f711SHong Zhang 
308d595f711SHong Zhang     /* add the k-th row into il and jl */
3097b056e98SHong Zhang     if (nzk > 1) {
310d595f711SHong Zhang       i     = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
311d595f711SHong Zhang       jl[k] = jl[i]; jl[i] = k;
312d595f711SHong Zhang       il[k] = ui[k] + 1;
313d595f711SHong Zhang     }
314d595f711SHong Zhang     ui_ptr[k] = current_space->array;
31526fbe8dcSKarl Rupp 
316d595f711SHong Zhang     current_space->array           += nzk;
317d595f711SHong Zhang     current_space->local_used      += nzk;
318d595f711SHong Zhang     current_space->local_remaining -= nzk;
319d595f711SHong Zhang 
320d595f711SHong Zhang     ui[k+1] = ui[k] + nzk;
321d595f711SHong Zhang   }
322d595f711SHong Zhang 
323d595f711SHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
324d595f711SHong Zhang   ierr = PetscFree4(ui_ptr,il,jl,cols);CHKERRQ(ierr);
325d595f711SHong Zhang 
326d595f711SHong Zhang   /* destroy list of free space and other temporary array(s) */
327854ce69bSBarry Smith   ierr = PetscMalloc1(ui[mbs]+1,&uj);CHKERRQ(ierr);
3287b056e98SHong Zhang   ierr = PetscFreeSpaceContiguous_Cholesky(&free_space,uj,mbs,ui,udiag);CHKERRQ(ierr); /* store matrix factor */
329d595f711SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
330d595f711SHong Zhang 
331d595f711SHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
332367daffbSBarry Smith   ierr = MatSeqSBAIJSetPreallocation(fact,bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr);
333d595f711SHong Zhang 
3347b056e98SHong Zhang   b               = (Mat_SeqSBAIJ*)fact->data;
335d595f711SHong Zhang   b->singlemalloc = PETSC_FALSE;
336d595f711SHong Zhang   b->free_a       = PETSC_TRUE;
337d595f711SHong Zhang   b->free_ij      = PETSC_TRUE;
33826fbe8dcSKarl Rupp 
339854ce69bSBarry Smith   ierr = PetscMalloc1(ui[mbs]+1,&b->a);CHKERRQ(ierr);
34026fbe8dcSKarl Rupp 
341d595f711SHong Zhang   b->j         = uj;
342d595f711SHong Zhang   b->i         = ui;
343c6d0d4f0SHong Zhang   b->diag      = udiag;
344c6d0d4f0SHong Zhang   b->free_diag = PETSC_TRUE;
345d595f711SHong Zhang   b->ilen      = 0;
346d595f711SHong Zhang   b->imax      = 0;
347d595f711SHong Zhang   b->row       = perm;
348d595f711SHong Zhang   b->icol      = perm;
34926fbe8dcSKarl Rupp 
350d595f711SHong Zhang   ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
3517b056e98SHong Zhang   ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
35226fbe8dcSKarl Rupp 
3537b056e98SHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
35426fbe8dcSKarl Rupp 
355854ce69bSBarry Smith   ierr = PetscMalloc1(mbs+1,&b->solve_work);CHKERRQ(ierr);
3563bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)fact,ui[mbs]*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
35726fbe8dcSKarl Rupp 
358d595f711SHong Zhang   b->maxnz = b->nz = ui[mbs];
359d595f711SHong Zhang 
36095b5ac22SHong Zhang   fact->info.factor_mallocs   = reallocs;
36195b5ac22SHong Zhang   fact->info.fill_ratio_given = fill;
362d595f711SHong Zhang   if (ai[mbs] != 0) {
36395b5ac22SHong Zhang     fact->info.fill_ratio_needed = ((PetscReal)ui[mbs])/ai[mbs];
364d595f711SHong Zhang   } else {
36595b5ac22SHong Zhang     fact->info.fill_ratio_needed = 0.0;
366d595f711SHong Zhang   }
36795b5ac22SHong Zhang #if defined(PETSC_USE_INFO)
36895b5ac22SHong Zhang   if (ai[mbs] != 0) {
36995b5ac22SHong Zhang     PetscReal af = fact->info.fill_ratio_needed;
37057622a8eSBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)fill,(double)af);CHKERRQ(ierr);
37157622a8eSBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);CHKERRQ(ierr);
37257622a8eSBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%g) for best performance.\n",(double)af);CHKERRQ(ierr);
37395b5ac22SHong Zhang   } else {
37495b5ac22SHong Zhang     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
37595b5ac22SHong Zhang   }
37695b5ac22SHong Zhang #endif
377c6d0d4f0SHong Zhang   fact->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering;
378d595f711SHong Zhang   PetscFunctionReturn(0);
379d595f711SHong Zhang }
380d595f711SHong Zhang 
381d595f711SHong Zhang #undef __FUNCT__
382d595f711SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqSBAIJ_inplace"
383d595f711SHong Zhang PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ_inplace(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
384d595f711SHong Zhang {
385d595f711SHong Zhang   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data;
386d595f711SHong Zhang   Mat_SeqSBAIJ       *b;
387d595f711SHong Zhang   PetscErrorCode     ierr;
388ace3abfcSBarry Smith   PetscBool          perm_identity,missing;
389d595f711SHong Zhang   PetscReal          fill = info->fill;
390d595f711SHong Zhang   const PetscInt     *rip,*ai,*aj;
391d595f711SHong Zhang   PetscInt           i,mbs=a->mbs,bs=A->rmap->bs,reallocs=0,prow,d;
392d595f711SHong Zhang   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
393d595f711SHong Zhang   PetscInt           nlnk,*lnk,ncols,*cols,*uj,**ui_ptr,*uj_ptr;
3940298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
395d595f711SHong Zhang   PetscBT            lnkbt;
39649b5e25fSSatish Balay 
39749b5e25fSSatish Balay   PetscFunctionBegin;
39858ebbce7SBarry Smith   ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr);
399e32f2f54SBarry Smith   if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
40058ebbce7SBarry Smith 
40110c27e3dSHong Zhang   /*
40210c27e3dSHong Zhang    This code originally uses Modified Sparse Row (MSR) storage
40310c27e3dSHong Zhang    (see page 85, "Iterative Methods ..." by Saad) for the output matrix B - bad choise!
40410c27e3dSHong Zhang    Then it is rewritten so the factor B takes seqsbaij format. However the associated
40510c27e3dSHong Zhang    MatCholeskyFactorNumeric_() have not been modified for the cases of bs>1 or !perm_identity,
40610c27e3dSHong Zhang    thus the original code in MSR format is still used for these cases.
40710c27e3dSHong Zhang    The code below should replace MatCholeskyFactorSymbolic_SeqSBAIJ_MSR() whenever
40810c27e3dSHong Zhang    MatCholeskyFactorNumeric_() is modified for using sbaij symbolic factor.
40910c27e3dSHong Zhang   */
410fff829cfSHong Zhang   if (bs > 1) {
411719d5645SBarry Smith     ierr = MatCholeskyFactorSymbolic_SeqSBAIJ_MSR(fact,A,perm,info);CHKERRQ(ierr);
412fff829cfSHong Zhang     PetscFunctionReturn(0);
413fff829cfSHong Zhang   }
41410c27e3dSHong Zhang 
41598a8e78dSHong Zhang   /* check whether perm is the identity mapping */
41698a8e78dSHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
41798a8e78dSHong Zhang 
418fff829cfSHong Zhang   if (perm_identity) {
419fff829cfSHong Zhang     a->permute = PETSC_FALSE;
42026fbe8dcSKarl Rupp 
421fff829cfSHong Zhang     ai = a->i; aj = a->j;
4226c4ed002SBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix reordering is not supported for sbaij matrix. Use aij format");
423fff829cfSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
424fff829cfSHong Zhang 
425fff829cfSHong Zhang   /* initialization */
426854ce69bSBarry Smith   ierr  = PetscMalloc1(mbs+1,&ui);CHKERRQ(ierr);
42798a8e78dSHong Zhang   ui[0] = 0;
42898a8e78dSHong Zhang 
42998a8e78dSHong Zhang   /* jl: linked list for storing indices of the pivot rows
4307625dc9aSHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
431dcca6d9dSJed Brown   ierr = PetscMalloc4(mbs,&ui_ptr,mbs,&il,mbs,&jl,mbs,&cols);CHKERRQ(ierr);
4327625dc9aSHong Zhang   for (i=0; i<mbs; i++) {
4337625dc9aSHong Zhang     jl[i] = mbs; il[i] = 0;
434fff829cfSHong Zhang   }
435fff829cfSHong Zhang 
43698a8e78dSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
4377625dc9aSHong Zhang   nlnk = mbs + 1;
4387625dc9aSHong Zhang   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
439fff829cfSHong Zhang 
4407625dc9aSHong Zhang   /* initial FreeSpace size is fill*(ai[mbs]+1) */
441f91af8c7SBarry Smith   ierr          = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,ai[mbs]+1),&free_space);CHKERRQ(ierr);
44298a8e78dSHong Zhang   current_space = free_space;
44398a8e78dSHong Zhang 
4447625dc9aSHong Zhang   for (k=0; k<mbs; k++) {  /* for each active row k */
44598a8e78dSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
44698a8e78dSHong Zhang     nzk   = 0;
44798a8e78dSHong Zhang     ncols = ai[rip[k]+1] - ai[rip[k]];
4487625dc9aSHong Zhang     for (j=0; j<ncols; j++) {
4497625dc9aSHong Zhang       i       = *(aj + ai[rip[k]] + j);
4507625dc9aSHong Zhang       cols[j] = rip[i];
4517625dc9aSHong Zhang     }
4527625dc9aSHong Zhang     ierr = PetscLLAdd(ncols,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
45398a8e78dSHong Zhang     nzk += nlnk;
45498a8e78dSHong Zhang 
45598a8e78dSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
45698a8e78dSHong Zhang     prow = jl[k]; /* 1st pivot row */
457fff829cfSHong Zhang 
458fff829cfSHong Zhang     while (prow < k) {
459fff829cfSHong Zhang       nextprow = jl[prow];
46098a8e78dSHong Zhang       /* merge prow into k-th row */
4617625dc9aSHong Zhang       jmin   = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
46298a8e78dSHong Zhang       jmax   = ui[prow+1];
46398a8e78dSHong Zhang       ncols  = jmax-jmin;
4647625dc9aSHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
4655a8e39fbSHong Zhang       ierr   = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
46698a8e78dSHong Zhang       nzk   += nlnk;
467fff829cfSHong Zhang 
46898a8e78dSHong Zhang       /* update il and jl for prow */
469fff829cfSHong Zhang       if (jmin < jmax) {
470fff829cfSHong Zhang         il[prow] = jmin;
47126fbe8dcSKarl Rupp 
4727625dc9aSHong Zhang         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
473fff829cfSHong Zhang       }
474fff829cfSHong Zhang       prow = nextprow;
475fff829cfSHong Zhang     }
476fff829cfSHong Zhang 
47798a8e78dSHong Zhang     /* if free space is not available, make more free space */
47898a8e78dSHong Zhang     if (current_space->local_remaining<nzk) {
4797625dc9aSHong Zhang       i    = mbs - k + 1; /* num of unfactored rows */
480f91af8c7SBarry Smith       i    = PetscMin(PetscIntMultTruncate(i,nzk), PetscIntMultTruncate(i,i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
481a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
48298a8e78dSHong Zhang       reallocs++;
48398a8e78dSHong Zhang     }
48498a8e78dSHong Zhang 
48598a8e78dSHong Zhang     /* copy data into free space, then initialize lnk */
4867625dc9aSHong Zhang     ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
48798a8e78dSHong Zhang 
48898a8e78dSHong Zhang     /* add the k-th row into il and jl */
48998a8e78dSHong Zhang     if (nzk-1 > 0) {
4907625dc9aSHong Zhang       i     = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
491fff829cfSHong Zhang       jl[k] = jl[i]; jl[i] = k;
49298a8e78dSHong Zhang       il[k] = ui[k] + 1;
493fff829cfSHong Zhang     }
4947625dc9aSHong Zhang     ui_ptr[k] = current_space->array;
49526fbe8dcSKarl Rupp 
49698a8e78dSHong Zhang     current_space->array           += nzk;
49798a8e78dSHong Zhang     current_space->local_used      += nzk;
49898a8e78dSHong Zhang     current_space->local_remaining -= nzk;
499fff829cfSHong Zhang 
50098a8e78dSHong Zhang     ui[k+1] = ui[k] + nzk;
501fff829cfSHong Zhang   }
502fff829cfSHong Zhang 
503fff829cfSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
504d8c74875SBarry Smith   ierr = PetscFree4(ui_ptr,il,jl,cols);CHKERRQ(ierr);
505fff829cfSHong Zhang 
50698a8e78dSHong Zhang   /* destroy list of free space and other temporary array(s) */
507854ce69bSBarry Smith   ierr = PetscMalloc1(ui[mbs]+1,&uj);CHKERRQ(ierr);
508a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
50998a8e78dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
510fff829cfSHong Zhang 
51198a8e78dSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
512367daffbSBarry Smith   ierr = MatSeqSBAIJSetPreallocation(fact,bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr);
51398a8e78dSHong Zhang 
51495b5ac22SHong Zhang   b               = (Mat_SeqSBAIJ*)fact->data;
515fff829cfSHong Zhang   b->singlemalloc = PETSC_FALSE;
516e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
517e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
51826fbe8dcSKarl Rupp 
519854ce69bSBarry Smith   ierr = PetscMalloc1(ui[mbs]+1,&b->a);CHKERRQ(ierr);
52026fbe8dcSKarl Rupp 
52198a8e78dSHong Zhang   b->j    = uj;
52298a8e78dSHong Zhang   b->i    = ui;
523fff829cfSHong Zhang   b->diag = 0;
524fff829cfSHong Zhang   b->ilen = 0;
525fff829cfSHong Zhang   b->imax = 0;
526fff829cfSHong Zhang   b->row  = perm;
52726fbe8dcSKarl Rupp 
528fff829cfSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
52926fbe8dcSKarl Rupp 
530fff829cfSHong Zhang   ierr     = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
531fff829cfSHong Zhang   b->icol  = perm;
532fff829cfSHong Zhang   ierr     = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
533854ce69bSBarry Smith   ierr     = PetscMalloc1(mbs+1,&b->solve_work);CHKERRQ(ierr);
5343bb1ff40SBarry Smith   ierr     = PetscLogObjectMemory((PetscObject)fact,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
5357625dc9aSHong Zhang   b->maxnz = b->nz = ui[mbs];
536fff829cfSHong Zhang 
53795b5ac22SHong Zhang   fact->info.factor_mallocs   = reallocs;
53895b5ac22SHong Zhang   fact->info.fill_ratio_given = fill;
5397625dc9aSHong Zhang   if (ai[mbs] != 0) {
54095b5ac22SHong Zhang     fact->info.fill_ratio_needed = ((PetscReal)ui[mbs])/ai[mbs];
541fff829cfSHong Zhang   } else {
54295b5ac22SHong Zhang     fact->info.fill_ratio_needed = 0.0;
543fff829cfSHong Zhang   }
54495b5ac22SHong Zhang #if defined(PETSC_USE_INFO)
54595b5ac22SHong Zhang   if (ai[mbs] != 0) {
54695b5ac22SHong Zhang     PetscReal af = fact->info.fill_ratio_needed;
54757622a8eSBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)fill,(double)af);CHKERRQ(ierr);
54857622a8eSBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);CHKERRQ(ierr);
54957622a8eSBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%g) for best performance.\n",(double)af);CHKERRQ(ierr);
55095b5ac22SHong Zhang   } else {
55195b5ac22SHong Zhang     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
55295b5ac22SHong Zhang   }
55395b5ac22SHong Zhang #endif
554d595f711SHong Zhang   ierr = MatSeqSBAIJSetNumericFactorization_inplace(fact,perm_identity);CHKERRQ(ierr);
555064503c5SHong Zhang   PetscFunctionReturn(0);
556064503c5SHong Zhang }
557d595f711SHong Zhang 
5584a2ae208SSatish Balay #undef __FUNCT__
5594a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N"
5600481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat C,Mat A,const MatFactorInfo *info)
56149b5e25fSSatish Balay {
5624c16a6a6SHong Zhang   Mat_SeqSBAIJ   *a   = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data;
5634c16a6a6SHong Zhang   IS             perm = b->row;
5646849ba73SBarry Smith   PetscErrorCode ierr;
5655d0c19d7SBarry Smith   const PetscInt *ai,*aj,*perm_ptr,mbs=a->mbs,*bi=b->i,*bj=b->j;
5665d0c19d7SBarry Smith   PetscInt       i,j;
5675d0c19d7SBarry Smith   PetscInt       *a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
5683bc0b13bSBarry Smith   PetscInt       bs  =A->rmap->bs,bs2 = a->bs2;
5694c16a6a6SHong Zhang   MatScalar      *ba = b->a,*aa,*ap,*dk,*uik;
5704c16a6a6SHong Zhang   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
57128de702eSHong Zhang   MatScalar      *work;
57213f74950SBarry Smith   PetscInt       *pivots;
5735f8bbccaSHong Zhang   PetscBool      allowzeropivot,zeropivotdetected;
5744c16a6a6SHong Zhang 
5754c16a6a6SHong Zhang   PetscFunctionBegin;
5764c16a6a6SHong Zhang   /* initialization */
5771795a4d1SJed Brown   ierr = PetscCalloc1(bs2*mbs,&rtmp);CHKERRQ(ierr);
578dcca6d9dSJed Brown   ierr = PetscMalloc2(mbs,&il,mbs,&jl);CHKERRQ(ierr);
5795f8bbccaSHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
5805f8bbccaSHong Zhang 
5816df5ee2eSHong Zhang   il[0] = 0;
5826df5ee2eSHong Zhang   for (i=0; i<mbs; i++) jl[i] = mbs;
5836df5ee2eSHong Zhang 
584dcca6d9dSJed Brown   ierr = PetscMalloc3(bs2,&dk,bs2,&uik,bs,&work);CHKERRQ(ierr);
585785e854fSJed Brown   ierr = PetscMalloc1(bs,&pivots);CHKERRQ(ierr);
5864c16a6a6SHong Zhang 
5874c16a6a6SHong Zhang   ierr = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr);
5884c16a6a6SHong Zhang 
5894c16a6a6SHong Zhang   /* check permutation */
5904c16a6a6SHong Zhang   if (!a->permute) {
5914c16a6a6SHong Zhang     ai = a->i; aj = a->j; aa = a->a;
5924c16a6a6SHong Zhang   } else {
5934c16a6a6SHong Zhang     ai   = a->inew; aj = a->jnew;
594785e854fSJed Brown     ierr = PetscMalloc1(bs2*ai[mbs],&aa);CHKERRQ(ierr);
5954c16a6a6SHong Zhang     ierr = PetscMemcpy(aa,a->a,bs2*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr);
596785e854fSJed Brown     ierr = PetscMalloc1(ai[mbs],&a2anew);CHKERRQ(ierr);
59713f74950SBarry Smith     ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(PetscInt));CHKERRQ(ierr);
5984c16a6a6SHong Zhang 
5994c16a6a6SHong Zhang     for (i=0; i<mbs; i++) {
6004c16a6a6SHong Zhang       jmin = ai[i]; jmax = ai[i+1];
6014c16a6a6SHong Zhang       for (j=jmin; j<jmax; j++) {
6024c16a6a6SHong Zhang         while (a2anew[j] != j) {
6034c16a6a6SHong Zhang           k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
6044c16a6a6SHong Zhang           for (k1=0; k1<bs2; k1++) {
6054c16a6a6SHong Zhang             dk[k1]       = aa[k*bs2+k1];
6064c16a6a6SHong Zhang             aa[k*bs2+k1] = aa[j*bs2+k1];
6074c16a6a6SHong Zhang             aa[j*bs2+k1] = dk[k1];
6084c16a6a6SHong Zhang           }
6094c16a6a6SHong Zhang         }
6104c16a6a6SHong Zhang         /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
6114c16a6a6SHong Zhang         if (i > aj[j]) {
6124c16a6a6SHong Zhang           /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */
6134c16a6a6SHong Zhang           ap = aa + j*bs2;                     /* ptr to the beginning of j-th block of aa */
6144c16a6a6SHong Zhang           for (k=0; k<bs2; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */
6154c16a6a6SHong Zhang           for (k=0; k<bs; k++) {               /* j-th block of aa <- dk^T */
6164c16a6a6SHong Zhang             for (k1=0; k1<bs; k1++) *ap++ = dk[k + bs*k1];
6174c16a6a6SHong Zhang           }
6184c16a6a6SHong Zhang         }
6194c16a6a6SHong Zhang       }
6204c16a6a6SHong Zhang     }
621323b833fSBarry Smith     ierr = PetscFree(a2anew);CHKERRQ(ierr);
6224c16a6a6SHong Zhang   }
6234c16a6a6SHong Zhang 
6244c16a6a6SHong Zhang   /* for each row k */
6254c16a6a6SHong Zhang   for (k = 0; k<mbs; k++) {
6264c16a6a6SHong Zhang 
6274c16a6a6SHong Zhang     /*initialize k-th row with elements nonzero in row perm(k) of A */
6284c16a6a6SHong Zhang     jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];
629057f5ba7SHong Zhang 
6304c16a6a6SHong Zhang     ap = aa + jmin*bs2;
6314c16a6a6SHong Zhang     for (j = jmin; j < jmax; j++) {
6324c16a6a6SHong Zhang       vj       = perm_ptr[aj[j]];   /* block col. index */
6334c16a6a6SHong Zhang       rtmp_ptr = rtmp + vj*bs2;
6344c16a6a6SHong Zhang       for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++;
6354c16a6a6SHong Zhang     }
6364c16a6a6SHong Zhang 
6374c16a6a6SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
6384c16a6a6SHong Zhang     ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
6394c16a6a6SHong Zhang     i    = jl[k]; /* first row to be added to k_th row  */
6404c16a6a6SHong Zhang 
641057f5ba7SHong Zhang     while (i < k) {
6424c16a6a6SHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
6434c16a6a6SHong Zhang 
6444c16a6a6SHong Zhang       /* compute multiplier */
6454c16a6a6SHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
6464c16a6a6SHong Zhang 
6474c16a6a6SHong Zhang       /* uik = -inv(Di)*U_bar(i,k) */
6484c16a6a6SHong Zhang       diag = ba + i*bs2;
6494c16a6a6SHong Zhang       u    = ba + ili*bs2;
6504c16a6a6SHong Zhang       ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
65196b95a6bSBarry Smith       PetscKernel_A_gets_A_minus_B_times_C(bs,uik,diag,u);
6524c16a6a6SHong Zhang 
6534c16a6a6SHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
65496b95a6bSBarry Smith       PetscKernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u);
6553bc0b13bSBarry Smith       ierr = PetscLogFlops(4.0*bs*bs2);CHKERRQ(ierr);
6564c16a6a6SHong Zhang 
6574c16a6a6SHong Zhang       /* update -U(i,k) */
6584c16a6a6SHong Zhang       ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
6594c16a6a6SHong Zhang 
6604c16a6a6SHong Zhang       /* add multiple of row i to k-th row ... */
6614c16a6a6SHong Zhang       jmin = ili + 1; jmax = bi[i+1];
6624c16a6a6SHong Zhang       if (jmin < jmax) {
6634c16a6a6SHong Zhang         for (j=jmin; j<jmax; j++) {
6644c16a6a6SHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j) */
6654c16a6a6SHong Zhang           rtmp_ptr = rtmp + bj[j]*bs2;
6664c16a6a6SHong Zhang           u        = ba + j*bs2;
66796b95a6bSBarry Smith           PetscKernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u);
6684c16a6a6SHong Zhang         }
6693bc0b13bSBarry Smith         ierr = PetscLogFlops(2.0*bs*bs2*(jmax-jmin));CHKERRQ(ierr);
6704c16a6a6SHong Zhang 
6714c16a6a6SHong Zhang         /* ... add i to row list for next nonzero entry */
6724c16a6a6SHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
6734c16a6a6SHong Zhang         j     = bj[jmin];
6744c16a6a6SHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
6754c16a6a6SHong Zhang       }
6764c16a6a6SHong Zhang       i = nexti;
6774c16a6a6SHong Zhang     }
6784c16a6a6SHong Zhang 
6794c16a6a6SHong Zhang     /* save nonzero entries in k-th row of U ... */
6804c16a6a6SHong Zhang 
6814c16a6a6SHong Zhang     /* invert diagonal block */
6824c16a6a6SHong Zhang     diag = ba+k*bs2;
6834c16a6a6SHong Zhang     ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr);
6845f8bbccaSHong Zhang 
6855f8bbccaSHong Zhang     ierr = PetscKernel_A_gets_inverse_A(bs,diag,pivots,work,allowzeropivot,&zeropivotdetected);CHKERRQ(ierr);
686*7b6c816cSBarry Smith     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
6874c16a6a6SHong Zhang 
6884c16a6a6SHong Zhang     jmin = bi[k]; jmax = bi[k+1];
6894c16a6a6SHong Zhang     if (jmin < jmax) {
6904c16a6a6SHong Zhang       for (j=jmin; j<jmax; j++) {
6914c16a6a6SHong Zhang         vj       = bj[j];      /* block col. index of U */
6924c16a6a6SHong Zhang         u        = ba + j*bs2;
6934c16a6a6SHong Zhang         rtmp_ptr = rtmp + vj*bs2;
6944c16a6a6SHong Zhang         for (k1=0; k1<bs2; k1++) {
6954c16a6a6SHong Zhang           *u++        = *rtmp_ptr;
6964c16a6a6SHong Zhang           *rtmp_ptr++ = 0.0;
6974c16a6a6SHong Zhang         }
6984c16a6a6SHong Zhang       }
6994c16a6a6SHong Zhang 
7004c16a6a6SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
7014c16a6a6SHong Zhang       il[k] = jmin;
7024c16a6a6SHong Zhang       i     = bj[jmin];
7034c16a6a6SHong Zhang       jl[k] = jl[i]; jl[i] = k;
7044c16a6a6SHong Zhang     }
7054c16a6a6SHong Zhang   }
7064c16a6a6SHong Zhang 
7074c16a6a6SHong Zhang   ierr = PetscFree(rtmp);CHKERRQ(ierr);
708d8c74875SBarry Smith   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
709d8c74875SBarry Smith   ierr = PetscFree3(dk,uik,work);CHKERRQ(ierr);
7104c16a6a6SHong Zhang   ierr = PetscFree(pivots);CHKERRQ(ierr);
7114c16a6a6SHong Zhang   if (a->permute) {
7124c16a6a6SHong Zhang     ierr = PetscFree(aa);CHKERRQ(ierr);
7134c16a6a6SHong Zhang   }
7144c16a6a6SHong Zhang 
7154c16a6a6SHong Zhang   ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr);
71626fbe8dcSKarl Rupp 
7174f79d315SHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_N_inplace;
7184f79d315SHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_N_inplace;
7194f79d315SHong Zhang   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_N_inplace;
7204f79d315SHong Zhang   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_N_inplace;
721db4efbfdSBarry Smith 
7224c16a6a6SHong Zhang   C->assembled    = PETSC_TRUE;
7234c16a6a6SHong Zhang   C->preallocated = PETSC_TRUE;
72426fbe8dcSKarl Rupp 
725efee365bSSatish Balay   ierr = PetscLogFlops(1.3333*bs*bs2*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
7264c16a6a6SHong Zhang   PetscFunctionReturn(0);
7274c16a6a6SHong Zhang }
728d4edadd4SHong Zhang 
7294a2ae208SSatish Balay #undef __FUNCT__
7304a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering"
7310481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info)
732671cb588SHong Zhang {
733671cb588SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data;
734dfbe8321SBarry Smith   PetscErrorCode ierr;
73513f74950SBarry Smith   PetscInt       i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
73613f74950SBarry Smith   PetscInt       *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
7373bc0b13bSBarry Smith   PetscInt       bs  =A->rmap->bs,bs2 = a->bs2;
738671cb588SHong Zhang   MatScalar      *ba = b->a,*aa,*ap,*dk,*uik;
739671cb588SHong Zhang   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
74028de702eSHong Zhang   MatScalar      *work;
74113f74950SBarry Smith   PetscInt       *pivots;
7425f8bbccaSHong Zhang   PetscBool      allowzeropivot,zeropivotdetected;
743671cb588SHong Zhang 
744671cb588SHong Zhang   PetscFunctionBegin;
7451795a4d1SJed Brown   ierr = PetscCalloc1(bs2*mbs,&rtmp);CHKERRQ(ierr);
746dcca6d9dSJed Brown   ierr = PetscMalloc2(mbs,&il,mbs,&jl);CHKERRQ(ierr);
7476df5ee2eSHong Zhang   il[0] = 0;
7486df5ee2eSHong Zhang   for (i=0; i<mbs; i++) jl[i] = mbs;
7496df5ee2eSHong Zhang 
750dcca6d9dSJed Brown   ierr = PetscMalloc3(bs2,&dk,bs2,&uik,bs,&work);CHKERRQ(ierr);
751785e854fSJed Brown   ierr = PetscMalloc1(bs,&pivots);CHKERRQ(ierr);
7525f8bbccaSHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
753671cb588SHong Zhang 
754671cb588SHong Zhang   ai = a->i; aj = a->j; aa = a->a;
755671cb588SHong Zhang 
756671cb588SHong Zhang   /* for each row k */
757671cb588SHong Zhang   for (k = 0; k<mbs; k++) {
758671cb588SHong Zhang 
759671cb588SHong Zhang     /*initialize k-th row with elements nonzero in row k of A */
760671cb588SHong Zhang     jmin = ai[k]; jmax = ai[k+1];
761671cb588SHong Zhang     ap   = aa + jmin*bs2;
762671cb588SHong Zhang     for (j = jmin; j < jmax; j++) {
763671cb588SHong Zhang       vj       = aj[j];   /* block col. index */
764671cb588SHong Zhang       rtmp_ptr = rtmp + vj*bs2;
765671cb588SHong Zhang       for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++;
766671cb588SHong Zhang     }
767671cb588SHong Zhang 
768671cb588SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
769671cb588SHong Zhang     ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
770671cb588SHong Zhang     i    = jl[k]; /* first row to be added to k_th row  */
771671cb588SHong Zhang 
772057f5ba7SHong Zhang     while (i < k) {
773671cb588SHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
774671cb588SHong Zhang 
775671cb588SHong Zhang       /* compute multiplier */
776671cb588SHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
777671cb588SHong Zhang 
778671cb588SHong Zhang       /* uik = -inv(Di)*U_bar(i,k) */
779671cb588SHong Zhang       diag = ba + i*bs2;
780671cb588SHong Zhang       u    = ba + ili*bs2;
781671cb588SHong Zhang       ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
78296b95a6bSBarry Smith       PetscKernel_A_gets_A_minus_B_times_C(bs,uik,diag,u);
783671cb588SHong Zhang 
784671cb588SHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
78596b95a6bSBarry Smith       PetscKernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u);
7863bc0b13bSBarry Smith       ierr = PetscLogFlops(2.0*bs*bs2);CHKERRQ(ierr);
787671cb588SHong Zhang 
788671cb588SHong Zhang       /* update -U(i,k) */
789671cb588SHong Zhang       ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
790671cb588SHong Zhang 
791671cb588SHong Zhang       /* add multiple of row i to k-th row ... */
792671cb588SHong Zhang       jmin = ili + 1; jmax = bi[i+1];
793671cb588SHong Zhang       if (jmin < jmax) {
794671cb588SHong Zhang         for (j=jmin; j<jmax; j++) {
795671cb588SHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j) */
796671cb588SHong Zhang           rtmp_ptr = rtmp + bj[j]*bs2;
797671cb588SHong Zhang           u        = ba + j*bs2;
79896b95a6bSBarry Smith           PetscKernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u);
799671cb588SHong Zhang         }
8003bc0b13bSBarry Smith         ierr = PetscLogFlops(2.0*bs*bs2*(jmax-jmin));CHKERRQ(ierr);
801671cb588SHong Zhang 
802671cb588SHong Zhang         /* ... add i to row list for next nonzero entry */
803671cb588SHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
804671cb588SHong Zhang         j     = bj[jmin];
805671cb588SHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
806671cb588SHong Zhang       }
807671cb588SHong Zhang       i = nexti;
808671cb588SHong Zhang     }
809671cb588SHong Zhang 
810671cb588SHong Zhang     /* save nonzero entries in k-th row of U ... */
811671cb588SHong Zhang 
812671cb588SHong Zhang     /* invert diagonal block */
813671cb588SHong Zhang     diag = ba+k*bs2;
814671cb588SHong Zhang     ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr);
8155f8bbccaSHong Zhang 
8165f8bbccaSHong Zhang     ierr = PetscKernel_A_gets_inverse_A(bs,diag,pivots,work,allowzeropivot,&zeropivotdetected);CHKERRQ(ierr);
817*7b6c816cSBarry Smith     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
818671cb588SHong Zhang 
819671cb588SHong Zhang     jmin = bi[k]; jmax = bi[k+1];
820671cb588SHong Zhang     if (jmin < jmax) {
821671cb588SHong Zhang       for (j=jmin; j<jmax; j++) {
822671cb588SHong Zhang         vj       = bj[j];      /* block col. index of U */
823671cb588SHong Zhang         u        = ba + j*bs2;
824671cb588SHong Zhang         rtmp_ptr = rtmp + vj*bs2;
825671cb588SHong Zhang         for (k1=0; k1<bs2; k1++) {
826671cb588SHong Zhang           *u++        = *rtmp_ptr;
827671cb588SHong Zhang           *rtmp_ptr++ = 0.0;
828671cb588SHong Zhang         }
829671cb588SHong Zhang       }
830671cb588SHong Zhang 
831671cb588SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
832671cb588SHong Zhang       il[k] = jmin;
833671cb588SHong Zhang       i     = bj[jmin];
834671cb588SHong Zhang       jl[k] = jl[i]; jl[i] = k;
835671cb588SHong Zhang     }
836671cb588SHong Zhang   }
837671cb588SHong Zhang 
838671cb588SHong Zhang   ierr = PetscFree(rtmp);CHKERRQ(ierr);
839d8c74875SBarry Smith   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
840d8c74875SBarry Smith   ierr = PetscFree3(dk,uik,work);CHKERRQ(ierr);
841671cb588SHong Zhang   ierr = PetscFree(pivots);CHKERRQ(ierr);
842671cb588SHong Zhang 
8434f79d315SHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_N_NaturalOrdering_inplace;
8444f79d315SHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_N_NaturalOrdering_inplace;
8454f79d315SHong Zhang   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_N_NaturalOrdering_inplace;
8464f79d315SHong Zhang   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_N_NaturalOrdering_inplace;
847671cb588SHong Zhang   C->assembled           = PETSC_TRUE;
848671cb588SHong Zhang   C->preallocated        = PETSC_TRUE;
84926fbe8dcSKarl Rupp 
850efee365bSSatish Balay   ierr = PetscLogFlops(1.3333*bs*bs2*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
851671cb588SHong Zhang   PetscFunctionReturn(0);
852671cb588SHong Zhang }
853671cb588SHong Zhang 
85449b5e25fSSatish Balay /*
855fcf159c0SHong Zhang     Numeric U^T*D*U factorization for SBAIJ format. Modified from SNF of YSMP.
856cc0c071aSHong Zhang     Version for blocks 2 by 2.
85749b5e25fSSatish Balay */
8584a2ae208SSatish Balay #undef __FUNCT__
8594a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2"
8600481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat C,Mat A,const MatFactorInfo *info)
86149b5e25fSSatish Balay {
86291602c66SHong Zhang   Mat_SeqSBAIJ   *a   = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data;
863cc0c071aSHong Zhang   IS             perm = b->row;
8646849ba73SBarry Smith   PetscErrorCode ierr;
8655d0c19d7SBarry Smith   const PetscInt *ai,*aj,*perm_ptr;
8665d0c19d7SBarry Smith   PetscInt       i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
8675d0c19d7SBarry Smith   PetscInt       *a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
868d8c74875SBarry Smith   MatScalar      *ba = b->a,*aa,*ap;
869d8c74875SBarry Smith   MatScalar      *u,*diag,*rtmp,*rtmp_ptr,dk[4],uik[4];
87014d2772eSHong Zhang   PetscReal      shift = info->shiftamount;
871a455e926SHong Zhang   PetscBool      allowzeropivot,zeropivotdetected;
87249b5e25fSSatish Balay 
87349b5e25fSSatish Balay   PetscFunctionBegin;
8740164db54SHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
8750164db54SHong Zhang 
87691602c66SHong Zhang   /* initialization */
87791602c66SHong Zhang   /* il and jl record the first nonzero element in each row of the accessing
87891602c66SHong Zhang      window U(0:k, k:mbs-1).
87991602c66SHong Zhang      jl:    list of rows to be added to uneliminated rows
88091602c66SHong Zhang             i>= k: jl(i) is the first row to be added to row i
88191602c66SHong Zhang             i<  k: jl(i) is the row following row i in some list of rows
88291602c66SHong Zhang             jl(i) = mbs indicates the end of a list
88391602c66SHong Zhang      il(i): points to the first nonzero element in columns k,...,mbs-1 of
88491602c66SHong Zhang             row i of U */
8851795a4d1SJed Brown   ierr = PetscCalloc1(4*mbs,&rtmp);CHKERRQ(ierr);
886dcca6d9dSJed Brown   ierr = PetscMalloc2(mbs,&il,mbs,&jl);CHKERRQ(ierr);
8876df5ee2eSHong Zhang   il[0] = 0;
8886df5ee2eSHong Zhang   for (i=0; i<mbs; i++) jl[i] = mbs;
8896df5ee2eSHong Zhang 
890cc0c071aSHong Zhang   ierr = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr);
891cc0c071aSHong Zhang 
892cc0c071aSHong Zhang   /* check permutation */
893cc0c071aSHong Zhang   if (!a->permute) {
894cc0c071aSHong Zhang     ai = a->i; aj = a->j; aa = a->a;
895cc0c071aSHong Zhang   } else {
896cc0c071aSHong Zhang     ai   = a->inew; aj = a->jnew;
897785e854fSJed Brown     ierr = PetscMalloc1(4*ai[mbs],&aa);CHKERRQ(ierr);
898cc0c071aSHong Zhang     ierr = PetscMemcpy(aa,a->a,4*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr);
899785e854fSJed Brown     ierr = PetscMalloc1(ai[mbs],&a2anew);CHKERRQ(ierr);
90013f74950SBarry Smith     ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(PetscInt));CHKERRQ(ierr);
901cc0c071aSHong Zhang 
902cc0c071aSHong Zhang     for (i=0; i<mbs; i++) {
903cc0c071aSHong Zhang       jmin = ai[i]; jmax = ai[i+1];
904cc0c071aSHong Zhang       for (j=jmin; j<jmax; j++) {
905cc0c071aSHong Zhang         while (a2anew[j] != j) {
906cc0c071aSHong Zhang           k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
907cc0c071aSHong Zhang           for (k1=0; k1<4; k1++) {
908cc0c071aSHong Zhang             dk[k1]     = aa[k*4+k1];
909cc0c071aSHong Zhang             aa[k*4+k1] = aa[j*4+k1];
910cc0c071aSHong Zhang             aa[j*4+k1] = dk[k1];
911cc0c071aSHong Zhang           }
912cc0c071aSHong Zhang         }
913cc0c071aSHong Zhang         /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
914cc0c071aSHong Zhang         if (i > aj[j]) {
915a1723e09SHong Zhang           /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */
916cc0c071aSHong Zhang           ap    = aa + j*4;  /* ptr to the beginning of the block */
917cc0c071aSHong Zhang           dk[1] = ap[1];     /* swap ap[1] and ap[2] */
918cc0c071aSHong Zhang           ap[1] = ap[2];
919cc0c071aSHong Zhang           ap[2] = dk[1];
920cc0c071aSHong Zhang         }
921cc0c071aSHong Zhang       }
922cc0c071aSHong Zhang     }
923ac355199SBarry Smith     ierr = PetscFree(a2anew);CHKERRQ(ierr);
924cc0c071aSHong Zhang   }
9253845f261SHong Zhang 
92691602c66SHong Zhang   /* for each row k */
92791602c66SHong Zhang   for (k = 0; k<mbs; k++) {
92891602c66SHong Zhang 
92991602c66SHong Zhang     /*initialize k-th row with elements nonzero in row perm(k) of A */
930cc0c071aSHong Zhang     jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];
931cc0c071aSHong Zhang     ap   = aa + jmin*4;
93291602c66SHong Zhang     for (j = jmin; j < jmax; j++) {
933cc0c071aSHong Zhang       vj       = perm_ptr[aj[j]];   /* block col. index */
934cc0c071aSHong Zhang       rtmp_ptr = rtmp + vj*4;
935cc0c071aSHong Zhang       for (i=0; i<4; i++) *rtmp_ptr++ = *ap++;
93691602c66SHong Zhang     }
93791602c66SHong Zhang 
93891602c66SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
939cc0c071aSHong Zhang     ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr);
94091602c66SHong Zhang     i    = jl[k]; /* first row to be added to k_th row  */
94191602c66SHong Zhang 
942057f5ba7SHong Zhang     while (i < k) {
94391602c66SHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
94491602c66SHong Zhang 
9453845f261SHong Zhang       /* compute multiplier */
94691602c66SHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
9473845f261SHong Zhang 
9483845f261SHong Zhang       /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */
949cc0c071aSHong Zhang       diag   = ba + i*4;
950cc0c071aSHong Zhang       u      = ba + ili*4;
951cc0c071aSHong Zhang       uik[0] = -(diag[0]*u[0] + diag[2]*u[1]);
952cc0c071aSHong Zhang       uik[1] = -(diag[1]*u[0] + diag[3]*u[1]);
953cc0c071aSHong Zhang       uik[2] = -(diag[0]*u[2] + diag[2]*u[3]);
954cc0c071aSHong Zhang       uik[3] = -(diag[1]*u[2] + diag[3]*u[3]);
9553845f261SHong Zhang 
9563845f261SHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */
957cc0c071aSHong Zhang       dk[0] += uik[0]*u[0] + uik[1]*u[1];
958cc0c071aSHong Zhang       dk[1] += uik[2]*u[0] + uik[3]*u[1];
959cc0c071aSHong Zhang       dk[2] += uik[0]*u[2] + uik[1]*u[3];
960cc0c071aSHong Zhang       dk[3] += uik[2]*u[2] + uik[3]*u[3];
9613845f261SHong Zhang 
962dc0b31edSSatish Balay       ierr = PetscLogFlops(16.0*2.0);CHKERRQ(ierr);
963187a9f4bSHong Zhang 
9643845f261SHong Zhang       /* update -U(i,k): ba[ili] = uik */
965cc0c071aSHong Zhang       ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr);
96691602c66SHong Zhang 
96791602c66SHong Zhang       /* add multiple of row i to k-th row ... */
96891602c66SHong Zhang       jmin = ili + 1; jmax = bi[i+1];
96991602c66SHong Zhang       if (jmin < jmax) {
9703845f261SHong Zhang         for (j=jmin; j<jmax; j++) {
9713845f261SHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */
972cc0c071aSHong Zhang           rtmp_ptr     = rtmp + bj[j]*4;
973cc0c071aSHong Zhang           u            = ba + j*4;
974cc0c071aSHong Zhang           rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1];
975cc0c071aSHong Zhang           rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1];
976cc0c071aSHong Zhang           rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3];
977cc0c071aSHong Zhang           rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3];
9783845f261SHong Zhang         }
979dc0b31edSSatish Balay         ierr = PetscLogFlops(16.0*(jmax-jmin));CHKERRQ(ierr);
9803845f261SHong Zhang 
98191602c66SHong Zhang         /* ... add i to row list for next nonzero entry */
98291602c66SHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
98391602c66SHong Zhang         j     = bj[jmin];
98491602c66SHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
98591602c66SHong Zhang       }
986a1723e09SHong Zhang       i = nexti;
98791602c66SHong Zhang     }
988cc0c071aSHong Zhang 
98991602c66SHong Zhang     /* save nonzero entries in k-th row of U ... */
9903845f261SHong Zhang 
991cc0c071aSHong Zhang     /* invert diagonal block */
992cc0c071aSHong Zhang     diag = ba+k*4;
993cc0c071aSHong Zhang     ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr);
994a455e926SHong Zhang     ierr = PetscKernel_A_gets_inverse_A_2(diag,shift,allowzeropivot,&zeropivotdetected);CHKERRQ(ierr);
995*7b6c816cSBarry Smith     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
9963845f261SHong Zhang 
99791602c66SHong Zhang     jmin = bi[k]; jmax = bi[k+1];
99891602c66SHong Zhang     if (jmin < jmax) {
99991602c66SHong Zhang       for (j=jmin; j<jmax; j++) {
1000cc0c071aSHong Zhang         vj       = bj[j];      /* block col. index of U */
1001cc0c071aSHong Zhang         u        = ba + j*4;
1002cc0c071aSHong Zhang         rtmp_ptr = rtmp + vj*4;
1003cc0c071aSHong Zhang         for (k1=0; k1<4; k1++) {
1004cc0c071aSHong Zhang           *u++        = *rtmp_ptr;
1005cc0c071aSHong Zhang           *rtmp_ptr++ = 0.0;
1006cc0c071aSHong Zhang         }
1007cc0c071aSHong Zhang       }
10083845f261SHong Zhang 
100991602c66SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
101091602c66SHong Zhang       il[k] = jmin;
101191602c66SHong Zhang       i     = bj[jmin];
101291602c66SHong Zhang       jl[k] = jl[i]; jl[i] = k;
101391602c66SHong Zhang     }
101491602c66SHong Zhang   }
10153845f261SHong Zhang 
101649b5e25fSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
1017d8c74875SBarry Smith   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
101891602c66SHong Zhang   if (a->permute) {
101991602c66SHong Zhang     ierr = PetscFree(aa);CHKERRQ(ierr);
102091602c66SHong Zhang   }
1021cc0c071aSHong Zhang   ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr);
102226fbe8dcSKarl Rupp 
10234f79d315SHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_2_inplace;
10244f79d315SHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_2_inplace;
102549b5e25fSSatish Balay   C->assembled           = PETSC_TRUE;
10265c0bcdfcSHong Zhang   C->preallocated        = PETSC_TRUE;
102726fbe8dcSKarl Rupp 
1028efee365bSSatish Balay   ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
102949b5e25fSSatish Balay   PetscFunctionReturn(0);
103049b5e25fSSatish Balay }
103191602c66SHong Zhang 
103249b5e25fSSatish Balay /*
103349b5e25fSSatish Balay       Version for when blocks are 2 by 2 Using natural ordering
103449b5e25fSSatish Balay */
10354a2ae208SSatish Balay #undef __FUNCT__
10364a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering"
10370481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info)
103849b5e25fSSatish Balay {
1039ab27746eSHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data;
1040dfbe8321SBarry Smith   PetscErrorCode ierr;
104113f74950SBarry Smith   PetscInt       i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
104213f74950SBarry Smith   PetscInt       *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
1043d8c74875SBarry Smith   MatScalar      *ba = b->a,*aa,*ap,dk[8],uik[8];
1044ab27746eSHong Zhang   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
104514d2772eSHong Zhang   PetscReal      shift = info->shiftamount;
1046a455e926SHong Zhang   PetscBool      allowzeropivot,zeropivotdetected;
104749b5e25fSSatish Balay 
104849b5e25fSSatish Balay   PetscFunctionBegin;
10490164db54SHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
10500164db54SHong Zhang 
1051ab27746eSHong Zhang   /* initialization */
1052ab27746eSHong Zhang   /* il and jl record the first nonzero element in each row of the accessing
1053ab27746eSHong Zhang      window U(0:k, k:mbs-1).
1054ab27746eSHong Zhang      jl:    list of rows to be added to uneliminated rows
1055ab27746eSHong Zhang             i>= k: jl(i) is the first row to be added to row i
1056ab27746eSHong Zhang             i<  k: jl(i) is the row following row i in some list of rows
1057ab27746eSHong Zhang             jl(i) = mbs indicates the end of a list
1058ab27746eSHong Zhang      il(i): points to the first nonzero element in columns k,...,mbs-1 of
1059ab27746eSHong Zhang             row i of U */
10601795a4d1SJed Brown   ierr = PetscCalloc1(4*mbs,&rtmp);CHKERRQ(ierr);
1061dcca6d9dSJed Brown   ierr = PetscMalloc2(mbs,&il,mbs,&jl);CHKERRQ(ierr);
10626df5ee2eSHong Zhang   il[0] = 0;
10636df5ee2eSHong Zhang   for (i=0; i<mbs; i++) jl[i] = mbs;
10646df5ee2eSHong Zhang 
1065ab27746eSHong Zhang   ai = a->i; aj = a->j; aa = a->a;
1066ab27746eSHong Zhang 
1067ab27746eSHong Zhang   /* for each row k */
1068ab27746eSHong Zhang   for (k = 0; k<mbs; k++) {
1069ab27746eSHong Zhang 
1070ab27746eSHong Zhang     /*initialize k-th row with elements nonzero in row k of A */
1071ab27746eSHong Zhang     jmin = ai[k]; jmax = ai[k+1];
1072ab27746eSHong Zhang     ap   = aa + jmin*4;
1073ab27746eSHong Zhang     for (j = jmin; j < jmax; j++) {
1074ab27746eSHong Zhang       vj       = aj[j];   /* block col. index */
1075ab27746eSHong Zhang       rtmp_ptr = rtmp + vj*4;
1076ab27746eSHong Zhang       for (i=0; i<4; i++) *rtmp_ptr++ = *ap++;
107749b5e25fSSatish Balay     }
1078ab27746eSHong Zhang 
1079ab27746eSHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
1080ab27746eSHong Zhang     ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr);
1081ab27746eSHong Zhang     i    = jl[k]; /* first row to be added to k_th row  */
1082ab27746eSHong Zhang 
1083057f5ba7SHong Zhang     while (i < k) {
1084ab27746eSHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
1085ab27746eSHong Zhang 
1086ab27746eSHong Zhang       /* compute multiplier */
1087ab27746eSHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
1088ab27746eSHong Zhang 
1089ab27746eSHong Zhang       /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */
1090ab27746eSHong Zhang       diag   = ba + i*4;
1091ab27746eSHong Zhang       u      = ba + ili*4;
1092ab27746eSHong Zhang       uik[0] = -(diag[0]*u[0] + diag[2]*u[1]);
1093ab27746eSHong Zhang       uik[1] = -(diag[1]*u[0] + diag[3]*u[1]);
1094ab27746eSHong Zhang       uik[2] = -(diag[0]*u[2] + diag[2]*u[3]);
1095ab27746eSHong Zhang       uik[3] = -(diag[1]*u[2] + diag[3]*u[3]);
1096ab27746eSHong Zhang 
1097ab27746eSHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */
1098ab27746eSHong Zhang       dk[0] += uik[0]*u[0] + uik[1]*u[1];
1099ab27746eSHong Zhang       dk[1] += uik[2]*u[0] + uik[3]*u[1];
1100ab27746eSHong Zhang       dk[2] += uik[0]*u[2] + uik[1]*u[3];
1101ab27746eSHong Zhang       dk[3] += uik[2]*u[2] + uik[3]*u[3];
1102ab27746eSHong Zhang 
1103dc0b31edSSatish Balay       ierr = PetscLogFlops(16.0*2.0);CHKERRQ(ierr);
1104187a9f4bSHong Zhang 
1105ab27746eSHong Zhang       /* update -U(i,k): ba[ili] = uik */
1106ab27746eSHong Zhang       ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr);
1107ab27746eSHong Zhang 
1108ab27746eSHong Zhang       /* add multiple of row i to k-th row ... */
1109ab27746eSHong Zhang       jmin = ili + 1; jmax = bi[i+1];
1110ab27746eSHong Zhang       if (jmin < jmax) {
1111ab27746eSHong Zhang         for (j=jmin; j<jmax; j++) {
1112ab27746eSHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */
1113ab27746eSHong Zhang           rtmp_ptr     = rtmp + bj[j]*4;
1114ab27746eSHong Zhang           u            = ba + j*4;
1115ab27746eSHong Zhang           rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1];
1116ab27746eSHong Zhang           rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1];
1117ab27746eSHong Zhang           rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3];
1118ab27746eSHong Zhang           rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3];
111949b5e25fSSatish Balay         }
1120dc0b31edSSatish Balay         ierr = PetscLogFlops(16.0*(jmax-jmin));CHKERRQ(ierr);
1121ab27746eSHong Zhang 
1122ab27746eSHong Zhang         /* ... add i to row list for next nonzero entry */
1123ab27746eSHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
1124ab27746eSHong Zhang         j     = bj[jmin];
1125ab27746eSHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
112649b5e25fSSatish Balay       }
1127ab27746eSHong Zhang       i = nexti;
112849b5e25fSSatish Balay     }
1129ab27746eSHong Zhang 
1130ab27746eSHong Zhang     /* save nonzero entries in k-th row of U ... */
1131ab27746eSHong Zhang 
113249b5e25fSSatish Balay     /* invert diagonal block */
1133ab27746eSHong Zhang     diag = ba+k*4;
1134ab27746eSHong Zhang     ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr);
1135a455e926SHong Zhang     ierr = PetscKernel_A_gets_inverse_A_2(diag,shift,allowzeropivot,&zeropivotdetected);CHKERRQ(ierr);
1136*7b6c816cSBarry Smith     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1137ab27746eSHong Zhang 
1138ab27746eSHong Zhang     jmin = bi[k]; jmax = bi[k+1];
1139ab27746eSHong Zhang     if (jmin < jmax) {
1140ab27746eSHong Zhang       for (j=jmin; j<jmax; j++) {
1141ab27746eSHong Zhang         vj       = bj[j];      /* block col. index of U */
1142ab27746eSHong Zhang         u        = ba + j*4;
1143ab27746eSHong Zhang         rtmp_ptr = rtmp + vj*4;
1144ab27746eSHong Zhang         for (k1=0; k1<4; k1++) {
1145ab27746eSHong Zhang           *u++        = *rtmp_ptr;
1146ab27746eSHong Zhang           *rtmp_ptr++ = 0.0;
1147ab27746eSHong Zhang         }
1148ab27746eSHong Zhang       }
1149ab27746eSHong Zhang 
1150ab27746eSHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
1151ab27746eSHong Zhang       il[k] = jmin;
1152ab27746eSHong Zhang       i     = bj[jmin];
1153ab27746eSHong Zhang       jl[k] = jl[i]; jl[i] = k;
1154ab27746eSHong Zhang     }
115549b5e25fSSatish Balay   }
115649b5e25fSSatish Balay 
115749b5e25fSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
1158d8c74875SBarry Smith   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
1159ab27746eSHong Zhang 
11604f79d315SHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_2_NaturalOrdering_inplace;
11614f79d315SHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_2_NaturalOrdering_inplace;
11624f79d315SHong Zhang   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_2_NaturalOrdering_inplace;
11634f79d315SHong Zhang   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_2_NaturalOrdering_inplace;
116449b5e25fSSatish Balay   C->assembled           = PETSC_TRUE;
11655c0bcdfcSHong Zhang   C->preallocated        = PETSC_TRUE;
116626fbe8dcSKarl Rupp 
1167efee365bSSatish Balay   ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
116849b5e25fSSatish Balay   PetscFunctionReturn(0);
116949b5e25fSSatish Balay }
117049b5e25fSSatish Balay 
117149b5e25fSSatish Balay /*
117298a8e78dSHong Zhang     Numeric U^T*D*U factorization for SBAIJ format.
117391602c66SHong Zhang     Version for blocks are 1 by 1.
117449b5e25fSSatish Balay */
11754a2ae208SSatish Balay #undef __FUNCT__
1176d595f711SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1_inplace"
1177d595f711SHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_inplace(Mat C,Mat A,const MatFactorInfo *info)
117849b5e25fSSatish Balay {
117949b5e25fSSatish Balay   Mat_SeqSBAIJ   *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ*)C->data;
118049b5e25fSSatish Balay   IS             ip=b->row;
11816849ba73SBarry Smith   PetscErrorCode ierr;
11825d0c19d7SBarry Smith   const PetscInt *ai,*aj,*rip;
11835d0c19d7SBarry Smith   PetscInt       *a2anew,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j,*bcol;
1184997a0949SHong Zhang   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
1185997a0949SHong Zhang   MatScalar      *rtmp,*ba=b->a,*bval,*aa,dk,uikdi;
11860e95ead3SHong Zhang   PetscReal      rs;
11870e95ead3SHong Zhang   FactorShiftCtx sctx;
118849b5e25fSSatish Balay 
118949b5e25fSSatish Balay   PetscFunctionBegin;
11900e95ead3SHong Zhang   /* MatPivotSetUp(): initialize shift context sctx */
11910e95ead3SHong Zhang   ierr = PetscMemzero(&sctx,sizeof(FactorShiftCtx));CHKERRQ(ierr);
11923cea4cbeSHong Zhang 
119349b5e25fSSatish Balay   ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr);
1194cb718733SHong Zhang   if (!a->permute) {
11952d007305SHong Zhang     ai = a->i; aj = a->j; aa = a->a;
11962d007305SHong Zhang   } else {
11972d007305SHong Zhang     ai     = a->inew; aj = a->jnew;
1198fff829cfSHong Zhang     nz     = ai[mbs];
1199785e854fSJed Brown     ierr   = PetscMalloc1(nz,&aa);CHKERRQ(ierr);
1200fff829cfSHong Zhang     a2anew = a->a2anew;
1201997a0949SHong Zhang     bval   = a->a;
1202fff829cfSHong Zhang     for (j=0; j<nz; j++) {
1203997a0949SHong Zhang       aa[a2anew[j]] = *(bval++);
12042d007305SHong Zhang     }
12052d007305SHong Zhang   }
12062d007305SHong Zhang 
120791602c66SHong Zhang   /* initialization */
120849b5e25fSSatish Balay   /* il and jl record the first nonzero element in each row of the accessing
120949b5e25fSSatish Balay      window U(0:k, k:mbs-1).
121049b5e25fSSatish Balay      jl:    list of rows to be added to uneliminated rows
121149b5e25fSSatish Balay             i>= k: jl(i) is the first row to be added to row i
121249b5e25fSSatish Balay             i<  k: jl(i) is the row following row i in some list of rows
121349b5e25fSSatish Balay             jl(i) = mbs indicates the end of a list
121449b5e25fSSatish Balay      il(i): points to the first nonzero element in columns k,...,mbs-1 of
121549b5e25fSSatish Balay             row i of U */
1216dcca6d9dSJed Brown   ierr = PetscMalloc3(mbs,&rtmp,mbs,&il,mbs,&jl);CHKERRQ(ierr);
1217997a0949SHong Zhang 
1218997a0949SHong Zhang   do {
121907b50cabSHong Zhang     sctx.newshift = PETSC_FALSE;
12206df5ee2eSHong Zhang     il[0] = 0;
122149b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
12226df5ee2eSHong Zhang       rtmp[i] = 0.0; jl[i] = mbs;
122349b5e25fSSatish Balay     }
1224997a0949SHong Zhang 
122549b5e25fSSatish Balay     for (k = 0; k<mbs; k++) {
1226997a0949SHong Zhang       /*initialize k-th row by the perm[k]-th row of A */
122749b5e25fSSatish Balay       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
12287701de02SHong Zhang       bval = ba + bi[k];
122949b5e25fSSatish Balay       for (j = jmin; j < jmax; j++) {
1230997a0949SHong Zhang         col       = rip[aj[j]];
1231997a0949SHong Zhang         rtmp[col] = aa[j];
12327701de02SHong Zhang         *bval++   = 0.0; /* for in-place factorization */
123349b5e25fSSatish Balay       }
12343cea4cbeSHong Zhang 
12353cea4cbeSHong Zhang       /* shift the diagonal of the matrix */
12363cea4cbeSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
123749b5e25fSSatish Balay 
123891602c66SHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
123949b5e25fSSatish Balay       dk = rtmp[k];
124049b5e25fSSatish Balay       i  = jl[k]; /* first row to be added to k_th row  */
124149b5e25fSSatish Balay 
1242057f5ba7SHong Zhang       while (i < k) {
124349b5e25fSSatish Balay         nexti = jl[i]; /* next row to be added to k_th row */
124449b5e25fSSatish Balay 
1245fff829cfSHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
124649b5e25fSSatish Balay         ili     = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
1247997a0949SHong Zhang         uikdi   = -ba[ili]*ba[bi[i]]; /* diagonal(k) */
124849b5e25fSSatish Balay         dk     += uikdi*ba[ili];
1249658e7b3eSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
125049b5e25fSSatish Balay 
1251997a0949SHong Zhang         /* add multiple of row i to k-th row */
125249b5e25fSSatish Balay         jmin = ili + 1; jmax = bi[i+1];
125349b5e25fSSatish Balay         if (jmin < jmax) {
125449b5e25fSSatish Balay           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
1255dc0b31edSSatish Balay           ierr = PetscLogFlops(2.0*(jmax-jmin));CHKERRQ(ierr);
1256187a9f4bSHong Zhang 
1257fff829cfSHong Zhang           /* update il and jl for row i */
1258fff829cfSHong Zhang           il[i] = jmin;
1259fff829cfSHong Zhang           j     = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
126049b5e25fSSatish Balay         }
1261ab27746eSHong Zhang         i = nexti;
126249b5e25fSSatish Balay       }
126349b5e25fSSatish Balay 
12643cea4cbeSHong Zhang       /* shift the diagonals when zero pivot is detected */
12653cea4cbeSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
12663cea4cbeSHong Zhang       rs   = 0.0;
1267997a0949SHong Zhang       jmin = bi[k]+1;
1268997a0949SHong Zhang       nz   = bi[k+1] - jmin;
1269997a0949SHong Zhang       if (nz) {
1270997a0949SHong Zhang         bcol = bj + jmin;
1271997a0949SHong Zhang         while (nz--) {
127289f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
127389f845c8SHong Zhang           bcol++;
1274997a0949SHong Zhang         }
1275997a0949SHong Zhang       }
12763cea4cbeSHong Zhang 
12773cea4cbeSHong Zhang       sctx.rs = rs;
12783cea4cbeSHong Zhang       sctx.pv = dk;
12794cccfbddSHong Zhang       ierr    = MatPivotCheck(C,A,info,&sctx,k);CHKERRQ(ierr);
128007b50cabSHong Zhang       if (sctx.newshift) break;    /* sctx.shift_amount is updated */
12810e95ead3SHong Zhang       dk = sctx.pv;
128249b5e25fSSatish Balay 
1283997a0949SHong Zhang       /* copy data into U(k,:) */
128498a8e78dSHong Zhang       ba[bi[k]] = 1.0/dk; /* U(k,k) */
1285fff829cfSHong Zhang       jmin      = bi[k]+1; jmax = bi[k+1];
128649b5e25fSSatish Balay       if (jmin < jmax) {
128749b5e25fSSatish Balay         for (j=jmin; j<jmax; j++) {
1288997a0949SHong Zhang           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
128949b5e25fSSatish Balay         }
129098a8e78dSHong Zhang         /* add the k-th row into il and jl */
129149b5e25fSSatish Balay         il[k] = jmin;
129298a8e78dSHong Zhang         i     = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
129349b5e25fSSatish Balay       }
129449b5e25fSSatish Balay     }
129507b50cabSHong Zhang   } while (sctx.newshift);
1296d8c74875SBarry Smith   ierr = PetscFree3(rtmp,il,jl);CHKERRQ(ierr);
129798a8e78dSHong Zhang   if (a->permute) {ierr = PetscFree(aa);CHKERRQ(ierr);}
129849b5e25fSSatish Balay 
129949b5e25fSSatish Balay   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
130026fbe8dcSKarl Rupp 
13010a3c351aSHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_1_inplace;
13024f79d315SHong Zhang   C->ops->solves         = MatSolves_SeqSBAIJ_1_inplace;
13030a3c351aSHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_1_inplace;
13040a3c351aSHong Zhang   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_1_inplace;
13050a3c351aSHong Zhang   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_1_inplace;
130649b5e25fSSatish Balay   C->assembled           = PETSC_TRUE;
13075c0bcdfcSHong Zhang   C->preallocated        = PETSC_TRUE;
130826fbe8dcSKarl Rupp 
1309d0f46423SBarry Smith   ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr);
13103cea4cbeSHong Zhang   if (sctx.nshift) {
1311f4db908eSBarry Smith     if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
131257622a8eSBarry Smith       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr);
1313f4db908eSBarry Smith     } else if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
131457622a8eSBarry Smith       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr);
1315997a0949SHong Zhang     }
1316997a0949SHong Zhang   }
131749b5e25fSSatish Balay   PetscFunctionReturn(0);
131849b5e25fSSatish Balay }
131949b5e25fSSatish Balay 
1320671cb588SHong Zhang /*
132180d3d5a7SHong Zhang   Version for when blocks are 1 by 1 Using natural ordering under new datastructure
132280d3d5a7SHong Zhang   Modified from MatCholeskyFactorNumeric_SeqAIJ()
1323671cb588SHong Zhang */
13244a2ae208SSatish Balay #undef __FUNCT__
13254a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering"
1326d595f711SHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info)
1327d595f711SHong Zhang {
1328d595f711SHong Zhang   Mat_SeqSBAIJ   *a=(Mat_SeqSBAIJ*)A->data;
13297b056e98SHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)B->data;
1330d595f711SHong Zhang   PetscErrorCode ierr;
1331d595f711SHong Zhang   PetscInt       i,j,mbs=A->rmap->n,*bi=b->i,*bj=b->j,*bdiag=b->diag,*bjtmp;
1332545dd064SHong Zhang   PetscInt       *ai=a->i,*aj=a->j,*ajtmp;
1333d595f711SHong Zhang   PetscInt       k,jmin,jmax,*c2r,*il,col,nexti,ili,nz;
1334d595f711SHong Zhang   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
1335d90ac83dSHong Zhang   FactorShiftCtx sctx;
1336d595f711SHong Zhang   PetscReal      rs;
1337d595f711SHong Zhang   MatScalar      d,*v;
1338d595f711SHong Zhang 
1339d595f711SHong Zhang   PetscFunctionBegin;
1340dcca6d9dSJed Brown   ierr = PetscMalloc3(mbs,&rtmp,mbs,&il,mbs,&c2r);CHKERRQ(ierr);
1341545dd064SHong Zhang 
1342d595f711SHong Zhang   /* MatPivotSetUp(): initialize shift context sctx */
1343d90ac83dSHong Zhang   ierr = PetscMemzero(&sctx,sizeof(FactorShiftCtx));CHKERRQ(ierr);
1344d595f711SHong Zhang 
1345f4db908eSBarry Smith   if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) { /* set sctx.shift_top=max{rs} */
1346d595f711SHong Zhang     sctx.shift_top = info->zeropivot;
134726fbe8dcSKarl Rupp 
1348545dd064SHong Zhang     ierr = PetscMemzero(rtmp,mbs*sizeof(MatScalar));CHKERRQ(ierr);
134926fbe8dcSKarl Rupp 
1350d595f711SHong Zhang     for (i=0; i<mbs; i++) {
1351d595f711SHong Zhang       /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
1352d595f711SHong Zhang       d        = (aa)[a->diag[i]];
1353545dd064SHong Zhang       rtmp[i] += -PetscRealPart(d);  /* diagonal entry */
1354545dd064SHong Zhang       ajtmp    = aj + ai[i] + 1;     /* exclude diagonal */
1355545dd064SHong Zhang       v        = aa + ai[i] + 1;
1356545dd064SHong Zhang       nz       = ai[i+1] - ai[i] - 1;
1357545dd064SHong Zhang       for (j=0; j<nz; j++) {
1358545dd064SHong Zhang         rtmp[i]        += PetscAbsScalar(v[j]);
1359545dd064SHong Zhang         rtmp[ajtmp[j]] += PetscAbsScalar(v[j]);
1360545dd064SHong Zhang       }
13610838c725SBarry Smith       if (PetscRealPart(rtmp[i]) > sctx.shift_top) sctx.shift_top = PetscRealPart(rtmp[i]);
1362d595f711SHong Zhang     }
1363d595f711SHong Zhang     sctx.shift_top *= 1.1;
1364d595f711SHong Zhang     sctx.nshift_max = 5;
1365d595f711SHong Zhang     sctx.shift_lo   = 0.;
1366d595f711SHong Zhang     sctx.shift_hi   = 1.;
1367d595f711SHong Zhang   }
1368d595f711SHong Zhang 
1369d595f711SHong Zhang   /* allocate working arrays
1370d595f711SHong Zhang      c2r: linked list, keep track of pivot rows for a given column. c2r[col]: head of the list for a given col
1371d595f711SHong 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
1372d595f711SHong Zhang   */
1373d595f711SHong Zhang   do {
137407b50cabSHong Zhang     sctx.newshift = PETSC_FALSE;
1375d595f711SHong Zhang 
1376d595f711SHong Zhang     for (i=0; i<mbs; i++) c2r[i] = mbs;
13779e95ef84SSatish Balay     if (mbs) il[0] = 0;
1378d595f711SHong Zhang 
1379d595f711SHong Zhang     for (k = 0; k<mbs; k++) {
1380d595f711SHong Zhang       /* zero rtmp */
1381d595f711SHong Zhang       nz    = bi[k+1] - bi[k];
1382d595f711SHong Zhang       bjtmp = bj + bi[k];
13837b056e98SHong Zhang       for (j=0; j<nz; j++) rtmp[bjtmp[j]] = 0.0;
1384d595f711SHong Zhang 
1385d595f711SHong Zhang       /* load in initial unfactored row */
1386d595f711SHong Zhang       bval = ba + bi[k];
1387d595f711SHong Zhang       jmin = ai[k]; jmax = ai[k+1];
1388d595f711SHong Zhang       for (j = jmin; j < jmax; j++) {
1389d595f711SHong Zhang         col       = aj[j];
1390d595f711SHong Zhang         rtmp[col] = aa[j];
1391d595f711SHong Zhang         *bval++   = 0.0; /* for in-place factorization */
1392d595f711SHong Zhang       }
1393d595f711SHong Zhang       /* shift the diagonal of the matrix: ZeropivotApply() */
1394d595f711SHong Zhang       rtmp[k] += sctx.shift_amount;  /* shift the diagonal of the matrix */
1395d595f711SHong Zhang 
1396d595f711SHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
1397d595f711SHong Zhang       dk = rtmp[k];
1398d595f711SHong Zhang       i  = c2r[k]; /* first row to be added to k_th row  */
1399d595f711SHong Zhang 
1400d595f711SHong Zhang       while (i < k) {
1401d595f711SHong Zhang         nexti = c2r[i]; /* next row to be added to k_th row */
1402d595f711SHong Zhang 
1403d595f711SHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
1404d595f711SHong Zhang         ili     = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
1405d595f711SHong Zhang         uikdi   = -ba[ili]*ba[bdiag[i]]; /* diagonal(k) */
1406d595f711SHong Zhang         dk     += uikdi*ba[ili]; /* update diag[k] */
1407d595f711SHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
1408d595f711SHong Zhang 
1409d595f711SHong Zhang         /* add multiple of row i to k-th row */
1410d595f711SHong Zhang         jmin = ili + 1; jmax = bi[i+1];
1411d595f711SHong Zhang         if (jmin < jmax) {
1412d595f711SHong Zhang           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
1413d595f711SHong Zhang           /* update il and c2r for row i */
1414d595f711SHong Zhang           il[i] = jmin;
1415d595f711SHong Zhang           j     = bj[jmin]; c2r[i] = c2r[j]; c2r[j] = i;
1416d595f711SHong Zhang         }
1417d595f711SHong Zhang         i = nexti;
1418d595f711SHong Zhang       }
1419d595f711SHong Zhang 
1420d595f711SHong Zhang       /* copy data into U(k,:) */
1421d595f711SHong Zhang       rs   = 0.0;
1422d595f711SHong Zhang       jmin = bi[k]; jmax = bi[k+1]-1;
1423d595f711SHong Zhang       if (jmin < jmax) {
1424d595f711SHong Zhang         for (j=jmin; j<jmax; j++) {
1425d595f711SHong Zhang           col = bj[j]; ba[j] = rtmp[col]; rs += PetscAbsScalar(ba[j]);
1426d595f711SHong Zhang         }
1427d595f711SHong Zhang         /* add the k-th row into il and c2r */
1428d595f711SHong Zhang         il[k] = jmin;
1429d595f711SHong Zhang         i     = bj[jmin]; c2r[k] = c2r[i]; c2r[i] = k;
1430d595f711SHong Zhang       }
1431d595f711SHong Zhang 
1432d595f711SHong Zhang       sctx.rs = rs;
1433d595f711SHong Zhang       sctx.pv = dk;
14344cccfbddSHong Zhang       ierr    = MatPivotCheck(B,A,info,&sctx,k);CHKERRQ(ierr);
143507b50cabSHong Zhang       if (sctx.newshift) break;
1436d595f711SHong Zhang       dk = sctx.pv;
1437d595f711SHong Zhang 
1438d595f711SHong Zhang       ba[bdiag[k]] = 1.0/dk; /* U(k,k) */
1439d595f711SHong Zhang     }
144007b50cabSHong Zhang   } while (sctx.newshift);
1441d595f711SHong Zhang 
1442d595f711SHong Zhang   ierr = PetscFree3(rtmp,il,c2r);CHKERRQ(ierr);
1443d595f711SHong Zhang 
1444d595f711SHong Zhang   B->ops->solve          = MatSolve_SeqSBAIJ_1_NaturalOrdering;
14459dcdb97aSHong Zhang   B->ops->solves         = MatSolves_SeqSBAIJ_1;
1446d595f711SHong Zhang   B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering;
144780d3d5a7SHong Zhang   B->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering;
144880d3d5a7SHong Zhang   B->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering;
1449d595f711SHong Zhang 
14507b056e98SHong Zhang   B->assembled    = PETSC_TRUE;
14517b056e98SHong Zhang   B->preallocated = PETSC_TRUE;
145226fbe8dcSKarl Rupp 
14537b056e98SHong Zhang   ierr = PetscLogFlops(B->rmap->n);CHKERRQ(ierr);
1454d595f711SHong Zhang 
1455d595f711SHong Zhang   /* MatPivotView() */
1456d595f711SHong Zhang   if (sctx.nshift) {
1457f4db908eSBarry Smith     if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
145857622a8eSBarry Smith       ierr = PetscInfo4(A,"number of shift_pd tries %D, shift_amount %g, diagonal shifted up by %e fraction top_value %e\n",sctx.nshift,(double)sctx.shift_amount,(double)sctx.shift_fraction,(double)sctx.shift_top);CHKERRQ(ierr);
1459f4db908eSBarry Smith     } else if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
146057622a8eSBarry Smith       ierr = PetscInfo2(A,"number of shift_nz tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr);
1461f4db908eSBarry Smith     } else if (info->shifttype == (PetscReal)MAT_SHIFT_INBLOCKS) {
146257622a8eSBarry Smith       ierr = PetscInfo2(A,"number of shift_inblocks applied %D, each shift_amount %g\n",sctx.nshift,(double)info->shiftamount);CHKERRQ(ierr);
1463d595f711SHong Zhang     }
1464d595f711SHong Zhang   }
1465d595f711SHong Zhang   PetscFunctionReturn(0);
1466d595f711SHong Zhang }
1467d595f711SHong Zhang 
1468d595f711SHong Zhang #undef __FUNCT__
1469d595f711SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering_inplace"
1470d595f711SHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering_inplace(Mat C,Mat A,const MatFactorInfo *info)
1471671cb588SHong Zhang {
1472671cb588SHong Zhang   Mat_SeqSBAIJ   *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ*)C->data;
1473dfbe8321SBarry Smith   PetscErrorCode ierr;
147413f74950SBarry Smith   PetscInt       i,j,mbs = a->mbs;
147513f74950SBarry Smith   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
14763cea4cbeSHong Zhang   PetscInt       k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz;
1477653a6975SHong Zhang   MatScalar      *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval;
147814d2772eSHong Zhang   PetscReal      rs;
14790e95ead3SHong Zhang   FactorShiftCtx sctx;
1480653a6975SHong Zhang 
1481653a6975SHong Zhang   PetscFunctionBegin;
14820e95ead3SHong Zhang   /* MatPivotSetUp(): initialize shift context sctx */
14830e95ead3SHong Zhang   ierr = PetscMemzero(&sctx,sizeof(FactorShiftCtx));CHKERRQ(ierr);
14840e95ead3SHong Zhang 
1485653a6975SHong Zhang   /* initialization */
1486653a6975SHong Zhang   /* il and jl record the first nonzero element in each row of the accessing
1487653a6975SHong Zhang      window U(0:k, k:mbs-1).
1488653a6975SHong Zhang      jl:    list of rows to be added to uneliminated rows
1489653a6975SHong Zhang             i>= k: jl(i) is the first row to be added to row i
1490653a6975SHong Zhang             i<  k: jl(i) is the row following row i in some list of rows
1491653a6975SHong Zhang             jl(i) = mbs indicates the end of a list
1492653a6975SHong Zhang      il(i): points to the first nonzero element in U(i,k:mbs-1)
1493653a6975SHong Zhang   */
1494785e854fSJed Brown   ierr = PetscMalloc1(mbs,&rtmp);CHKERRQ(ierr);
1495dcca6d9dSJed Brown   ierr = PetscMalloc2(mbs,&il,mbs,&jl);CHKERRQ(ierr);
1496b00f7748SHong Zhang 
1497b00f7748SHong Zhang   do {
149807b50cabSHong Zhang     sctx.newshift = PETSC_FALSE;
14996df5ee2eSHong Zhang     il[0] = 0;
1500653a6975SHong Zhang     for (i=0; i<mbs; i++) {
15016df5ee2eSHong Zhang       rtmp[i] = 0.0; jl[i] = mbs;
1502653a6975SHong Zhang     }
1503653a6975SHong Zhang 
1504997a0949SHong Zhang     for (k = 0; k<mbs; k++) {
1505653a6975SHong Zhang       /*initialize k-th row with elements nonzero in row perm(k) of A */
1506653a6975SHong Zhang       nz   = ai[k+1] - ai[k];
1507653a6975SHong Zhang       acol = aj + ai[k];
1508653a6975SHong Zhang       aval = aa + ai[k];
1509653a6975SHong Zhang       bval = ba + bi[k];
1510653a6975SHong Zhang       while (nz--) {
1511653a6975SHong Zhang         rtmp[*acol++] = *aval++;
1512653a6975SHong Zhang         *bval++       = 0.0; /* for in-place factorization */
1513653a6975SHong Zhang       }
15143cea4cbeSHong Zhang 
15153cea4cbeSHong Zhang       /* shift the diagonal of the matrix */
15163cea4cbeSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
1517653a6975SHong Zhang 
1518653a6975SHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
1519653a6975SHong Zhang       dk = rtmp[k];
1520653a6975SHong Zhang       i  = jl[k]; /* first row to be added to k_th row  */
1521653a6975SHong Zhang 
1522653a6975SHong Zhang       while (i < k) {
1523653a6975SHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
1524653a6975SHong Zhang         /* compute multiplier, update D(k) and U(i,k) */
1525653a6975SHong Zhang         ili     = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
1526653a6975SHong Zhang         uikdi   = -ba[ili]*ba[bi[i]];
1527653a6975SHong Zhang         dk     += uikdi*ba[ili];
1528653a6975SHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
1529653a6975SHong Zhang 
1530653a6975SHong Zhang         /* add multiple of row i to k-th row ... */
1531653a6975SHong Zhang         jmin = ili + 1;
1532653a6975SHong Zhang         nz   = bi[i+1] - jmin;
1533653a6975SHong Zhang         if (nz > 0) {
1534653a6975SHong Zhang           bcol = bj + jmin;
1535653a6975SHong Zhang           bval = ba + jmin;
1536dc0b31edSSatish Balay           ierr = PetscLogFlops(2.0*nz);CHKERRQ(ierr);
1537653a6975SHong Zhang           while (nz--) rtmp[*bcol++] += uikdi*(*bval++);
1538187a9f4bSHong Zhang 
1539997a0949SHong Zhang           /* update il and jl for i-th row */
1540997a0949SHong Zhang           il[i] = jmin;
1541997a0949SHong Zhang           j     = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
1542653a6975SHong Zhang         }
1543653a6975SHong Zhang         i = nexti;
1544653a6975SHong Zhang       }
1545653a6975SHong Zhang 
15463cea4cbeSHong Zhang       /* shift the diagonals when zero pivot is detected */
15473cea4cbeSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
15483cea4cbeSHong Zhang       rs   = 0.0;
154921c26570Svictorle       jmin = bi[k]+1;
155021c26570Svictorle       nz   = bi[k+1] - jmin;
155121c26570Svictorle       if (nz) {
155221c26570Svictorle         bcol = bj + jmin;
155321c26570Svictorle         while (nz--) {
155489f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
155589f845c8SHong Zhang           bcol++;
155621c26570Svictorle         }
155721c26570Svictorle       }
15583cea4cbeSHong Zhang 
15593cea4cbeSHong Zhang       sctx.rs = rs;
15603cea4cbeSHong Zhang       sctx.pv = dk;
15614cccfbddSHong Zhang       ierr    = MatPivotCheck(C,A,info,&sctx,k);CHKERRQ(ierr);
156207b50cabSHong Zhang       if (sctx.newshift) break;    /* sctx.shift_amount is updated */
15630e95ead3SHong Zhang       dk = sctx.pv;
1564653a6975SHong Zhang 
1565997a0949SHong Zhang       /* copy data into U(k,:) */
1566653a6975SHong Zhang       ba[bi[k]] = 1.0/dk;
1567653a6975SHong Zhang       jmin      = bi[k]+1;
1568653a6975SHong Zhang       nz        = bi[k+1] - jmin;
1569653a6975SHong Zhang       if (nz) {
1570653a6975SHong Zhang         bcol = bj + jmin;
1571653a6975SHong Zhang         bval = ba + jmin;
1572653a6975SHong Zhang         while (nz--) {
1573653a6975SHong Zhang           *bval++       = rtmp[*bcol];
1574653a6975SHong Zhang           rtmp[*bcol++] = 0.0;
1575653a6975SHong Zhang         }
1576997a0949SHong Zhang         /* add k-th row into il and jl */
1577653a6975SHong Zhang         il[k] = jmin;
1578997a0949SHong Zhang         i     = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
1579653a6975SHong Zhang       }
1580b00f7748SHong Zhang     } /* end of for (k = 0; k<mbs; k++) */
158107b50cabSHong Zhang   } while (sctx.newshift);
1582653a6975SHong Zhang   ierr = PetscFree(rtmp);CHKERRQ(ierr);
1583d8c74875SBarry Smith   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
1584653a6975SHong Zhang 
15850a3c351aSHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
15864f79d315SHong Zhang   C->ops->solves         = MatSolves_SeqSBAIJ_1_inplace;
15870a3c351aSHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
15880a3c351aSHong Zhang   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
15890a3c351aSHong Zhang   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
1590db4efbfdSBarry Smith 
1591653a6975SHong Zhang   C->assembled    = PETSC_TRUE;
1592653a6975SHong Zhang   C->preallocated = PETSC_TRUE;
159326fbe8dcSKarl Rupp 
1594d0f46423SBarry Smith   ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr);
15953cea4cbeSHong Zhang   if (sctx.nshift) {
1596f4db908eSBarry Smith     if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
159757622a8eSBarry Smith       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr);
1598f4db908eSBarry Smith     } else if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
159957622a8eSBarry Smith       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr);
1600b00f7748SHong Zhang     }
1601fb10cecfSBarry Smith   }
1602653a6975SHong Zhang   PetscFunctionReturn(0);
1603653a6975SHong Zhang }
1604653a6975SHong Zhang 
1605653a6975SHong Zhang #undef __FUNCT__
16064a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactor_SeqSBAIJ"
16070481f469SBarry Smith PetscErrorCode MatCholeskyFactor_SeqSBAIJ(Mat A,IS perm,const MatFactorInfo *info)
160849b5e25fSSatish Balay {
1609dfbe8321SBarry Smith   PetscErrorCode ierr;
161049b5e25fSSatish Balay   Mat            C;
161149b5e25fSSatish Balay 
161249b5e25fSSatish Balay   PetscFunctionBegin;
1613719d5645SBarry Smith   ierr = MatGetFactor(A,"petsc",MAT_FACTOR_CHOLESKY,&C);CHKERRQ(ierr);
1614719d5645SBarry Smith   ierr = MatCholeskyFactorSymbolic(C,A,perm,info);CHKERRQ(ierr);
1615719d5645SBarry Smith   ierr = MatCholeskyFactorNumeric(C,A,info);CHKERRQ(ierr);
161626fbe8dcSKarl Rupp 
1617db4efbfdSBarry Smith   A->ops->solve          = C->ops->solve;
1618db4efbfdSBarry Smith   A->ops->solvetranspose = C->ops->solvetranspose;
161926fbe8dcSKarl Rupp 
162028be2f97SBarry Smith   ierr = MatHeaderMerge(A,&C);CHKERRQ(ierr);
162149b5e25fSSatish Balay   PetscFunctionReturn(0);
162249b5e25fSSatish Balay }
162349b5e25fSSatish Balay 
162449b5e25fSSatish Balay 
1625