xref: /petsc/src/mat/impls/sbaij/seq/sbaijfact.c (revision 9dcdb97a1f5de371b1590f3f46a2e86791a0ad4d)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
249b5e25fSSatish Balay 
37c4f633dSBarry Smith #include "../src/mat/impls/baij/seq/baij.h"
47c4f633dSBarry Smith #include "../src/mat/impls/sbaij/seq/sbaij.h"
5c60f0209SBarry Smith #include "../src/mat/blockinvert.h"
67c4f633dSBarry Smith #include "petscis.h"
749b5e25fSSatish Balay 
88dc9108eSHong Zhang #if !defined(PETSC_USE_COMPLEX)
95f9f512dSHong Zhang /*
105f9f512dSHong Zhang   input:
11c037c3f7SHong Zhang    F -- numeric factor
125f9f512dSHong Zhang   output:
13c037c3f7SHong Zhang    nneg, nzero, npos: matrix inertia
145f9f512dSHong Zhang */
155f9f512dSHong Zhang 
165f9f512dSHong Zhang #undef __FUNCT__
175f9f512dSHong Zhang #define __FUNCT__ "MatGetInertia_SeqSBAIJ"
1813f74950SBarry Smith PetscErrorCode MatGetInertia_SeqSBAIJ(Mat F,PetscInt *nneig,PetscInt *nzero,PetscInt *npos)
195f9f512dSHong Zhang {
20638f5ce0SDinesh Kaushik   Mat_SeqSBAIJ *fact_ptr = (Mat_SeqSBAIJ*)F->data;
21dd6ea824SBarry Smith   MatScalar    *dd = fact_ptr->a;
22d0f46423SBarry Smith   PetscInt     mbs=fact_ptr->mbs,bs=F->rmap->bs,i,nneig_tmp,npos_tmp,*fi = fact_ptr->i;
235f9f512dSHong Zhang 
245f9f512dSHong Zhang   PetscFunctionBegin;
2577431f27SBarry Smith   if (bs != 1) SETERRQ1(PETSC_ERR_SUP,"No support for bs: %D >1 yet",bs);
26eeeff2ecSHong Zhang   nneig_tmp = 0; npos_tmp = 0;
27eeeff2ecSHong Zhang   for (i=0; i<mbs; i++){
28eeeff2ecSHong Zhang     if (PetscRealPart(dd[*fi]) > 0.0){
29eeeff2ecSHong Zhang       npos_tmp++;
30eeeff2ecSHong Zhang     } else if (PetscRealPart(dd[*fi]) < 0.0){
31eeeff2ecSHong Zhang       nneig_tmp++;
325f9f512dSHong Zhang     }
33eeeff2ecSHong Zhang     fi++;
343e0d88b5SBarry Smith   }
35eeeff2ecSHong Zhang   if (nneig) *nneig = nneig_tmp;
36eeeff2ecSHong Zhang   if (npos)  *npos  = npos_tmp;
37eeeff2ecSHong Zhang   if (nzero) *nzero = mbs - nneig_tmp - npos_tmp;
38eeeff2ecSHong Zhang 
395f9f512dSHong Zhang   PetscFunctionReturn(0);
405f9f512dSHong Zhang }
418dc9108eSHong Zhang #endif /* !defined(PETSC_USE_COMPLEX) */
425f9f512dSHong Zhang 
435f9f512dSHong Zhang /*
445f9f512dSHong Zhang   Symbolic U^T*D*U factorization for SBAIJ format. Modified from SSF of YSMP.
4510c27e3dSHong Zhang   Use Modified Sparse Row (MSR) storage for u and ju. See page 85, "Iterative Methods ..." by Saad.
465f9f512dSHong Zhang */
4710c27e3dSHong Zhang #undef __FUNCT__
4810c27e3dSHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqSBAIJ_MSR"
490481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ_MSR(Mat F,Mat A,IS perm,const MatFactorInfo *info)
5010c27e3dSHong Zhang {
5110c27e3dSHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b;
5210c27e3dSHong Zhang   PetscErrorCode ierr;
535d0c19d7SBarry Smith   const PetscInt *rip,*ai,*aj;
545d0c19d7SBarry Smith   PetscInt       i,mbs = a->mbs,*jutmp,bs = A->rmap->bs,bs2=a->bs2;
5510c27e3dSHong Zhang   PetscInt       m,reallocs = 0,prow;
5610c27e3dSHong Zhang   PetscInt       *jl,*q,jmin,jmax,juidx,nzk,qm,*iu,*ju,k,j,vj,umax,maxadd;
5710c27e3dSHong Zhang   PetscReal      f = info->fill;
5810c27e3dSHong Zhang   PetscTruth     perm_identity;
5910c27e3dSHong Zhang 
6010c27e3dSHong Zhang   PetscFunctionBegin;
6110c27e3dSHong Zhang   /* check whether perm is the identity mapping */
6210c27e3dSHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
6310c27e3dSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
6410c27e3dSHong Zhang 
6510c27e3dSHong Zhang   if (perm_identity){ /* without permutation */
6610c27e3dSHong Zhang     a->permute = PETSC_FALSE;
6710c27e3dSHong Zhang     ai = a->i; aj = a->j;
6810c27e3dSHong Zhang   } else {            /* non-trivial permutation */
6910c27e3dSHong Zhang     a->permute = PETSC_TRUE;
7010c27e3dSHong Zhang     ierr = MatReorderingSeqSBAIJ(A,perm);CHKERRQ(ierr);
7110c27e3dSHong Zhang     ai = a->inew; aj = a->jnew;
7210c27e3dSHong Zhang   }
7310c27e3dSHong Zhang 
7410c27e3dSHong Zhang   /* initialization */
7510c27e3dSHong Zhang   ierr  = PetscMalloc((mbs+1)*sizeof(PetscInt),&iu);CHKERRQ(ierr);
7610c27e3dSHong Zhang   umax  = (PetscInt)(f*ai[mbs] + 1); umax += mbs + 1;
7710c27e3dSHong Zhang   ierr  = PetscMalloc(umax*sizeof(PetscInt),&ju);CHKERRQ(ierr);
7810c27e3dSHong Zhang   iu[0] = mbs+1;
7910c27e3dSHong Zhang   juidx = mbs + 1; /* index for ju */
80d8c74875SBarry Smith   /* jl linked list for pivot row -- linked list for col index */
81d8c74875SBarry Smith   ierr  = PetscMalloc2(mbs,PetscInt,&jl,mbs,PetscInt,&q);CHKERRQ(ierr);
8210c27e3dSHong Zhang   for (i=0; i<mbs; i++){
8310c27e3dSHong Zhang     jl[i] = mbs;
8410c27e3dSHong Zhang     q[i] = 0;
8510c27e3dSHong Zhang   }
8610c27e3dSHong Zhang 
8710c27e3dSHong Zhang   /* for each row k */
8810c27e3dSHong Zhang   for (k=0; k<mbs; k++){
8910c27e3dSHong Zhang     for (i=0; i<mbs; i++) q[i] = 0;  /* to be removed! */
9010c27e3dSHong Zhang     nzk  = 0; /* num. of nz blocks in k-th block row with diagonal block excluded */
9110c27e3dSHong Zhang     q[k] = mbs;
9210c27e3dSHong Zhang     /* initialize nonzero structure of k-th row to row rip[k] of A */
9310c27e3dSHong Zhang     jmin = ai[rip[k]] +1; /* exclude diag[k] */
9410c27e3dSHong Zhang     jmax = ai[rip[k]+1];
9510c27e3dSHong Zhang     for (j=jmin; j<jmax; j++){
9610c27e3dSHong Zhang       vj = rip[aj[j]]; /* col. value */
9710c27e3dSHong Zhang       if(vj > k){
9810c27e3dSHong Zhang         qm = k;
9910c27e3dSHong Zhang         do {
10010c27e3dSHong Zhang           m  = qm; qm = q[m];
10110c27e3dSHong Zhang         } while(qm < vj);
10210c27e3dSHong Zhang         if (qm == vj) {
10310c27e3dSHong Zhang           SETERRQ(PETSC_ERR_PLIB,"Duplicate entry in A\n");
10410c27e3dSHong Zhang         }
10510c27e3dSHong Zhang         nzk++;
10610c27e3dSHong Zhang         q[m]  = vj;
10710c27e3dSHong Zhang         q[vj] = qm;
10810c27e3dSHong Zhang       } /* if(vj > k) */
10910c27e3dSHong Zhang     } /* for (j=jmin; j<jmax; j++) */
11010c27e3dSHong Zhang 
11110c27e3dSHong Zhang     /* modify nonzero structure of k-th row by computing fill-in
11210c27e3dSHong Zhang        for each row i to be merged in */
11310c27e3dSHong Zhang     prow = k;
11410c27e3dSHong Zhang     prow = jl[prow]; /* next pivot row (== mbs for symbolic factorization) */
11510c27e3dSHong Zhang 
11610c27e3dSHong Zhang     while (prow < k){
11710c27e3dSHong Zhang       /* merge row prow into k-th row */
11810c27e3dSHong Zhang       jmin = iu[prow] + 1; jmax = iu[prow+1];
11910c27e3dSHong Zhang       qm = k;
12010c27e3dSHong Zhang       for (j=jmin; j<jmax; j++){
12110c27e3dSHong Zhang         vj = ju[j];
12210c27e3dSHong Zhang         do {
12310c27e3dSHong Zhang           m = qm; qm = q[m];
12410c27e3dSHong Zhang         } while (qm < vj);
12510c27e3dSHong Zhang         if (qm != vj){
12610c27e3dSHong Zhang          nzk++; q[m] = vj; q[vj] = qm; qm = vj;
12710c27e3dSHong Zhang         }
12810c27e3dSHong Zhang       }
12910c27e3dSHong Zhang       prow = jl[prow]; /* next pivot row */
13010c27e3dSHong Zhang     }
13110c27e3dSHong Zhang 
13210c27e3dSHong Zhang     /* add k to row list for first nonzero element in k-th row */
13310c27e3dSHong Zhang     if (nzk > 0){
13410c27e3dSHong Zhang       i = q[k]; /* col value of first nonzero element in U(k, k+1:mbs-1) */
13510c27e3dSHong Zhang       jl[k] = jl[i]; jl[i] = k;
13610c27e3dSHong Zhang     }
13710c27e3dSHong Zhang     iu[k+1] = iu[k] + nzk;
13810c27e3dSHong Zhang 
13910c27e3dSHong Zhang     /* allocate more space to ju if needed */
14010c27e3dSHong Zhang     if (iu[k+1] > umax) {
14110c27e3dSHong Zhang       /* estimate how much additional space we will need */
14210c27e3dSHong Zhang       /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */
14310c27e3dSHong Zhang       /* just double the memory each time */
14410c27e3dSHong Zhang       maxadd = umax;
14510c27e3dSHong Zhang       if (maxadd < nzk) maxadd = (mbs-k)*(nzk+1)/2;
14610c27e3dSHong Zhang       umax += maxadd;
14710c27e3dSHong Zhang 
14810c27e3dSHong Zhang       /* allocate a longer ju */
14910c27e3dSHong Zhang       ierr = PetscMalloc(umax*sizeof(PetscInt),&jutmp);CHKERRQ(ierr);
15010c27e3dSHong Zhang       ierr = PetscMemcpy(jutmp,ju,iu[k]*sizeof(PetscInt));CHKERRQ(ierr);
15110c27e3dSHong Zhang       ierr = PetscFree(ju);CHKERRQ(ierr);
15210c27e3dSHong Zhang       ju   = jutmp;
15310c27e3dSHong Zhang       reallocs++; /* count how many times we realloc */
15410c27e3dSHong Zhang     }
15510c27e3dSHong Zhang 
15610c27e3dSHong Zhang     /* save nonzero structure of k-th row in ju */
15710c27e3dSHong Zhang     i=k;
15810c27e3dSHong Zhang     while (nzk --) {
15910c27e3dSHong Zhang       i           = q[i];
16010c27e3dSHong Zhang       ju[juidx++] = i;
16110c27e3dSHong Zhang     }
16210c27e3dSHong Zhang   }
16310c27e3dSHong Zhang 
1646cf91177SBarry Smith #if defined(PETSC_USE_INFO)
16510c27e3dSHong Zhang   if (ai[mbs] != 0) {
16610c27e3dSHong Zhang     PetscReal af = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]);
167ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr);
168ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
169ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr);
170ae15b995SBarry Smith     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
17110c27e3dSHong Zhang   } else {
172ae15b995SBarry Smith     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
17310c27e3dSHong Zhang   }
17463ba0a88SBarry Smith #endif
17510c27e3dSHong Zhang 
17610c27e3dSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
177d8c74875SBarry Smith   ierr = PetscFree2(jl,q);CHKERRQ(ierr);
17810c27e3dSHong Zhang 
17910c27e3dSHong Zhang   /* put together the new matrix */
180719d5645SBarry Smith   ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(F,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
18110c27e3dSHong Zhang 
182719d5645SBarry Smith   /* ierr = PetscLogObjectParent(B,iperm);CHKERRQ(ierr); */
183719d5645SBarry Smith   b = (Mat_SeqSBAIJ*)(F)->data;
18410c27e3dSHong Zhang   b->singlemalloc = PETSC_FALSE;
185e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
186e6b907acSBarry Smith   b->free_ij       = PETSC_TRUE;
18710c27e3dSHong Zhang   ierr = PetscMalloc((iu[mbs]+1)*sizeof(MatScalar)*bs2,&b->a);CHKERRQ(ierr);
18810c27e3dSHong Zhang   b->j    = ju;
18910c27e3dSHong Zhang   b->i    = iu;
19010c27e3dSHong Zhang   b->diag = 0;
19110c27e3dSHong Zhang   b->ilen = 0;
19210c27e3dSHong Zhang   b->imax = 0;
19310c27e3dSHong Zhang   b->row  = perm;
19410c27e3dSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
19510c27e3dSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
19610c27e3dSHong Zhang   b->icol = perm;
19710c27e3dSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
19810c27e3dSHong Zhang   ierr    = PetscMalloc((bs*mbs+bs)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
19910c27e3dSHong Zhang   /* In b structure:  Free imax, ilen, old a, old j.
20010c27e3dSHong Zhang      Allocate idnew, solve_work, new a, new j */
201719d5645SBarry Smith   ierr = PetscLogObjectMemory(F,(iu[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
20210c27e3dSHong Zhang   b->maxnz = b->nz = iu[mbs];
20310c27e3dSHong Zhang 
204719d5645SBarry Smith   (F)->info.factor_mallocs    = reallocs;
205719d5645SBarry Smith   (F)->info.fill_ratio_given  = f;
20610c27e3dSHong Zhang   if (ai[mbs] != 0) {
207719d5645SBarry Smith     (F)->info.fill_ratio_needed = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]);
20810c27e3dSHong Zhang   } else {
209719d5645SBarry Smith     (F)->info.fill_ratio_needed = 0.0;
21010c27e3dSHong Zhang   }
211d595f711SHong Zhang   ierr = MatSeqSBAIJSetNumericFactorization_inplace(F,perm_identity);CHKERRQ(ierr);
21210c27e3dSHong Zhang   PetscFunctionReturn(0);
21310c27e3dSHong Zhang }
21410c27e3dSHong Zhang /*
21510c27e3dSHong Zhang     Symbolic U^T*D*U factorization for SBAIJ format.
21610c27e3dSHong Zhang */
21798a8e78dSHong Zhang #include "petscbt.h"
2187c4f633dSBarry Smith #include "../src/mat/utils/freespace.h"
2194a2ae208SSatish Balay #undef __FUNCT__
2204a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqSBAIJ"
2210481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
22249b5e25fSSatish Balay {
22398a8e78dSHong Zhang   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data;
22498a8e78dSHong Zhang   Mat_SeqSBAIJ       *b;
2256849ba73SBarry Smith   PetscErrorCode     ierr;
22658ebbce7SBarry Smith   PetscTruth         perm_identity,missing;
22798a8e78dSHong Zhang   PetscReal          fill = info->fill;
2287b056e98SHong Zhang   const PetscInt     *rip,*ai=a->i,*aj=a->j;
2295d0c19d7SBarry Smith   PetscInt           i,mbs=a->mbs,bs=A->rmap->bs,reallocs=0,prow,d;
23098a8e78dSHong Zhang   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
231c6d0d4f0SHong Zhang   PetscInt           nlnk,*lnk,ncols,*cols,*uj,**ui_ptr,*uj_ptr,*udiag;
232a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
23398a8e78dSHong Zhang   PetscBT            lnkbt;
234d595f711SHong Zhang   PetscTruth         olddatastruct=PETSC_FALSE;
235d595f711SHong Zhang 
236d595f711SHong Zhang   PetscFunctionBegin;
237d595f711SHong Zhang     ierr = PetscOptionsGetTruth(PETSC_NULL,"-cholesky_old",&olddatastruct,PETSC_NULL);CHKERRQ(ierr);
238d595f711SHong Zhang   if (olddatastruct || bs>1 ){
239d595f711SHong Zhang     ierr = MatCholeskyFactorSymbolic_SeqSBAIJ_inplace(fact,A,perm,info);CHKERRQ(ierr);
240d595f711SHong Zhang     PetscFunctionReturn(0);
241d595f711SHong Zhang   }
242d595f711SHong Zhang 
243c6d0d4f0SHong Zhang   if (A->rmap->n != A->cmap->n) SETERRQ2(PETSC_ERR_ARG_WRONG,"Must be square matrix, rows %D columns %D",A->rmap->n,A->cmap->n);
244d595f711SHong Zhang   ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr);
245d595f711SHong Zhang   if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
246d595f711SHong Zhang 
247d595f711SHong Zhang   /* check whether perm is the identity mapping */
248d595f711SHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
249c6d0d4f0SHong Zhang   if (!perm_identity) SETERRQ(PETSC_ERR_SUP,"Matrix reordering is not supported for sbaij matrix. Use aij format");
250d595f711SHong Zhang   a->permute = PETSC_FALSE;
251d595f711SHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
252d595f711SHong Zhang 
253d595f711SHong Zhang   /* initialization */
254d595f711SHong Zhang   ierr  = PetscMalloc((mbs+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
255c6d0d4f0SHong Zhang   ierr  = PetscMalloc((mbs+1)*sizeof(PetscInt),&udiag);CHKERRQ(ierr);
256d595f711SHong Zhang   ui[0] = 0;
257d595f711SHong Zhang 
258d595f711SHong Zhang   /* jl: linked list for storing indices of the pivot rows
259d595f711SHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
260d595f711SHong Zhang   ierr = PetscMalloc4(mbs,PetscInt*,&ui_ptr,mbs,PetscInt,&il,mbs,PetscInt,&jl,mbs,PetscInt,&cols);CHKERRQ(ierr);
261d595f711SHong Zhang   for (i=0; i<mbs; i++){
262d595f711SHong Zhang     jl[i] = mbs; il[i] = 0;
263d595f711SHong Zhang   }
264d595f711SHong Zhang 
265d595f711SHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
266d595f711SHong Zhang   nlnk = mbs + 1;
267d595f711SHong Zhang   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
268d595f711SHong Zhang 
269d595f711SHong Zhang   /* initial FreeSpace size is fill*(ai[mbs]+1) */
270d595f711SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[mbs]+1)),&free_space);CHKERRQ(ierr);
271d595f711SHong Zhang   current_space = free_space;
272d595f711SHong Zhang 
273d595f711SHong Zhang   for (k=0; k<mbs; k++){  /* for each active row k */
274d595f711SHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
275d595f711SHong Zhang     nzk   = 0;
276c6d0d4f0SHong Zhang     ncols = ai[k+1] - ai[k];
277c6d0d4f0SHong Zhang     if (!ncols) SETERRQ1(PETSC_ERR_MAT_CH_ZRPVT,"Empty row %D in matrix ",k);
278d595f711SHong Zhang     for (j=0; j<ncols; j++){
279c6d0d4f0SHong Zhang       i = *(aj + ai[k] + j);
280c6d0d4f0SHong Zhang       cols[j] = i;
281d595f711SHong Zhang     }
282d595f711SHong Zhang     ierr = PetscLLAdd(ncols,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
283d595f711SHong Zhang     nzk += nlnk;
284d595f711SHong Zhang 
285d595f711SHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
286d595f711SHong Zhang     prow = jl[k]; /* 1st pivot row */
287d595f711SHong Zhang 
288d595f711SHong Zhang     while (prow < k){
289d595f711SHong Zhang       nextprow = jl[prow];
290d595f711SHong Zhang       /* merge prow into k-th row */
291d595f711SHong Zhang       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
292d595f711SHong Zhang       jmax = ui[prow+1];
293d595f711SHong Zhang       ncols = jmax-jmin;
294d595f711SHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
295d595f711SHong Zhang       ierr = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
296d595f711SHong Zhang       nzk += nlnk;
297d595f711SHong Zhang 
298d595f711SHong Zhang       /* update il and jl for prow */
299d595f711SHong Zhang       if (jmin < jmax){
300d595f711SHong Zhang         il[prow] = jmin;
301d595f711SHong Zhang         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
302d595f711SHong Zhang       }
303d595f711SHong Zhang       prow = nextprow;
304d595f711SHong Zhang     }
305d595f711SHong Zhang 
306d595f711SHong Zhang     /* if free space is not available, make more free space */
307d595f711SHong Zhang     if (current_space->local_remaining<nzk) {
308d595f711SHong Zhang       i  = mbs - k + 1; /* num of unfactored rows */
3097b056e98SHong Zhang       i *= PetscMin(nzk, i-1); /* i*nzk, i*(i-1): estimated and max additional space needed */
310d595f711SHong Zhang       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
311d595f711SHong Zhang       reallocs++;
312d595f711SHong Zhang     }
313d595f711SHong Zhang 
314d595f711SHong Zhang     /* copy data into free space, then initialize lnk */
315d595f711SHong Zhang     ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
316d595f711SHong Zhang 
317d595f711SHong Zhang     /* add the k-th row into il and jl */
3187b056e98SHong Zhang     if (nzk > 1){
319d595f711SHong Zhang       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
320d595f711SHong Zhang       jl[k] = jl[i]; jl[i] = k;
321d595f711SHong Zhang       il[k] = ui[k] + 1;
322d595f711SHong Zhang     }
323d595f711SHong Zhang     ui_ptr[k] = current_space->array;
324d595f711SHong Zhang     current_space->array           += nzk;
325d595f711SHong Zhang     current_space->local_used      += nzk;
326d595f711SHong Zhang     current_space->local_remaining -= nzk;
327d595f711SHong Zhang 
328d595f711SHong Zhang     ui[k+1] = ui[k] + nzk;
329d595f711SHong Zhang   }
330d595f711SHong Zhang 
331d595f711SHong Zhang #if defined(PETSC_USE_INFO)
332d595f711SHong Zhang   if (ai[mbs] != 0) {
333d595f711SHong Zhang     PetscReal af = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
334d595f711SHong Zhang     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
335d595f711SHong Zhang     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
336d595f711SHong Zhang     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
337d595f711SHong Zhang   } else {
338d595f711SHong Zhang     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
339d595f711SHong Zhang   }
340d595f711SHong Zhang #endif
341d595f711SHong Zhang 
342d595f711SHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
343d595f711SHong Zhang   ierr = PetscFree4(ui_ptr,il,jl,cols);CHKERRQ(ierr);
344d595f711SHong Zhang 
345d595f711SHong Zhang   /* destroy list of free space and other temporary array(s) */
346d595f711SHong Zhang   ierr = PetscMalloc((ui[mbs]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
3477b056e98SHong Zhang   ierr = PetscFreeSpaceContiguous_Cholesky(&free_space,uj,mbs,ui,udiag);CHKERRQ(ierr); /* store matrix factor */
348d595f711SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
349d595f711SHong Zhang 
350d595f711SHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
351d595f711SHong Zhang   ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(fact,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
352d595f711SHong Zhang 
3537b056e98SHong Zhang   b = (Mat_SeqSBAIJ*)fact->data;
354d595f711SHong Zhang   b->singlemalloc = PETSC_FALSE;
355d595f711SHong Zhang   b->free_a       = PETSC_TRUE;
356d595f711SHong Zhang   b->free_ij      = PETSC_TRUE;
357d595f711SHong Zhang   ierr = PetscMalloc((ui[mbs]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
358d595f711SHong Zhang   b->j    = uj;
359d595f711SHong Zhang   b->i    = ui;
360c6d0d4f0SHong Zhang   b->diag = udiag;
361c6d0d4f0SHong Zhang   b->free_diag = PETSC_TRUE;
362d595f711SHong Zhang   b->ilen = 0;
363d595f711SHong Zhang   b->imax = 0;
364d595f711SHong Zhang   b->row  = perm;
365d595f711SHong Zhang   b->icol = perm;
366d595f711SHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
3677b056e98SHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
3687b056e98SHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
369d595f711SHong Zhang   ierr    = PetscMalloc((mbs+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
3707b056e98SHong Zhang   ierr    = PetscLogObjectMemory(fact,ui[mbs]*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
371d595f711SHong Zhang   b->maxnz = b->nz = ui[mbs];
372d595f711SHong Zhang 
373d595f711SHong Zhang   (fact)->info.factor_mallocs    = reallocs;
374d595f711SHong Zhang   (fact)->info.fill_ratio_given  = fill;
375d595f711SHong Zhang   if (ai[mbs] != 0) {
376d595f711SHong Zhang     (fact)->info.fill_ratio_needed = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
377d595f711SHong Zhang   } else {
378d595f711SHong Zhang     (fact)->info.fill_ratio_needed = 0.0;
379d595f711SHong Zhang   }
380c6d0d4f0SHong Zhang   fact->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering;
381d595f711SHong Zhang   PetscFunctionReturn(0);
382d595f711SHong Zhang }
383d595f711SHong Zhang 
384d595f711SHong Zhang #undef __FUNCT__
385d595f711SHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqSBAIJ_inplace"
386d595f711SHong Zhang PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ_inplace(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
387d595f711SHong Zhang {
388d595f711SHong Zhang   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data;
389d595f711SHong Zhang   Mat_SeqSBAIJ       *b;
390d595f711SHong Zhang   PetscErrorCode     ierr;
391d595f711SHong Zhang   PetscTruth         perm_identity,missing;
392d595f711SHong Zhang   PetscReal          fill = info->fill;
393d595f711SHong Zhang   const PetscInt     *rip,*ai,*aj;
394d595f711SHong Zhang   PetscInt           i,mbs=a->mbs,bs=A->rmap->bs,reallocs=0,prow,d;
395d595f711SHong Zhang   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
396d595f711SHong Zhang   PetscInt           nlnk,*lnk,ncols,*cols,*uj,**ui_ptr,*uj_ptr;
397d595f711SHong Zhang   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
398d595f711SHong Zhang   PetscBT            lnkbt;
39949b5e25fSSatish Balay 
40049b5e25fSSatish Balay   PetscFunctionBegin;
40158ebbce7SBarry Smith   ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr);
40258ebbce7SBarry Smith   if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
40358ebbce7SBarry Smith 
40410c27e3dSHong Zhang   /*
40510c27e3dSHong Zhang    This code originally uses Modified Sparse Row (MSR) storage
40610c27e3dSHong Zhang    (see page 85, "Iterative Methods ..." by Saad) for the output matrix B - bad choise!
40710c27e3dSHong Zhang    Then it is rewritten so the factor B takes seqsbaij format. However the associated
40810c27e3dSHong Zhang    MatCholeskyFactorNumeric_() have not been modified for the cases of bs>1 or !perm_identity,
40910c27e3dSHong Zhang    thus the original code in MSR format is still used for these cases.
41010c27e3dSHong Zhang    The code below should replace MatCholeskyFactorSymbolic_SeqSBAIJ_MSR() whenever
41110c27e3dSHong Zhang    MatCholeskyFactorNumeric_() is modified for using sbaij symbolic factor.
41210c27e3dSHong Zhang   */
413fff829cfSHong Zhang   if (bs > 1){
414719d5645SBarry Smith     ierr = MatCholeskyFactorSymbolic_SeqSBAIJ_MSR(fact,A,perm,info);CHKERRQ(ierr);
415fff829cfSHong Zhang     PetscFunctionReturn(0);
416fff829cfSHong Zhang   }
41710c27e3dSHong Zhang 
41898a8e78dSHong Zhang   /* check whether perm is the identity mapping */
41998a8e78dSHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
42098a8e78dSHong Zhang 
421fff829cfSHong Zhang   if (perm_identity){
422fff829cfSHong Zhang     a->permute = PETSC_FALSE;
423fff829cfSHong Zhang     ai = a->i; aj = a->j;
424fff829cfSHong Zhang   } else {
425d33dd558SHong Zhang     SETERRQ(PETSC_ERR_SUP,"Matrix reordering is not supported for sbaij matrix. Use aij format");
426d33dd558SHong Zhang     /* There are bugs for reordeing. Needs further work.
427d33dd558SHong Zhang        MatReordering for sbaij cannot be efficient. User should use aij formt! */
428fff829cfSHong Zhang     a->permute = PETSC_TRUE;
429fff829cfSHong Zhang     ierr = MatReorderingSeqSBAIJ(A,perm);CHKERRQ(ierr);
430fff829cfSHong Zhang     ai = a->inew; aj = a->jnew;
431fff829cfSHong Zhang   }
432fff829cfSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
433fff829cfSHong Zhang 
434fff829cfSHong Zhang   /* initialization */
4357625dc9aSHong Zhang   ierr  = PetscMalloc((mbs+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
43698a8e78dSHong Zhang   ui[0] = 0;
43798a8e78dSHong Zhang 
43898a8e78dSHong Zhang   /* jl: linked list for storing indices of the pivot rows
4397625dc9aSHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
440d8c74875SBarry Smith   ierr = PetscMalloc4(mbs,PetscInt*,&ui_ptr,mbs,PetscInt,&il,mbs,PetscInt,&jl,mbs,PetscInt,&cols);CHKERRQ(ierr);
4417625dc9aSHong Zhang   for (i=0; i<mbs; i++){
4427625dc9aSHong Zhang     jl[i] = mbs; il[i] = 0;
443fff829cfSHong Zhang   }
444fff829cfSHong Zhang 
44598a8e78dSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
4467625dc9aSHong Zhang   nlnk = mbs + 1;
4477625dc9aSHong Zhang   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
448fff829cfSHong Zhang 
4497625dc9aSHong Zhang   /* initial FreeSpace size is fill*(ai[mbs]+1) */
450a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[mbs]+1)),&free_space);CHKERRQ(ierr);
45198a8e78dSHong Zhang   current_space = free_space;
45298a8e78dSHong Zhang 
4537625dc9aSHong Zhang   for (k=0; k<mbs; k++){  /* for each active row k */
45498a8e78dSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
45598a8e78dSHong Zhang     nzk   = 0;
45698a8e78dSHong Zhang     ncols = ai[rip[k]+1] - ai[rip[k]];
4577625dc9aSHong Zhang     for (j=0; j<ncols; j++){
4587625dc9aSHong Zhang       i = *(aj + ai[rip[k]] + j);
4597625dc9aSHong Zhang       cols[j] = rip[i];
4607625dc9aSHong Zhang     }
4617625dc9aSHong Zhang     ierr = PetscLLAdd(ncols,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
46298a8e78dSHong Zhang     nzk += nlnk;
46398a8e78dSHong Zhang 
46498a8e78dSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
46598a8e78dSHong Zhang     prow = jl[k]; /* 1st pivot row */
466fff829cfSHong Zhang 
467fff829cfSHong Zhang     while (prow < k){
468fff829cfSHong Zhang       nextprow = jl[prow];
46998a8e78dSHong Zhang       /* merge prow into k-th row */
4707625dc9aSHong Zhang       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
47198a8e78dSHong Zhang       jmax = ui[prow+1];
47298a8e78dSHong Zhang       ncols = jmax-jmin;
4737625dc9aSHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
4745a8e39fbSHong Zhang       ierr = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
47598a8e78dSHong Zhang       nzk += nlnk;
476fff829cfSHong Zhang 
47798a8e78dSHong Zhang       /* update il and jl for prow */
478fff829cfSHong Zhang       if (jmin < jmax){
479fff829cfSHong Zhang         il[prow] = jmin;
4807625dc9aSHong Zhang         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
481fff829cfSHong Zhang       }
482fff829cfSHong Zhang       prow = nextprow;
483fff829cfSHong Zhang     }
484fff829cfSHong Zhang 
48598a8e78dSHong Zhang     /* if free space is not available, make more free space */
48698a8e78dSHong Zhang     if (current_space->local_remaining<nzk) {
4877625dc9aSHong Zhang       i = mbs - k + 1; /* num of unfactored rows */
48898a8e78dSHong Zhang       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
489a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
49098a8e78dSHong Zhang       reallocs++;
49198a8e78dSHong Zhang     }
49298a8e78dSHong Zhang 
49398a8e78dSHong Zhang     /* copy data into free space, then initialize lnk */
4947625dc9aSHong Zhang     ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
49598a8e78dSHong Zhang 
49698a8e78dSHong Zhang     /* add the k-th row into il and jl */
49798a8e78dSHong Zhang     if (nzk-1 > 0){
4987625dc9aSHong Zhang       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
499fff829cfSHong Zhang       jl[k] = jl[i]; jl[i] = k;
50098a8e78dSHong Zhang       il[k] = ui[k] + 1;
501fff829cfSHong Zhang     }
5027625dc9aSHong Zhang     ui_ptr[k] = current_space->array;
50398a8e78dSHong Zhang     current_space->array           += nzk;
50498a8e78dSHong Zhang     current_space->local_used      += nzk;
50598a8e78dSHong Zhang     current_space->local_remaining -= nzk;
506fff829cfSHong Zhang 
50798a8e78dSHong Zhang     ui[k+1] = ui[k] + nzk;
508fff829cfSHong Zhang   }
509fff829cfSHong Zhang 
5106cf91177SBarry Smith #if defined(PETSC_USE_INFO)
5117625dc9aSHong Zhang   if (ai[mbs] != 0) {
5127625dc9aSHong Zhang     PetscReal af = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
513ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
514ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
515ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
516fff829cfSHong Zhang   } else {
517ae15b995SBarry Smith     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
518fff829cfSHong Zhang   }
51963ba0a88SBarry Smith #endif
520fff829cfSHong Zhang 
521fff829cfSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
522d8c74875SBarry Smith   ierr = PetscFree4(ui_ptr,il,jl,cols);CHKERRQ(ierr);
523fff829cfSHong Zhang 
52498a8e78dSHong Zhang   /* destroy list of free space and other temporary array(s) */
5257625dc9aSHong Zhang   ierr = PetscMalloc((ui[mbs]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
526a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
52798a8e78dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
528fff829cfSHong Zhang 
52998a8e78dSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
530719d5645SBarry Smith   ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(fact,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
53198a8e78dSHong Zhang 
532719d5645SBarry Smith   b = (Mat_SeqSBAIJ*)(fact)->data;
533fff829cfSHong Zhang   b->singlemalloc = PETSC_FALSE;
534e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
535e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
5367625dc9aSHong Zhang   ierr = PetscMalloc((ui[mbs]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
53798a8e78dSHong Zhang   b->j    = uj;
53898a8e78dSHong Zhang   b->i    = ui;
539fff829cfSHong Zhang   b->diag = 0;
540fff829cfSHong Zhang   b->ilen = 0;
541fff829cfSHong Zhang   b->imax = 0;
542fff829cfSHong Zhang   b->row  = perm;
543fff829cfSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
544fff829cfSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
545fff829cfSHong Zhang   b->icol = perm;
546fff829cfSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
5477625dc9aSHong Zhang   ierr    = PetscMalloc((mbs+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
548719d5645SBarry Smith   ierr    = PetscLogObjectMemory(fact,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
5497625dc9aSHong Zhang   b->maxnz = b->nz = ui[mbs];
550fff829cfSHong Zhang 
551719d5645SBarry Smith   (fact)->info.factor_mallocs    = reallocs;
552719d5645SBarry Smith   (fact)->info.fill_ratio_given  = fill;
5537625dc9aSHong Zhang   if (ai[mbs] != 0) {
554719d5645SBarry Smith     (fact)->info.fill_ratio_needed = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
555fff829cfSHong Zhang   } else {
556719d5645SBarry Smith     (fact)->info.fill_ratio_needed = 0.0;
557fff829cfSHong Zhang   }
558d595f711SHong Zhang   ierr = MatSeqSBAIJSetNumericFactorization_inplace(fact,perm_identity);CHKERRQ(ierr);
559064503c5SHong Zhang   PetscFunctionReturn(0);
560064503c5SHong Zhang }
561d595f711SHong Zhang 
5624a2ae208SSatish Balay #undef __FUNCT__
5634a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N"
5640481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat C,Mat A,const MatFactorInfo *info)
56549b5e25fSSatish Balay {
5664c16a6a6SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
5674c16a6a6SHong Zhang   IS             perm = b->row;
5686849ba73SBarry Smith   PetscErrorCode ierr;
5695d0c19d7SBarry Smith   const PetscInt *ai,*aj,*perm_ptr,mbs=a->mbs,*bi=b->i,*bj=b->j;
5705d0c19d7SBarry Smith   PetscInt       i,j;
5715d0c19d7SBarry Smith   PetscInt       *a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
572d0f46423SBarry Smith   PetscInt       bs=A->rmap->bs,bs2 = a->bs2,bslog = 0;
5734c16a6a6SHong Zhang   MatScalar      *ba = b->a,*aa,*ap,*dk,*uik;
5744c16a6a6SHong Zhang   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
57528de702eSHong Zhang   MatScalar      *work;
57613f74950SBarry Smith   PetscInt       *pivots;
5774c16a6a6SHong Zhang 
5784c16a6a6SHong Zhang   PetscFunctionBegin;
5794c16a6a6SHong Zhang   /* initialization */
58082502324SSatish Balay   ierr = PetscMalloc(bs2*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
5814c16a6a6SHong Zhang   ierr = PetscMemzero(rtmp,bs2*mbs*sizeof(MatScalar));CHKERRQ(ierr);
582d8c74875SBarry Smith   ierr = PetscMalloc2(mbs,PetscInt,&il,mbs,PetscInt,&jl);CHKERRQ(ierr);
5834c16a6a6SHong Zhang   for (i=0; i<mbs; i++) {
5844c16a6a6SHong Zhang     jl[i] = mbs; il[0] = 0;
5854c16a6a6SHong Zhang   }
586d8c74875SBarry Smith   ierr = PetscMalloc3(bs2,MatScalar,&dk,bs2,MatScalar,&uik,bs,MatScalar,&work);CHKERRQ(ierr);
58713f74950SBarry Smith   ierr = PetscMalloc(bs*sizeof(PetscInt),&pivots);CHKERRQ(ierr);
5884c16a6a6SHong Zhang 
5894c16a6a6SHong Zhang   ierr  = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr);
5904c16a6a6SHong Zhang 
5914c16a6a6SHong Zhang   /* check permutation */
5924c16a6a6SHong Zhang   if (!a->permute){
5934c16a6a6SHong Zhang     ai = a->i; aj = a->j; aa = a->a;
5944c16a6a6SHong Zhang   } else {
5954c16a6a6SHong Zhang     ai   = a->inew; aj = a->jnew;
59682502324SSatish Balay     ierr = PetscMalloc(bs2*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr);
5974c16a6a6SHong Zhang     ierr = PetscMemcpy(aa,a->a,bs2*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr);
59813f74950SBarry Smith     ierr = PetscMalloc(ai[mbs]*sizeof(PetscInt),&a2anew);CHKERRQ(ierr);
59913f74950SBarry Smith     ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(PetscInt));CHKERRQ(ierr);
6004c16a6a6SHong Zhang 
601187a9f4bSHong Zhang     /* flops in while loop */
602187a9f4bSHong Zhang     bslog = 2*bs*bs2;
603187a9f4bSHong Zhang 
6044c16a6a6SHong Zhang     for (i=0; i<mbs; i++){
6054c16a6a6SHong Zhang       jmin = ai[i]; jmax = ai[i+1];
6064c16a6a6SHong Zhang       for (j=jmin; j<jmax; j++){
6074c16a6a6SHong Zhang         while (a2anew[j] != j){
6084c16a6a6SHong Zhang           k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
6094c16a6a6SHong Zhang           for (k1=0; k1<bs2; k1++){
6104c16a6a6SHong Zhang             dk[k1]       = aa[k*bs2+k1];
6114c16a6a6SHong Zhang             aa[k*bs2+k1] = aa[j*bs2+k1];
6124c16a6a6SHong Zhang             aa[j*bs2+k1] = dk[k1];
6134c16a6a6SHong Zhang           }
6144c16a6a6SHong Zhang         }
6154c16a6a6SHong Zhang         /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
6164c16a6a6SHong Zhang         if (i > aj[j]){
6174c16a6a6SHong Zhang           /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */
6184c16a6a6SHong Zhang           ap = aa + j*bs2;                     /* ptr to the beginning of j-th block of aa */
6194c16a6a6SHong Zhang           for (k=0; k<bs2; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */
6204c16a6a6SHong Zhang           for (k=0; k<bs; k++){               /* j-th block of aa <- dk^T */
6214c16a6a6SHong Zhang             for (k1=0; k1<bs; k1++) *ap++ = dk[k + bs*k1];
6224c16a6a6SHong Zhang           }
6234c16a6a6SHong Zhang         }
6244c16a6a6SHong Zhang       }
6254c16a6a6SHong Zhang     }
626323b833fSBarry Smith     ierr = PetscFree(a2anew);CHKERRQ(ierr);
6274c16a6a6SHong Zhang   }
6284c16a6a6SHong Zhang 
6294c16a6a6SHong Zhang   /* for each row k */
6304c16a6a6SHong Zhang   for (k = 0; k<mbs; k++){
6314c16a6a6SHong Zhang 
6324c16a6a6SHong Zhang     /*initialize k-th row with elements nonzero in row perm(k) of A */
6334c16a6a6SHong Zhang     jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];
634057f5ba7SHong Zhang 
6354c16a6a6SHong Zhang     ap = aa + jmin*bs2;
6364c16a6a6SHong Zhang     for (j = jmin; j < jmax; j++){
6374c16a6a6SHong Zhang       vj = perm_ptr[aj[j]];         /* block col. index */
6384c16a6a6SHong Zhang       rtmp_ptr = rtmp + vj*bs2;
6394c16a6a6SHong Zhang       for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++;
6404c16a6a6SHong Zhang     }
6414c16a6a6SHong Zhang 
6424c16a6a6SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
6434c16a6a6SHong Zhang     ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
6444c16a6a6SHong Zhang     i = jl[k]; /* first row to be added to k_th row  */
6454c16a6a6SHong Zhang 
646057f5ba7SHong Zhang     while (i < k){
6474c16a6a6SHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
6484c16a6a6SHong Zhang 
6494c16a6a6SHong Zhang       /* compute multiplier */
6504c16a6a6SHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
6514c16a6a6SHong Zhang 
6524c16a6a6SHong Zhang       /* uik = -inv(Di)*U_bar(i,k) */
6534c16a6a6SHong Zhang       diag = ba + i*bs2;
6544c16a6a6SHong Zhang       u    = ba + ili*bs2;
6554c16a6a6SHong Zhang       ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
6564c16a6a6SHong Zhang       Kernel_A_gets_A_minus_B_times_C(bs,uik,diag,u);
6574c16a6a6SHong Zhang 
6584c16a6a6SHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
6594c16a6a6SHong Zhang       Kernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u);
660dc0b31edSSatish Balay       ierr = PetscLogFlops(bslog*2.0);CHKERRQ(ierr);
6614c16a6a6SHong Zhang 
6624c16a6a6SHong Zhang       /* update -U(i,k) */
6634c16a6a6SHong Zhang       ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
6644c16a6a6SHong Zhang 
6654c16a6a6SHong Zhang       /* add multiple of row i to k-th row ... */
6664c16a6a6SHong Zhang       jmin = ili + 1; jmax = bi[i+1];
6674c16a6a6SHong Zhang       if (jmin < jmax){
6684c16a6a6SHong Zhang         for (j=jmin; j<jmax; j++) {
6694c16a6a6SHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j) */
6704c16a6a6SHong Zhang           rtmp_ptr = rtmp + bj[j]*bs2;
6714c16a6a6SHong Zhang           u = ba + j*bs2;
6724c16a6a6SHong Zhang           Kernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u);
6734c16a6a6SHong Zhang         }
674187a9f4bSHong Zhang         ierr = PetscLogFlops(bslog*(jmax-jmin));CHKERRQ(ierr);
6754c16a6a6SHong Zhang 
6764c16a6a6SHong Zhang         /* ... add i to row list for next nonzero entry */
6774c16a6a6SHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
6784c16a6a6SHong Zhang         j     = bj[jmin];
6794c16a6a6SHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
6804c16a6a6SHong Zhang       }
6814c16a6a6SHong Zhang       i = nexti;
6824c16a6a6SHong Zhang     }
6834c16a6a6SHong Zhang 
6844c16a6a6SHong Zhang     /* save nonzero entries in k-th row of U ... */
6854c16a6a6SHong Zhang 
6864c16a6a6SHong Zhang     /* invert diagonal block */
6874c16a6a6SHong Zhang     diag = ba+k*bs2;
6884c16a6a6SHong Zhang     ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr);
689d230e6fdSBarry Smith     ierr = Kernel_A_gets_inverse_A(bs,diag,pivots,work);CHKERRQ(ierr);
6904c16a6a6SHong Zhang 
6914c16a6a6SHong Zhang     jmin = bi[k]; jmax = bi[k+1];
6924c16a6a6SHong Zhang     if (jmin < jmax) {
6934c16a6a6SHong Zhang       for (j=jmin; j<jmax; j++){
6944c16a6a6SHong Zhang          vj = bj[j];           /* block col. index of U */
6954c16a6a6SHong Zhang          u   = ba + j*bs2;
6964c16a6a6SHong Zhang          rtmp_ptr = rtmp + vj*bs2;
6974c16a6a6SHong Zhang          for (k1=0; k1<bs2; k1++){
6984c16a6a6SHong Zhang            *u++        = *rtmp_ptr;
6994c16a6a6SHong Zhang            *rtmp_ptr++ = 0.0;
7004c16a6a6SHong Zhang          }
7014c16a6a6SHong Zhang       }
7024c16a6a6SHong Zhang 
7034c16a6a6SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
7044c16a6a6SHong Zhang       il[k] = jmin;
7054c16a6a6SHong Zhang       i     = bj[jmin];
7064c16a6a6SHong Zhang       jl[k] = jl[i]; jl[i] = k;
7074c16a6a6SHong Zhang     }
7084c16a6a6SHong Zhang   }
7094c16a6a6SHong Zhang 
7104c16a6a6SHong Zhang   ierr = PetscFree(rtmp);CHKERRQ(ierr);
711d8c74875SBarry Smith   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
712d8c74875SBarry Smith   ierr = PetscFree3(dk,uik,work);CHKERRQ(ierr);
7134c16a6a6SHong Zhang   ierr = PetscFree(pivots);CHKERRQ(ierr);
7144c16a6a6SHong Zhang   if (a->permute){
7154c16a6a6SHong Zhang     ierr = PetscFree(aa);CHKERRQ(ierr);
7164c16a6a6SHong Zhang   }
7174c16a6a6SHong Zhang 
7184c16a6a6SHong Zhang   ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr);
7194f79d315SHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_N_inplace;
7204f79d315SHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_N_inplace;
7214f79d315SHong Zhang   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_N_inplace;
7224f79d315SHong Zhang   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_N_inplace;
723db4efbfdSBarry Smith 
7244c16a6a6SHong Zhang   C->assembled    = PETSC_TRUE;
7254c16a6a6SHong Zhang   C->preallocated = PETSC_TRUE;
726efee365bSSatish Balay   ierr = PetscLogFlops(1.3333*bs*bs2*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
7274c16a6a6SHong Zhang   PetscFunctionReturn(0);
7284c16a6a6SHong Zhang }
729d4edadd4SHong Zhang 
7304a2ae208SSatish Balay #undef __FUNCT__
7314a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering"
7320481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info)
733671cb588SHong Zhang {
734671cb588SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
735dfbe8321SBarry Smith   PetscErrorCode ierr;
73613f74950SBarry Smith   PetscInt       i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
73713f74950SBarry Smith   PetscInt       *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
738d0f46423SBarry Smith   PetscInt       bs=A->rmap->bs,bs2 = a->bs2,bslog;
739671cb588SHong Zhang   MatScalar      *ba = b->a,*aa,*ap,*dk,*uik;
740671cb588SHong Zhang   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
74128de702eSHong Zhang   MatScalar      *work;
74213f74950SBarry Smith   PetscInt       *pivots;
743671cb588SHong Zhang 
744671cb588SHong Zhang   PetscFunctionBegin;
74582502324SSatish Balay   ierr = PetscMalloc(bs2*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
746671cb588SHong Zhang   ierr = PetscMemzero(rtmp,bs2*mbs*sizeof(MatScalar));CHKERRQ(ierr);
747d8c74875SBarry Smith   ierr = PetscMalloc2(mbs,PetscInt,&il,mbs,PetscInt,&jl);CHKERRQ(ierr);
748671cb588SHong Zhang   for (i=0; i<mbs; i++) {
749671cb588SHong Zhang     jl[i] = mbs; il[0] = 0;
750671cb588SHong Zhang   }
751d8c74875SBarry Smith   ierr = PetscMalloc3(bs2,MatScalar,&dk,bs2,MatScalar,&uik,bs,MatScalar,&work);CHKERRQ(ierr);
75213f74950SBarry Smith   ierr = PetscMalloc(bs*sizeof(PetscInt),&pivots);CHKERRQ(ierr);
753671cb588SHong Zhang 
754671cb588SHong Zhang   ai = a->i; aj = a->j; aa = a->a;
755671cb588SHong Zhang 
756187a9f4bSHong Zhang   /* flops in while loop */
757187a9f4bSHong Zhang   bslog = 2*bs*bs2;
758187a9f4bSHong Zhang 
759671cb588SHong Zhang   /* for each row k */
760671cb588SHong Zhang   for (k = 0; k<mbs; k++){
761671cb588SHong Zhang 
762671cb588SHong Zhang     /*initialize k-th row with elements nonzero in row k of A */
763671cb588SHong Zhang     jmin = ai[k]; jmax = ai[k+1];
764671cb588SHong Zhang     ap = aa + jmin*bs2;
765671cb588SHong Zhang     for (j = jmin; j < jmax; j++){
766671cb588SHong Zhang       vj = aj[j];         /* block col. index */
767671cb588SHong Zhang       rtmp_ptr = rtmp + vj*bs2;
768671cb588SHong Zhang       for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++;
769671cb588SHong Zhang     }
770671cb588SHong Zhang 
771671cb588SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
772671cb588SHong Zhang     ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
773671cb588SHong Zhang     i = jl[k]; /* first row to be added to k_th row  */
774671cb588SHong Zhang 
775057f5ba7SHong Zhang     while (i < k){
776671cb588SHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
777671cb588SHong Zhang 
778671cb588SHong Zhang       /* compute multiplier */
779671cb588SHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
780671cb588SHong Zhang 
781671cb588SHong Zhang       /* uik = -inv(Di)*U_bar(i,k) */
782671cb588SHong Zhang       diag = ba + i*bs2;
783671cb588SHong Zhang       u    = ba + ili*bs2;
784671cb588SHong Zhang       ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
785671cb588SHong Zhang       Kernel_A_gets_A_minus_B_times_C(bs,uik,diag,u);
786671cb588SHong Zhang 
787671cb588SHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
788671cb588SHong Zhang       Kernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u);
789dc0b31edSSatish Balay       ierr = PetscLogFlops(bslog*2.0);CHKERRQ(ierr);
790671cb588SHong Zhang 
791671cb588SHong Zhang       /* update -U(i,k) */
792671cb588SHong Zhang       ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
793671cb588SHong Zhang 
794671cb588SHong Zhang       /* add multiple of row i to k-th row ... */
795671cb588SHong Zhang       jmin = ili + 1; jmax = bi[i+1];
796671cb588SHong Zhang       if (jmin < jmax){
797671cb588SHong Zhang         for (j=jmin; j<jmax; j++) {
798671cb588SHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j) */
799671cb588SHong Zhang           rtmp_ptr = rtmp + bj[j]*bs2;
800671cb588SHong Zhang           u = ba + j*bs2;
801671cb588SHong Zhang           Kernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u);
802671cb588SHong Zhang         }
803187a9f4bSHong Zhang         ierr = PetscLogFlops(bslog*(jmax-jmin));CHKERRQ(ierr);
804671cb588SHong Zhang 
805671cb588SHong Zhang         /* ... add i to row list for next nonzero entry */
806671cb588SHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
807671cb588SHong Zhang         j     = bj[jmin];
808671cb588SHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
809671cb588SHong Zhang       }
810671cb588SHong Zhang       i = nexti;
811671cb588SHong Zhang     }
812671cb588SHong Zhang 
813671cb588SHong Zhang     /* save nonzero entries in k-th row of U ... */
814671cb588SHong Zhang 
815671cb588SHong Zhang     /* invert diagonal block */
816671cb588SHong Zhang     diag = ba+k*bs2;
817671cb588SHong Zhang     ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr);
818d230e6fdSBarry Smith     ierr = Kernel_A_gets_inverse_A(bs,diag,pivots,work);CHKERRQ(ierr);
819671cb588SHong Zhang 
820671cb588SHong Zhang     jmin = bi[k]; jmax = bi[k+1];
821671cb588SHong Zhang     if (jmin < jmax) {
822671cb588SHong Zhang       for (j=jmin; j<jmax; j++){
823671cb588SHong Zhang          vj = bj[j];           /* block col. index of U */
824671cb588SHong Zhang          u   = ba + j*bs2;
825671cb588SHong Zhang          rtmp_ptr = rtmp + vj*bs2;
826671cb588SHong Zhang          for (k1=0; k1<bs2; k1++){
827671cb588SHong Zhang            *u++        = *rtmp_ptr;
828671cb588SHong Zhang            *rtmp_ptr++ = 0.0;
829671cb588SHong Zhang          }
830671cb588SHong Zhang       }
831671cb588SHong Zhang 
832671cb588SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
833671cb588SHong Zhang       il[k] = jmin;
834671cb588SHong Zhang       i     = bj[jmin];
835671cb588SHong Zhang       jl[k] = jl[i]; jl[i] = k;
836671cb588SHong Zhang     }
837671cb588SHong Zhang   }
838671cb588SHong Zhang 
839671cb588SHong Zhang   ierr = PetscFree(rtmp);CHKERRQ(ierr);
840d8c74875SBarry Smith   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
841d8c74875SBarry Smith   ierr = PetscFree3(dk,uik,work);CHKERRQ(ierr);
842671cb588SHong Zhang   ierr = PetscFree(pivots);CHKERRQ(ierr);
843671cb588SHong Zhang 
8444f79d315SHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_N_NaturalOrdering_inplace;
8454f79d315SHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_N_NaturalOrdering_inplace;
8464f79d315SHong Zhang   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_N_NaturalOrdering_inplace;
8474f79d315SHong Zhang   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_N_NaturalOrdering_inplace;
848671cb588SHong Zhang   C->assembled = PETSC_TRUE;
849671cb588SHong Zhang   C->preallocated = PETSC_TRUE;
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];
87062bba022SBarry Smith   PetscReal      shift = info->shiftinblocks;
87149b5e25fSSatish Balay 
87249b5e25fSSatish Balay   PetscFunctionBegin;
87391602c66SHong Zhang   /* initialization */
87491602c66SHong Zhang   /* il and jl record the first nonzero element in each row of the accessing
87591602c66SHong Zhang      window U(0:k, k:mbs-1).
87691602c66SHong Zhang      jl:    list of rows to be added to uneliminated rows
87791602c66SHong Zhang             i>= k: jl(i) is the first row to be added to row i
87891602c66SHong Zhang             i<  k: jl(i) is the row following row i in some list of rows
87991602c66SHong Zhang             jl(i) = mbs indicates the end of a list
88091602c66SHong Zhang      il(i): points to the first nonzero element in columns k,...,mbs-1 of
88191602c66SHong Zhang             row i of U */
88282502324SSatish Balay   ierr = PetscMalloc(4*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
883cc0c071aSHong Zhang   ierr = PetscMemzero(rtmp,4*mbs*sizeof(MatScalar));CHKERRQ(ierr);
884d8c74875SBarry Smith   ierr = PetscMalloc2(mbs,PetscInt,&il,mbs,PetscInt,&jl);CHKERRQ(ierr);
88591602c66SHong Zhang   for (i=0; i<mbs; i++) {
8863845f261SHong Zhang     jl[i] = mbs; il[0] = 0;
88791602c66SHong Zhang   }
888cc0c071aSHong Zhang   ierr = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr);
889cc0c071aSHong Zhang 
890cc0c071aSHong Zhang   /* check permutation */
891cc0c071aSHong Zhang   if (!a->permute){
892cc0c071aSHong Zhang     ai = a->i; aj = a->j; aa = a->a;
893cc0c071aSHong Zhang   } else {
894cc0c071aSHong Zhang     ai   = a->inew; aj = a->jnew;
89582502324SSatish Balay     ierr = PetscMalloc(4*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr);
896cc0c071aSHong Zhang     ierr = PetscMemcpy(aa,a->a,4*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr);
89713f74950SBarry Smith     ierr = PetscMalloc(ai[mbs]*sizeof(PetscInt),&a2anew);CHKERRQ(ierr);
89813f74950SBarry Smith     ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(PetscInt));CHKERRQ(ierr);
899cc0c071aSHong Zhang 
900cc0c071aSHong Zhang     for (i=0; i<mbs; i++){
901cc0c071aSHong Zhang       jmin = ai[i]; jmax = ai[i+1];
902cc0c071aSHong Zhang       for (j=jmin; j<jmax; j++){
903cc0c071aSHong Zhang         while (a2anew[j] != j){
904cc0c071aSHong Zhang           k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
905cc0c071aSHong Zhang           for (k1=0; k1<4; k1++){
906cc0c071aSHong Zhang             dk[k1]       = aa[k*4+k1];
907cc0c071aSHong Zhang             aa[k*4+k1] = aa[j*4+k1];
908cc0c071aSHong Zhang             aa[j*4+k1] = dk[k1];
909cc0c071aSHong Zhang           }
910cc0c071aSHong Zhang         }
911cc0c071aSHong Zhang         /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
912cc0c071aSHong Zhang         if (i > aj[j]){
913a1723e09SHong Zhang           /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */
914cc0c071aSHong Zhang           ap = aa + j*4;     /* ptr to the beginning of the block */
915cc0c071aSHong Zhang           dk[1] = ap[1];     /* swap ap[1] and ap[2] */
916cc0c071aSHong Zhang           ap[1] = ap[2];
917cc0c071aSHong Zhang           ap[2] = dk[1];
918cc0c071aSHong Zhang         }
919cc0c071aSHong Zhang       }
920cc0c071aSHong Zhang     }
921ac355199SBarry Smith     ierr = PetscFree(a2anew);CHKERRQ(ierr);
922cc0c071aSHong Zhang   }
9233845f261SHong Zhang 
92491602c66SHong Zhang   /* for each row k */
92591602c66SHong Zhang   for (k = 0; k<mbs; k++){
92691602c66SHong Zhang 
92791602c66SHong Zhang     /*initialize k-th row with elements nonzero in row perm(k) of A */
928cc0c071aSHong Zhang     jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];
929cc0c071aSHong Zhang     ap = aa + jmin*4;
93091602c66SHong Zhang     for (j = jmin; j < jmax; j++){
931cc0c071aSHong Zhang       vj = perm_ptr[aj[j]];         /* block col. index */
932cc0c071aSHong Zhang       rtmp_ptr = rtmp + vj*4;
933cc0c071aSHong Zhang       for (i=0; i<4; i++) *rtmp_ptr++ = *ap++;
93491602c66SHong Zhang     }
93591602c66SHong Zhang 
93691602c66SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
937cc0c071aSHong Zhang     ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr);
93891602c66SHong Zhang     i = jl[k]; /* first row to be added to k_th row  */
93991602c66SHong Zhang 
940057f5ba7SHong Zhang     while (i < k){
94191602c66SHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
94291602c66SHong Zhang 
9433845f261SHong Zhang       /* compute multiplier */
94491602c66SHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
9453845f261SHong Zhang 
9463845f261SHong Zhang       /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */
947cc0c071aSHong Zhang       diag = ba + i*4;
948cc0c071aSHong Zhang       u    = ba + ili*4;
949cc0c071aSHong Zhang       uik[0] = -(diag[0]*u[0] + diag[2]*u[1]);
950cc0c071aSHong Zhang       uik[1] = -(diag[1]*u[0] + diag[3]*u[1]);
951cc0c071aSHong Zhang       uik[2] = -(diag[0]*u[2] + diag[2]*u[3]);
952cc0c071aSHong Zhang       uik[3] = -(diag[1]*u[2] + diag[3]*u[3]);
9533845f261SHong Zhang 
9543845f261SHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */
955cc0c071aSHong Zhang       dk[0] += uik[0]*u[0] + uik[1]*u[1];
956cc0c071aSHong Zhang       dk[1] += uik[2]*u[0] + uik[3]*u[1];
957cc0c071aSHong Zhang       dk[2] += uik[0]*u[2] + uik[1]*u[3];
958cc0c071aSHong Zhang       dk[3] += uik[2]*u[2] + uik[3]*u[3];
9593845f261SHong Zhang 
960dc0b31edSSatish Balay       ierr = PetscLogFlops(16.0*2.0);CHKERRQ(ierr);
961187a9f4bSHong Zhang 
9623845f261SHong Zhang       /* update -U(i,k): ba[ili] = uik */
963cc0c071aSHong Zhang       ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr);
96491602c66SHong Zhang 
96591602c66SHong Zhang       /* add multiple of row i to k-th row ... */
96691602c66SHong Zhang       jmin = ili + 1; jmax = bi[i+1];
96791602c66SHong Zhang       if (jmin < jmax){
9683845f261SHong Zhang         for (j=jmin; j<jmax; j++) {
9693845f261SHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */
970cc0c071aSHong Zhang           rtmp_ptr = rtmp + bj[j]*4;
971cc0c071aSHong Zhang           u = ba + j*4;
972cc0c071aSHong Zhang           rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1];
973cc0c071aSHong Zhang           rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1];
974cc0c071aSHong Zhang           rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3];
975cc0c071aSHong Zhang           rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3];
9763845f261SHong Zhang         }
977dc0b31edSSatish Balay         ierr = PetscLogFlops(16.0*(jmax-jmin));CHKERRQ(ierr);
9783845f261SHong Zhang 
97991602c66SHong Zhang         /* ... add i to row list for next nonzero entry */
98091602c66SHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
98191602c66SHong Zhang         j     = bj[jmin];
98291602c66SHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
98391602c66SHong Zhang       }
984a1723e09SHong Zhang       i = nexti;
98591602c66SHong Zhang     }
986cc0c071aSHong Zhang 
98791602c66SHong Zhang     /* save nonzero entries in k-th row of U ... */
9883845f261SHong Zhang 
989cc0c071aSHong Zhang     /* invert diagonal block */
990cc0c071aSHong Zhang     diag = ba+k*4;
991cc0c071aSHong Zhang     ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr);
99262bba022SBarry Smith     ierr = Kernel_A_gets_inverse_A_2(diag,shift);CHKERRQ(ierr);
9933845f261SHong Zhang 
99491602c66SHong Zhang     jmin = bi[k]; jmax = bi[k+1];
99591602c66SHong Zhang     if (jmin < jmax) {
99691602c66SHong Zhang       for (j=jmin; j<jmax; j++){
997cc0c071aSHong Zhang          vj = bj[j];           /* block col. index of U */
998cc0c071aSHong Zhang          u   = ba + j*4;
999cc0c071aSHong Zhang          rtmp_ptr = rtmp + vj*4;
1000cc0c071aSHong Zhang          for (k1=0; k1<4; k1++){
1001cc0c071aSHong Zhang            *u++        = *rtmp_ptr;
1002cc0c071aSHong Zhang            *rtmp_ptr++ = 0.0;
1003cc0c071aSHong Zhang          }
1004cc0c071aSHong Zhang       }
10053845f261SHong Zhang 
100691602c66SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
100791602c66SHong Zhang       il[k] = jmin;
100891602c66SHong Zhang       i     = bj[jmin];
100991602c66SHong Zhang       jl[k] = jl[i]; jl[i] = k;
101091602c66SHong Zhang     }
101191602c66SHong Zhang   }
10123845f261SHong Zhang 
101349b5e25fSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
1014d8c74875SBarry Smith   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
101591602c66SHong Zhang   if (a->permute) {
101691602c66SHong Zhang     ierr = PetscFree(aa);CHKERRQ(ierr);
101791602c66SHong Zhang   }
1018cc0c071aSHong Zhang   ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr);
10194f79d315SHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_2_inplace;
10204f79d315SHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_2_inplace;
102149b5e25fSSatish Balay   C->assembled = PETSC_TRUE;
10225c0bcdfcSHong Zhang   C->preallocated = PETSC_TRUE;
1023efee365bSSatish Balay   ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
102449b5e25fSSatish Balay   PetscFunctionReturn(0);
102549b5e25fSSatish Balay }
102691602c66SHong Zhang 
102749b5e25fSSatish Balay /*
102849b5e25fSSatish Balay       Version for when blocks are 2 by 2 Using natural ordering
102949b5e25fSSatish Balay */
10304a2ae208SSatish Balay #undef __FUNCT__
10314a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering"
10320481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info)
103349b5e25fSSatish Balay {
1034ab27746eSHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
1035dfbe8321SBarry Smith   PetscErrorCode ierr;
103613f74950SBarry Smith   PetscInt       i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
103713f74950SBarry Smith   PetscInt       *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
1038d8c74875SBarry Smith   MatScalar      *ba = b->a,*aa,*ap,dk[8],uik[8];
1039ab27746eSHong Zhang   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
104062bba022SBarry Smith   PetscReal      shift = info->shiftinblocks;
104149b5e25fSSatish Balay 
104249b5e25fSSatish Balay   PetscFunctionBegin;
1043ab27746eSHong Zhang   /* initialization */
1044ab27746eSHong Zhang   /* il and jl record the first nonzero element in each row of the accessing
1045ab27746eSHong Zhang      window U(0:k, k:mbs-1).
1046ab27746eSHong Zhang      jl:    list of rows to be added to uneliminated rows
1047ab27746eSHong Zhang             i>= k: jl(i) is the first row to be added to row i
1048ab27746eSHong Zhang             i<  k: jl(i) is the row following row i in some list of rows
1049ab27746eSHong Zhang             jl(i) = mbs indicates the end of a list
1050ab27746eSHong Zhang      il(i): points to the first nonzero element in columns k,...,mbs-1 of
1051ab27746eSHong Zhang             row i of U */
105282502324SSatish Balay   ierr = PetscMalloc(4*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
1053ab27746eSHong Zhang   ierr = PetscMemzero(rtmp,4*mbs*sizeof(MatScalar));CHKERRQ(ierr);
1054d8c74875SBarry Smith   ierr = PetscMalloc2(mbs,PetscInt,&il,mbs,PetscInt,&jl);CHKERRQ(ierr);
1055ab27746eSHong Zhang   for (i=0; i<mbs; i++) {
1056ab27746eSHong Zhang     jl[i] = mbs; il[0] = 0;
105749b5e25fSSatish Balay   }
1058ab27746eSHong Zhang   ai = a->i; aj = a->j; aa = a->a;
1059ab27746eSHong Zhang 
1060ab27746eSHong Zhang   /* for each row k */
1061ab27746eSHong Zhang   for (k = 0; k<mbs; k++){
1062ab27746eSHong Zhang 
1063ab27746eSHong Zhang     /*initialize k-th row with elements nonzero in row k of A */
1064ab27746eSHong Zhang     jmin = ai[k]; jmax = ai[k+1];
1065ab27746eSHong Zhang     ap = aa + jmin*4;
1066ab27746eSHong Zhang     for (j = jmin; j < jmax; j++){
1067ab27746eSHong Zhang       vj = aj[j];         /* block col. index */
1068ab27746eSHong Zhang       rtmp_ptr = rtmp + vj*4;
1069ab27746eSHong Zhang       for (i=0; i<4; i++) *rtmp_ptr++ = *ap++;
107049b5e25fSSatish Balay     }
1071ab27746eSHong Zhang 
1072ab27746eSHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
1073ab27746eSHong Zhang     ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr);
1074ab27746eSHong Zhang     i = jl[k]; /* first row to be added to k_th row  */
1075ab27746eSHong Zhang 
1076057f5ba7SHong Zhang     while (i < k){
1077ab27746eSHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
1078ab27746eSHong Zhang 
1079ab27746eSHong Zhang       /* compute multiplier */
1080ab27746eSHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
1081ab27746eSHong Zhang 
1082ab27746eSHong Zhang       /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */
1083ab27746eSHong Zhang       diag = ba + i*4;
1084ab27746eSHong Zhang       u    = ba + ili*4;
1085ab27746eSHong Zhang       uik[0] = -(diag[0]*u[0] + diag[2]*u[1]);
1086ab27746eSHong Zhang       uik[1] = -(diag[1]*u[0] + diag[3]*u[1]);
1087ab27746eSHong Zhang       uik[2] = -(diag[0]*u[2] + diag[2]*u[3]);
1088ab27746eSHong Zhang       uik[3] = -(diag[1]*u[2] + diag[3]*u[3]);
1089ab27746eSHong Zhang 
1090ab27746eSHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */
1091ab27746eSHong Zhang       dk[0] += uik[0]*u[0] + uik[1]*u[1];
1092ab27746eSHong Zhang       dk[1] += uik[2]*u[0] + uik[3]*u[1];
1093ab27746eSHong Zhang       dk[2] += uik[0]*u[2] + uik[1]*u[3];
1094ab27746eSHong Zhang       dk[3] += uik[2]*u[2] + uik[3]*u[3];
1095ab27746eSHong Zhang 
1096dc0b31edSSatish Balay       ierr = PetscLogFlops(16.0*2.0);CHKERRQ(ierr);
1097187a9f4bSHong Zhang 
1098ab27746eSHong Zhang       /* update -U(i,k): ba[ili] = uik */
1099ab27746eSHong Zhang       ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr);
1100ab27746eSHong Zhang 
1101ab27746eSHong Zhang       /* add multiple of row i to k-th row ... */
1102ab27746eSHong Zhang       jmin = ili + 1; jmax = bi[i+1];
1103ab27746eSHong Zhang       if (jmin < jmax){
1104ab27746eSHong Zhang         for (j=jmin; j<jmax; j++) {
1105ab27746eSHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */
1106ab27746eSHong Zhang           rtmp_ptr = rtmp + bj[j]*4;
1107ab27746eSHong Zhang           u = ba + j*4;
1108ab27746eSHong Zhang           rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1];
1109ab27746eSHong Zhang           rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1];
1110ab27746eSHong Zhang           rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3];
1111ab27746eSHong Zhang           rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3];
111249b5e25fSSatish Balay         }
1113dc0b31edSSatish Balay         ierr = PetscLogFlops(16.0*(jmax-jmin));CHKERRQ(ierr);
1114ab27746eSHong Zhang 
1115ab27746eSHong Zhang         /* ... add i to row list for next nonzero entry */
1116ab27746eSHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
1117ab27746eSHong Zhang         j     = bj[jmin];
1118ab27746eSHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
111949b5e25fSSatish Balay       }
1120ab27746eSHong Zhang       i = nexti;
112149b5e25fSSatish Balay     }
1122ab27746eSHong Zhang 
1123ab27746eSHong Zhang     /* save nonzero entries in k-th row of U ... */
1124ab27746eSHong Zhang 
112549b5e25fSSatish Balay     /* invert diagonal block */
1126ab27746eSHong Zhang     diag = ba+k*4;
1127ab27746eSHong Zhang     ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr);
112862bba022SBarry Smith     ierr = Kernel_A_gets_inverse_A_2(diag,shift);CHKERRQ(ierr);
1129ab27746eSHong Zhang 
1130ab27746eSHong Zhang     jmin = bi[k]; jmax = bi[k+1];
1131ab27746eSHong Zhang     if (jmin < jmax) {
1132ab27746eSHong Zhang       for (j=jmin; j<jmax; j++){
1133ab27746eSHong Zhang          vj = bj[j];           /* block col. index of U */
1134ab27746eSHong Zhang          u   = ba + j*4;
1135ab27746eSHong Zhang          rtmp_ptr = rtmp + vj*4;
1136ab27746eSHong Zhang          for (k1=0; k1<4; k1++){
1137ab27746eSHong Zhang            *u++        = *rtmp_ptr;
1138ab27746eSHong Zhang            *rtmp_ptr++ = 0.0;
1139ab27746eSHong Zhang          }
1140ab27746eSHong Zhang       }
1141ab27746eSHong Zhang 
1142ab27746eSHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
1143ab27746eSHong Zhang       il[k] = jmin;
1144ab27746eSHong Zhang       i     = bj[jmin];
1145ab27746eSHong Zhang       jl[k] = jl[i]; jl[i] = k;
1146ab27746eSHong Zhang     }
114749b5e25fSSatish Balay   }
114849b5e25fSSatish Balay 
114949b5e25fSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
1150d8c74875SBarry Smith   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
1151ab27746eSHong Zhang 
11524f79d315SHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_2_NaturalOrdering_inplace;
11534f79d315SHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_2_NaturalOrdering_inplace;
11544f79d315SHong Zhang   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_2_NaturalOrdering_inplace;
11554f79d315SHong Zhang   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_2_NaturalOrdering_inplace;
115649b5e25fSSatish Balay   C->assembled = PETSC_TRUE;
11575c0bcdfcSHong Zhang   C->preallocated = PETSC_TRUE;
1158efee365bSSatish Balay   ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
115949b5e25fSSatish Balay   PetscFunctionReturn(0);
116049b5e25fSSatish Balay }
116149b5e25fSSatish Balay 
116249b5e25fSSatish Balay /*
116398a8e78dSHong Zhang     Numeric U^T*D*U factorization for SBAIJ format.
116491602c66SHong Zhang     Version for blocks are 1 by 1.
116549b5e25fSSatish Balay */
11664a2ae208SSatish Balay #undef __FUNCT__
1167d595f711SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1_inplace"
1168d595f711SHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_inplace(Mat C,Mat A,const MatFactorInfo *info)
116949b5e25fSSatish Balay {
117049b5e25fSSatish Balay   Mat_SeqSBAIJ   *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ *)C->data;
117149b5e25fSSatish Balay   IS             ip=b->row;
11726849ba73SBarry Smith   PetscErrorCode ierr;
11735d0c19d7SBarry Smith   const PetscInt *ai,*aj,*rip;
11745d0c19d7SBarry Smith   PetscInt       *a2anew,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j,*bcol;
1175997a0949SHong Zhang   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
1176997a0949SHong Zhang   MatScalar      *rtmp,*ba=b->a,*bval,*aa,dk,uikdi;
11773cea4cbeSHong Zhang   PetscReal      zeropivot,rs,shiftnz;
1178fbf22428SSatish Balay   PetscReal      shiftpd;
11793cea4cbeSHong Zhang   ChShift_Ctx    sctx;
11803cea4cbeSHong Zhang   PetscInt       newshift;
118149b5e25fSSatish Balay 
118249b5e25fSSatish Balay   PetscFunctionBegin;
11833cea4cbeSHong Zhang   /* initialization */
11843cea4cbeSHong Zhang   shiftnz   = info->shiftnz;
11853cea4cbeSHong Zhang   shiftpd   = info->shiftpd;
11863cea4cbeSHong Zhang   zeropivot = info->zeropivot;
11873cea4cbeSHong Zhang 
118849b5e25fSSatish Balay   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
1189cb718733SHong Zhang   if (!a->permute){
11902d007305SHong Zhang     ai = a->i; aj = a->j; aa = a->a;
11912d007305SHong Zhang   } else {
11922d007305SHong Zhang     ai = a->inew; aj = a->jnew;
1193fff829cfSHong Zhang     nz = ai[mbs];
1194fff829cfSHong Zhang     ierr = PetscMalloc(nz*sizeof(MatScalar),&aa);CHKERRQ(ierr);
1195fff829cfSHong Zhang     a2anew = a->a2anew;
1196997a0949SHong Zhang     bval   = a->a;
1197fff829cfSHong Zhang     for (j=0; j<nz; j++){
1198997a0949SHong Zhang       aa[a2anew[j]] = *(bval++);
11992d007305SHong Zhang     }
12002d007305SHong Zhang   }
12012d007305SHong Zhang 
120291602c66SHong Zhang   /* initialization */
120349b5e25fSSatish Balay   /* il and jl record the first nonzero element in each row of the accessing
120449b5e25fSSatish Balay      window U(0:k, k:mbs-1).
120549b5e25fSSatish Balay      jl:    list of rows to be added to uneliminated rows
120649b5e25fSSatish Balay             i>= k: jl(i) is the first row to be added to row i
120749b5e25fSSatish Balay             i<  k: jl(i) is the row following row i in some list of rows
120849b5e25fSSatish Balay             jl(i) = mbs indicates the end of a list
120949b5e25fSSatish Balay      il(i): points to the first nonzero element in columns k,...,mbs-1 of
121049b5e25fSSatish Balay             row i of U */
1211d8c74875SBarry Smith   ierr = PetscMalloc3(mbs,MatScalar,&rtmp,mbs,PetscInt,&il,mbs,PetscInt,&jl);CHKERRQ(ierr);
1212997a0949SHong Zhang 
12133cea4cbeSHong Zhang   sctx.shift_amount = 0;
12143cea4cbeSHong Zhang   sctx.nshift       = 0;
1215997a0949SHong Zhang   do {
12163cea4cbeSHong Zhang     sctx.chshift = PETSC_FALSE;
121749b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
121849b5e25fSSatish Balay       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
121949b5e25fSSatish Balay     }
1220997a0949SHong Zhang 
122149b5e25fSSatish Balay     for (k = 0; k<mbs; k++){
1222997a0949SHong Zhang       /*initialize k-th row by the perm[k]-th row of A */
122349b5e25fSSatish Balay       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
12247701de02SHong Zhang       bval = ba + bi[k];
122549b5e25fSSatish Balay       for (j = jmin; j < jmax; j++){
1226997a0949SHong Zhang         col = rip[aj[j]];
1227997a0949SHong Zhang         rtmp[col] = aa[j];
12287701de02SHong Zhang         *bval++  = 0.0; /* for in-place factorization */
122949b5e25fSSatish Balay       }
12303cea4cbeSHong Zhang 
12313cea4cbeSHong Zhang       /* shift the diagonal of the matrix */
12323cea4cbeSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
123349b5e25fSSatish Balay 
123491602c66SHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
123549b5e25fSSatish Balay       dk = rtmp[k];
123649b5e25fSSatish Balay       i = jl[k]; /* first row to be added to k_th row  */
123749b5e25fSSatish Balay 
1238057f5ba7SHong Zhang       while (i < k){
123949b5e25fSSatish Balay         nexti = jl[i]; /* next row to be added to k_th row */
124049b5e25fSSatish Balay 
1241fff829cfSHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
124249b5e25fSSatish Balay         ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
1243997a0949SHong Zhang         uikdi = - ba[ili]*ba[bi[i]];  /* diagonal(k) */
124449b5e25fSSatish Balay         dk += uikdi*ba[ili];
1245658e7b3eSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
124649b5e25fSSatish Balay 
1247997a0949SHong Zhang         /* add multiple of row i to k-th row */
124849b5e25fSSatish Balay         jmin = ili + 1; jmax = bi[i+1];
124949b5e25fSSatish Balay         if (jmin < jmax){
125049b5e25fSSatish Balay           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
1251dc0b31edSSatish Balay           ierr = PetscLogFlops(2.0*(jmax-jmin));CHKERRQ(ierr);
1252187a9f4bSHong Zhang 
1253fff829cfSHong Zhang           /* update il and jl for row i */
1254fff829cfSHong Zhang           il[i] = jmin;
1255fff829cfSHong Zhang           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
125649b5e25fSSatish Balay         }
1257ab27746eSHong Zhang         i = nexti;
125849b5e25fSSatish Balay       }
125949b5e25fSSatish Balay 
12603cea4cbeSHong Zhang       /* shift the diagonals when zero pivot is detected */
12613cea4cbeSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
12623cea4cbeSHong Zhang       rs   = 0.0;
1263997a0949SHong Zhang       jmin = bi[k]+1;
1264997a0949SHong Zhang       nz   = bi[k+1] - jmin;
1265997a0949SHong Zhang       if (nz){
1266997a0949SHong Zhang         bcol = bj + jmin;
1267997a0949SHong Zhang         while (nz--){
126889f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
126989f845c8SHong Zhang           bcol++;
1270997a0949SHong Zhang         }
1271997a0949SHong Zhang       }
12723cea4cbeSHong Zhang 
12733cea4cbeSHong Zhang       sctx.rs = rs;
12743cea4cbeSHong Zhang       sctx.pv = dk;
127545938b79SHong Zhang       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
127645938b79SHong Zhang       if (newshift == 1) break;    /* sctx.shift_amount is updated */
127749b5e25fSSatish Balay 
1278997a0949SHong Zhang       /* copy data into U(k,:) */
127998a8e78dSHong Zhang       ba[bi[k]] = 1.0/dk; /* U(k,k) */
1280fff829cfSHong Zhang       jmin = bi[k]+1; jmax = bi[k+1];
128149b5e25fSSatish Balay       if (jmin < jmax) {
128249b5e25fSSatish Balay         for (j=jmin; j<jmax; j++){
1283997a0949SHong Zhang           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
128449b5e25fSSatish Balay         }
128598a8e78dSHong Zhang         /* add the k-th row into il and jl */
128649b5e25fSSatish Balay         il[k] = jmin;
128798a8e78dSHong Zhang         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
128849b5e25fSSatish Balay       }
128949b5e25fSSatish Balay     }
12903cea4cbeSHong Zhang   } while (sctx.chshift);
1291d8c74875SBarry Smith   ierr = PetscFree3(rtmp,il,jl);CHKERRQ(ierr);
129298a8e78dSHong Zhang   if (a->permute){ierr = PetscFree(aa);CHKERRQ(ierr);}
129349b5e25fSSatish Balay 
129449b5e25fSSatish Balay   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
12950a3c351aSHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_1_inplace;
12964f79d315SHong Zhang   C->ops->solves         = MatSolves_SeqSBAIJ_1_inplace;
12970a3c351aSHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_1_inplace;
12980a3c351aSHong Zhang   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_1_inplace;
12990a3c351aSHong Zhang   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_1_inplace;
130049b5e25fSSatish Balay   C->assembled    = PETSC_TRUE;
13015c0bcdfcSHong Zhang   C->preallocated = PETSC_TRUE;
1302d0f46423SBarry Smith   ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr);
13033cea4cbeSHong Zhang   if (sctx.nshift){
13043cea4cbeSHong Zhang     if (shiftnz) {
13051e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
13063cea4cbeSHong Zhang     } else if (shiftpd) {
13071e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
1308997a0949SHong Zhang     }
1309997a0949SHong Zhang   }
131049b5e25fSSatish Balay   PetscFunctionReturn(0);
131149b5e25fSSatish Balay }
131249b5e25fSSatish Balay 
1313671cb588SHong Zhang /*
1314d595f711SHong Zhang   Version for when blocks are 1 by 1 Using natural ordering - modified from MatCholeskyFactorNumeric_SeqAIJ()
1315671cb588SHong Zhang */
13164a2ae208SSatish Balay #undef __FUNCT__
13174a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering"
1318d595f711SHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info)
1319d595f711SHong Zhang {
1320d595f711SHong Zhang   Mat_SeqSBAIJ   *a=(Mat_SeqSBAIJ*)A->data;
13217b056e98SHong Zhang   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)B->data;
1322d595f711SHong Zhang   PetscErrorCode ierr;
1323d595f711SHong Zhang   PetscInt       i,j,mbs=A->rmap->n,*bi=b->i,*bj=b->j,*bdiag=b->diag,*bjtmp;
1324d595f711SHong Zhang   PetscInt       *ai=a->i,*aj=a->j;
1325d595f711SHong Zhang   PetscInt       k,jmin,jmax,*c2r,*il,col,nexti,ili,nz;
1326d595f711SHong Zhang   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
1327d595f711SHong Zhang 
1328d595f711SHong Zhang   LUShift_Ctx    sctx;
1329d595f711SHong Zhang   PetscInt       newshift;
1330d595f711SHong Zhang   PetscReal      rs;
1331d595f711SHong Zhang   MatScalar      d,*v;
1332d595f711SHong Zhang 
1333d595f711SHong Zhang   PetscFunctionBegin;
1334d595f711SHong Zhang   /* MatPivotSetUp(): initialize shift context sctx */
1335d595f711SHong Zhang   ierr = PetscMemzero(&sctx,sizeof(LUShift_Ctx));CHKERRQ(ierr);
1336d595f711SHong Zhang 
1337d595f711SHong Zhang   /* if both shift schemes are chosen by user, only use info->shiftpd */
1338d595f711SHong Zhang   if (info->shiftpd) { /* set sctx.shift_top=max{rs} */
1339d595f711SHong Zhang     sctx.shift_top = info->zeropivot;
1340d595f711SHong Zhang     for (i=0; i<mbs; i++) {
1341d595f711SHong Zhang       /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
1342d595f711SHong Zhang       d  = (aa)[a->diag[i]];
1343d595f711SHong Zhang       rs = -PetscAbsScalar(d) - PetscRealPart(d);
1344d595f711SHong Zhang       v  = aa+ai[i];
1345d595f711SHong Zhang       nz = ai[i+1] - ai[i];
1346d595f711SHong Zhang       for (j=0; j<nz; j++)
1347d595f711SHong Zhang 	rs += PetscAbsScalar(v[j]);
1348d595f711SHong Zhang       if (rs>sctx.shift_top) sctx.shift_top = rs;
1349d595f711SHong Zhang     }
1350d595f711SHong Zhang     sctx.shift_top   *= 1.1;
1351d595f711SHong Zhang     sctx.nshift_max   = 5;
1352d595f711SHong Zhang     sctx.shift_lo     = 0.;
1353d595f711SHong Zhang     sctx.shift_hi     = 1.;
1354d595f711SHong Zhang   }
1355d595f711SHong Zhang 
1356d595f711SHong Zhang   /* allocate working arrays
1357d595f711SHong Zhang      c2r: linked list, keep track of pivot rows for a given column. c2r[col]: head of the list for a given col
1358d595f711SHong 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
1359d595f711SHong Zhang   */
1360d595f711SHong Zhang   ierr = PetscMalloc3(mbs,MatScalar,&rtmp,mbs,PetscInt,&il,mbs,PetscInt,&c2r);CHKERRQ(ierr);
1361d595f711SHong Zhang 
1362d595f711SHong Zhang   do {
1363d595f711SHong Zhang     sctx.lushift = PETSC_FALSE;
1364d595f711SHong Zhang 
1365d595f711SHong Zhang     for (i=0; i<mbs; i++)  c2r[i] = mbs;
1366d595f711SHong Zhang     il[0] = 0;
1367d595f711SHong Zhang 
1368d595f711SHong Zhang     for (k = 0; k<mbs; k++){
1369d595f711SHong Zhang       /* zero rtmp */
1370d595f711SHong Zhang       nz = bi[k+1] - bi[k];
1371d595f711SHong Zhang       bjtmp = bj + bi[k];
13727b056e98SHong Zhang       for (j=0; j<nz; j++) rtmp[bjtmp[j]] = 0.0;
1373d595f711SHong Zhang 
1374d595f711SHong Zhang       /* load in initial unfactored row */
1375d595f711SHong Zhang       bval = ba + bi[k];
1376d595f711SHong Zhang       jmin = ai[k]; jmax = ai[k+1];
1377d595f711SHong Zhang       for (j = jmin; j < jmax; j++){
1378d595f711SHong Zhang         col = aj[j];
1379d595f711SHong Zhang         rtmp[col] = aa[j];
1380d595f711SHong Zhang         *bval++   = 0.0; /* for in-place factorization */
1381d595f711SHong Zhang       }
1382d595f711SHong Zhang       /* shift the diagonal of the matrix: ZeropivotApply() */
1383d595f711SHong Zhang       rtmp[k] += sctx.shift_amount;  /* shift the diagonal of the matrix */
1384d595f711SHong Zhang 
1385d595f711SHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
1386d595f711SHong Zhang       dk = rtmp[k];
1387d595f711SHong Zhang       i  = c2r[k]; /* first row to be added to k_th row  */
1388d595f711SHong Zhang 
1389d595f711SHong Zhang       while (i < k){
1390d595f711SHong Zhang         nexti = c2r[i]; /* next row to be added to k_th row */
1391d595f711SHong Zhang 
1392d595f711SHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
1393d595f711SHong Zhang         ili   = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
1394d595f711SHong Zhang         uikdi = - ba[ili]*ba[bdiag[i]];  /* diagonal(k) */
1395d595f711SHong Zhang         dk   += uikdi*ba[ili]; /* update diag[k] */
1396d595f711SHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
1397d595f711SHong Zhang 
1398d595f711SHong Zhang         /* add multiple of row i to k-th row */
1399d595f711SHong Zhang         jmin = ili + 1; jmax = bi[i+1];
1400d595f711SHong Zhang         if (jmin < jmax){
1401d595f711SHong Zhang           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
1402d595f711SHong Zhang           /* update il and c2r for row i */
1403d595f711SHong Zhang           il[i] = jmin;
1404d595f711SHong Zhang           j = bj[jmin]; c2r[i] = c2r[j]; c2r[j] = i;
1405d595f711SHong Zhang         }
1406d595f711SHong Zhang         i = nexti;
1407d595f711SHong Zhang       }
1408d595f711SHong Zhang 
1409d595f711SHong Zhang       /* copy data into U(k,:) */
1410d595f711SHong Zhang       rs   = 0.0;
1411d595f711SHong Zhang       jmin = bi[k]; jmax = bi[k+1]-1;
1412d595f711SHong Zhang       if (jmin < jmax) {
1413d595f711SHong Zhang         for (j=jmin; j<jmax; j++){
1414d595f711SHong Zhang           col = bj[j]; ba[j] = rtmp[col]; rs += PetscAbsScalar(ba[j]);
1415d595f711SHong Zhang         }
1416d595f711SHong Zhang         /* add the k-th row into il and c2r */
1417d595f711SHong Zhang         il[k] = jmin;
1418d595f711SHong Zhang         i = bj[jmin]; c2r[k] = c2r[i]; c2r[i] = k;
1419d595f711SHong Zhang       }
1420d595f711SHong Zhang 
1421d595f711SHong Zhang       /* MatPivotCheck() */
1422d595f711SHong Zhang       sctx.rs  = rs;
1423d595f711SHong Zhang       sctx.pv  = dk;
1424d595f711SHong Zhang       if (info->shiftnz){
1425d595f711SHong Zhang         ierr = MatPivotCheck_nz(info,sctx,k,newshift);CHKERRQ(ierr);
1426d595f711SHong Zhang       } else if (info->shiftpd){
1427d595f711SHong Zhang         ierr = MatPivotCheck_pd(info,sctx,k,newshift);CHKERRQ(ierr);
1428d595f711SHong Zhang       } else if (info->shiftinblocks){
1429d595f711SHong Zhang         ierr = MatPivotCheck_inblocks(info,sctx,k,newshift);CHKERRQ(ierr);
1430d595f711SHong Zhang       } else {
1431d595f711SHong Zhang         ierr = MatPivotCheck_none(info,sctx,k,newshift);CHKERRQ(ierr);
1432d595f711SHong Zhang       }
1433d595f711SHong Zhang       dk = sctx.pv;
1434d595f711SHong Zhang       if (newshift == 1) break;
1435d595f711SHong Zhang 
1436d595f711SHong Zhang       ba[bdiag[k]] = 1.0/dk; /* U(k,k) */
1437d595f711SHong Zhang     }
1438d595f711SHong Zhang   } while (sctx.lushift);
1439d595f711SHong Zhang 
1440d595f711SHong Zhang   ierr = PetscFree3(rtmp,il,c2r);CHKERRQ(ierr);
1441d595f711SHong Zhang 
1442d595f711SHong Zhang   B->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1443*9dcdb97aSHong Zhang   B->ops->solves          = MatSolves_SeqSBAIJ_1;
1444d595f711SHong Zhang   B->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1445d595f711SHong Zhang   B->ops->forwardsolve    = 0;
1446d595f711SHong Zhang   B->ops->backwardsolve   = 0;
1447d595f711SHong Zhang 
14487b056e98SHong Zhang   B->assembled    = PETSC_TRUE;
14497b056e98SHong Zhang   B->preallocated = PETSC_TRUE;
14507b056e98SHong Zhang   ierr = PetscLogFlops(B->rmap->n);CHKERRQ(ierr);
1451d595f711SHong Zhang 
1452d595f711SHong Zhang   /* MatPivotView() */
1453d595f711SHong Zhang   if (sctx.nshift){
1454d595f711SHong Zhang     if (info->shiftpd) {
1455d595f711SHong Zhang       ierr = PetscInfo4(A,"number of shift_pd tries %D, shift_amount %G, diagonal shifted up by %e fraction top_value %e\n",sctx.nshift,sctx.shift_amount,sctx.shift_fraction,sctx.shift_top);CHKERRQ(ierr);
1456d595f711SHong Zhang     } else if (info->shiftnz) {
1457d595f711SHong Zhang       ierr = PetscInfo2(A,"number of shift_nz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
1458d595f711SHong Zhang     } else if (info->shiftinblocks){
1459d595f711SHong Zhang       ierr = PetscInfo2(A,"number of shift_inblocks applied %D, each shift_amount %G\n",sctx.nshift,info->shiftinblocks);CHKERRQ(ierr);
1460d595f711SHong Zhang     }
1461d595f711SHong Zhang   }
1462c6d0d4f0SHong Zhang   /*
1463c6d0d4f0SHong Zhang   for (k = 0; k<mbs; k++){
1464c6d0d4f0SHong Zhang     printf(" row %d, diag %d ",k,bdiag[k]);
1465c6d0d4f0SHong Zhang     nz = bi[k+1] - bi[k];
1466c6d0d4f0SHong Zhang     bjtmp = bj + bi[k];
1467c6d0d4f0SHong Zhang     bval  = ba + bi[k];
1468c6d0d4f0SHong Zhang       for (j=0; j<nz; j++) {
1469c6d0d4f0SHong Zhang         printf(" (%d, %g),",bjtmp[j],bval[j]);
1470c6d0d4f0SHong Zhang       }
1471c6d0d4f0SHong Zhang       printf(" \n");
1472c6d0d4f0SHong Zhang   }
1473c6d0d4f0SHong Zhang   */
1474d595f711SHong Zhang   PetscFunctionReturn(0);
1475d595f711SHong Zhang }
1476d595f711SHong Zhang 
1477d595f711SHong Zhang #undef __FUNCT__
1478d595f711SHong Zhang #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering_inplace"
1479d595f711SHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering_inplace(Mat C,Mat A,const MatFactorInfo *info)
1480671cb588SHong Zhang {
1481671cb588SHong Zhang   Mat_SeqSBAIJ   *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ *)C->data;
1482dfbe8321SBarry Smith   PetscErrorCode ierr;
148313f74950SBarry Smith   PetscInt       i,j,mbs = a->mbs;
148413f74950SBarry Smith   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
14853cea4cbeSHong Zhang   PetscInt       k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz;
1486653a6975SHong Zhang   MatScalar      *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval;
14873cea4cbeSHong Zhang   PetscReal      zeropivot,rs,shiftnz;
1488fbf22428SSatish Balay   PetscReal      shiftpd;
14893cea4cbeSHong Zhang   ChShift_Ctx    sctx;
14903cea4cbeSHong Zhang   PetscInt       newshift;
1491653a6975SHong Zhang 
1492653a6975SHong Zhang   PetscFunctionBegin;
1493653a6975SHong Zhang   /* initialization */
14943cea4cbeSHong Zhang   shiftnz   = info->shiftnz;
14953cea4cbeSHong Zhang   shiftpd   = info->shiftpd;
14963cea4cbeSHong Zhang   zeropivot = info->zeropivot;
14973cea4cbeSHong Zhang 
1498653a6975SHong Zhang   /* il and jl record the first nonzero element in each row of the accessing
1499653a6975SHong Zhang      window U(0:k, k:mbs-1).
1500653a6975SHong Zhang      jl:    list of rows to be added to uneliminated rows
1501653a6975SHong Zhang             i>= k: jl(i) is the first row to be added to row i
1502653a6975SHong Zhang             i<  k: jl(i) is the row following row i in some list of rows
1503653a6975SHong Zhang             jl(i) = mbs indicates the end of a list
1504653a6975SHong Zhang      il(i): points to the first nonzero element in U(i,k:mbs-1)
1505653a6975SHong Zhang   */
1506653a6975SHong Zhang   ierr = PetscMalloc(mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
1507d8c74875SBarry Smith   ierr = PetscMalloc2(mbs,PetscInt,&il,mbs,PetscInt,&jl);CHKERRQ(ierr);
1508b00f7748SHong Zhang 
15093cea4cbeSHong Zhang   sctx.shift_amount = 0;
15103cea4cbeSHong Zhang   sctx.nshift       = 0;
1511b00f7748SHong Zhang   do {
15123cea4cbeSHong Zhang     sctx.chshift = PETSC_FALSE;
1513653a6975SHong Zhang     for (i=0; i<mbs; i++) {
1514653a6975SHong Zhang       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
1515653a6975SHong Zhang     }
1516653a6975SHong Zhang 
1517997a0949SHong Zhang     for (k = 0; k<mbs; k++){
1518653a6975SHong Zhang       /*initialize k-th row with elements nonzero in row perm(k) of A */
1519653a6975SHong Zhang       nz   = ai[k+1] - ai[k];
1520653a6975SHong Zhang       acol = aj + ai[k];
1521653a6975SHong Zhang       aval = aa + ai[k];
1522653a6975SHong Zhang       bval = ba + bi[k];
1523653a6975SHong Zhang       while (nz -- ){
1524653a6975SHong Zhang         rtmp[*acol++] = *aval++;
1525653a6975SHong Zhang         *bval++       = 0.0; /* for in-place factorization */
1526653a6975SHong Zhang       }
15273cea4cbeSHong Zhang 
15283cea4cbeSHong Zhang       /* shift the diagonal of the matrix */
15293cea4cbeSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
1530653a6975SHong Zhang 
1531653a6975SHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
1532653a6975SHong Zhang       dk = rtmp[k];
1533653a6975SHong Zhang       i  = jl[k]; /* first row to be added to k_th row  */
1534653a6975SHong Zhang 
1535653a6975SHong Zhang       while (i < k){
1536653a6975SHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
1537653a6975SHong Zhang         /* compute multiplier, update D(k) and U(i,k) */
1538653a6975SHong Zhang         ili   = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
1539653a6975SHong Zhang         uikdi = - ba[ili]*ba[bi[i]];
1540653a6975SHong Zhang         dk   += uikdi*ba[ili];
1541653a6975SHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
1542653a6975SHong Zhang 
1543653a6975SHong Zhang         /* add multiple of row i to k-th row ... */
1544653a6975SHong Zhang         jmin = ili + 1;
1545653a6975SHong Zhang         nz   = bi[i+1] - jmin;
1546653a6975SHong Zhang         if (nz > 0){
1547653a6975SHong Zhang           bcol = bj + jmin;
1548653a6975SHong Zhang           bval = ba + jmin;
1549dc0b31edSSatish Balay           ierr = PetscLogFlops(2.0*nz);CHKERRQ(ierr);
1550653a6975SHong Zhang           while (nz --) rtmp[*bcol++] += uikdi*(*bval++);
1551187a9f4bSHong Zhang 
1552997a0949SHong Zhang           /* update il and jl for i-th row */
1553997a0949SHong Zhang           il[i] = jmin;
1554997a0949SHong Zhang           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
1555653a6975SHong Zhang         }
1556653a6975SHong Zhang         i = nexti;
1557653a6975SHong Zhang       }
1558653a6975SHong Zhang 
15593cea4cbeSHong Zhang       /* shift the diagonals when zero pivot is detected */
15603cea4cbeSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
15613cea4cbeSHong Zhang       rs   = 0.0;
156221c26570Svictorle       jmin = bi[k]+1;
156321c26570Svictorle       nz   = bi[k+1] - jmin;
156421c26570Svictorle       if (nz){
156521c26570Svictorle         bcol = bj + jmin;
156621c26570Svictorle         while (nz--){
156789f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
156889f845c8SHong Zhang           bcol++;
156921c26570Svictorle         }
157021c26570Svictorle       }
15713cea4cbeSHong Zhang 
15723cea4cbeSHong Zhang       sctx.rs = rs;
15733cea4cbeSHong Zhang       sctx.pv = dk;
157445938b79SHong Zhang       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
157545938b79SHong Zhang       if (newshift == 1) break;    /* sctx.shift_amount is updated */
1576653a6975SHong Zhang 
1577997a0949SHong Zhang       /* copy data into U(k,:) */
1578653a6975SHong Zhang       ba[bi[k]] = 1.0/dk;
1579653a6975SHong Zhang       jmin      = bi[k]+1;
1580653a6975SHong Zhang       nz        = bi[k+1] - jmin;
1581653a6975SHong Zhang       if (nz){
1582653a6975SHong Zhang         bcol = bj + jmin;
1583653a6975SHong Zhang         bval = ba + jmin;
1584653a6975SHong Zhang         while (nz--){
1585653a6975SHong Zhang           *bval++       = rtmp[*bcol];
1586653a6975SHong Zhang           rtmp[*bcol++] = 0.0;
1587653a6975SHong Zhang         }
1588997a0949SHong Zhang         /* add k-th row into il and jl */
1589653a6975SHong Zhang         il[k] = jmin;
1590997a0949SHong Zhang         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
1591653a6975SHong Zhang       }
1592b00f7748SHong Zhang     } /* end of for (k = 0; k<mbs; k++) */
15933cea4cbeSHong Zhang   } while (sctx.chshift);
1594653a6975SHong Zhang   ierr = PetscFree(rtmp);CHKERRQ(ierr);
1595d8c74875SBarry Smith   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
1596653a6975SHong Zhang 
15970a3c351aSHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
15984f79d315SHong Zhang   C->ops->solves         = MatSolves_SeqSBAIJ_1_inplace;
15990a3c351aSHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
16000a3c351aSHong Zhang   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
16010a3c351aSHong Zhang   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
1602db4efbfdSBarry Smith 
1603653a6975SHong Zhang   C->assembled    = PETSC_TRUE;
1604653a6975SHong Zhang   C->preallocated = PETSC_TRUE;
1605d0f46423SBarry Smith   ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr);
16063cea4cbeSHong Zhang   if (sctx.nshift){
16073cea4cbeSHong Zhang     if (shiftnz) {
16081e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
16093cea4cbeSHong Zhang     } else if (shiftpd) {
16101e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
1611b00f7748SHong Zhang     }
1612fb10cecfSBarry Smith   }
1613653a6975SHong Zhang   PetscFunctionReturn(0);
1614653a6975SHong Zhang }
1615653a6975SHong Zhang 
1616653a6975SHong Zhang #undef __FUNCT__
16174a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactor_SeqSBAIJ"
16180481f469SBarry Smith PetscErrorCode MatCholeskyFactor_SeqSBAIJ(Mat A,IS perm,const MatFactorInfo *info)
161949b5e25fSSatish Balay {
1620dfbe8321SBarry Smith   PetscErrorCode ierr;
162149b5e25fSSatish Balay   Mat            C;
162249b5e25fSSatish Balay 
162349b5e25fSSatish Balay   PetscFunctionBegin;
1624719d5645SBarry Smith   ierr = MatGetFactor(A,"petsc",MAT_FACTOR_CHOLESKY,&C);CHKERRQ(ierr);
1625719d5645SBarry Smith   ierr = MatCholeskyFactorSymbolic(C,A,perm,info);CHKERRQ(ierr);
1626719d5645SBarry Smith   ierr = MatCholeskyFactorNumeric(C,A,info);CHKERRQ(ierr);
1627db4efbfdSBarry Smith   A->ops->solve            = C->ops->solve;
1628db4efbfdSBarry Smith   A->ops->solvetranspose   = C->ops->solvetranspose;
16294445ddedSHong Zhang   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
163049b5e25fSSatish Balay   PetscFunctionReturn(0);
163149b5e25fSSatish Balay }
163249b5e25fSSatish Balay 
163349b5e25fSSatish Balay 
1634