xref: /petsc/src/mat/impls/sbaij/seq/sbaijfact.c (revision 0a3c351a51b8afeba0e28ebe856db6cb3362db0e)
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 
43db4efbfdSBarry Smith extern PetscErrorCode MatSeqSBAIJSetNumericFactorization(Mat,PetscTruth);
44db4efbfdSBarry Smith 
455f9f512dSHong Zhang /*
465f9f512dSHong Zhang   Symbolic U^T*D*U factorization for SBAIJ format. Modified from SSF of YSMP.
4710c27e3dSHong Zhang   Use Modified Sparse Row (MSR) storage for u and ju. See page 85, "Iterative Methods ..." by Saad.
485f9f512dSHong Zhang */
4910c27e3dSHong Zhang #undef __FUNCT__
5010c27e3dSHong Zhang #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqSBAIJ_MSR"
510481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ_MSR(Mat F,Mat A,IS perm,const MatFactorInfo *info)
5210c27e3dSHong Zhang {
5310c27e3dSHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b;
5410c27e3dSHong Zhang   PetscErrorCode ierr;
555d0c19d7SBarry Smith   const PetscInt *rip,*ai,*aj;
565d0c19d7SBarry Smith   PetscInt       i,mbs = a->mbs,*jutmp,bs = A->rmap->bs,bs2=a->bs2;
5710c27e3dSHong Zhang   PetscInt       m,reallocs = 0,prow;
5810c27e3dSHong Zhang   PetscInt       *jl,*q,jmin,jmax,juidx,nzk,qm,*iu,*ju,k,j,vj,umax,maxadd;
5910c27e3dSHong Zhang   PetscReal      f = info->fill;
6010c27e3dSHong Zhang   PetscTruth     perm_identity;
6110c27e3dSHong Zhang 
6210c27e3dSHong Zhang   PetscFunctionBegin;
6310c27e3dSHong Zhang   /* check whether perm is the identity mapping */
6410c27e3dSHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
6510c27e3dSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
6610c27e3dSHong Zhang 
6710c27e3dSHong Zhang   if (perm_identity){ /* without permutation */
6810c27e3dSHong Zhang     a->permute = PETSC_FALSE;
6910c27e3dSHong Zhang     ai = a->i; aj = a->j;
7010c27e3dSHong Zhang   } else {            /* non-trivial permutation */
7110c27e3dSHong Zhang     a->permute = PETSC_TRUE;
7210c27e3dSHong Zhang     ierr = MatReorderingSeqSBAIJ(A,perm);CHKERRQ(ierr);
7310c27e3dSHong Zhang     ai = a->inew; aj = a->jnew;
7410c27e3dSHong Zhang   }
7510c27e3dSHong Zhang 
7610c27e3dSHong Zhang   /* initialization */
7710c27e3dSHong Zhang   ierr  = PetscMalloc((mbs+1)*sizeof(PetscInt),&iu);CHKERRQ(ierr);
7810c27e3dSHong Zhang   umax  = (PetscInt)(f*ai[mbs] + 1); umax += mbs + 1;
7910c27e3dSHong Zhang   ierr  = PetscMalloc(umax*sizeof(PetscInt),&ju);CHKERRQ(ierr);
8010c27e3dSHong Zhang   iu[0] = mbs+1;
8110c27e3dSHong Zhang   juidx = mbs + 1; /* index for ju */
82d8c74875SBarry Smith   /* jl linked list for pivot row -- linked list for col index */
83d8c74875SBarry Smith   ierr  = PetscMalloc2(mbs,PetscInt,&jl,mbs,PetscInt,&q);CHKERRQ(ierr);
8410c27e3dSHong Zhang   for (i=0; i<mbs; i++){
8510c27e3dSHong Zhang     jl[i] = mbs;
8610c27e3dSHong Zhang     q[i] = 0;
8710c27e3dSHong Zhang   }
8810c27e3dSHong Zhang 
8910c27e3dSHong Zhang   /* for each row k */
9010c27e3dSHong Zhang   for (k=0; k<mbs; k++){
9110c27e3dSHong Zhang     for (i=0; i<mbs; i++) q[i] = 0;  /* to be removed! */
9210c27e3dSHong Zhang     nzk  = 0; /* num. of nz blocks in k-th block row with diagonal block excluded */
9310c27e3dSHong Zhang     q[k] = mbs;
9410c27e3dSHong Zhang     /* initialize nonzero structure of k-th row to row rip[k] of A */
9510c27e3dSHong Zhang     jmin = ai[rip[k]] +1; /* exclude diag[k] */
9610c27e3dSHong Zhang     jmax = ai[rip[k]+1];
9710c27e3dSHong Zhang     for (j=jmin; j<jmax; j++){
9810c27e3dSHong Zhang       vj = rip[aj[j]]; /* col. value */
9910c27e3dSHong Zhang       if(vj > k){
10010c27e3dSHong Zhang         qm = k;
10110c27e3dSHong Zhang         do {
10210c27e3dSHong Zhang           m  = qm; qm = q[m];
10310c27e3dSHong Zhang         } while(qm < vj);
10410c27e3dSHong Zhang         if (qm == vj) {
10510c27e3dSHong Zhang           SETERRQ(PETSC_ERR_PLIB,"Duplicate entry in A\n");
10610c27e3dSHong Zhang         }
10710c27e3dSHong Zhang         nzk++;
10810c27e3dSHong Zhang         q[m]  = vj;
10910c27e3dSHong Zhang         q[vj] = qm;
11010c27e3dSHong Zhang       } /* if(vj > k) */
11110c27e3dSHong Zhang     } /* for (j=jmin; j<jmax; j++) */
11210c27e3dSHong Zhang 
11310c27e3dSHong Zhang     /* modify nonzero structure of k-th row by computing fill-in
11410c27e3dSHong Zhang        for each row i to be merged in */
11510c27e3dSHong Zhang     prow = k;
11610c27e3dSHong Zhang     prow = jl[prow]; /* next pivot row (== mbs for symbolic factorization) */
11710c27e3dSHong Zhang 
11810c27e3dSHong Zhang     while (prow < k){
11910c27e3dSHong Zhang       /* merge row prow into k-th row */
12010c27e3dSHong Zhang       jmin = iu[prow] + 1; jmax = iu[prow+1];
12110c27e3dSHong Zhang       qm = k;
12210c27e3dSHong Zhang       for (j=jmin; j<jmax; j++){
12310c27e3dSHong Zhang         vj = ju[j];
12410c27e3dSHong Zhang         do {
12510c27e3dSHong Zhang           m = qm; qm = q[m];
12610c27e3dSHong Zhang         } while (qm < vj);
12710c27e3dSHong Zhang         if (qm != vj){
12810c27e3dSHong Zhang          nzk++; q[m] = vj; q[vj] = qm; qm = vj;
12910c27e3dSHong Zhang         }
13010c27e3dSHong Zhang       }
13110c27e3dSHong Zhang       prow = jl[prow]; /* next pivot row */
13210c27e3dSHong Zhang     }
13310c27e3dSHong Zhang 
13410c27e3dSHong Zhang     /* add k to row list for first nonzero element in k-th row */
13510c27e3dSHong Zhang     if (nzk > 0){
13610c27e3dSHong Zhang       i = q[k]; /* col value of first nonzero element in U(k, k+1:mbs-1) */
13710c27e3dSHong Zhang       jl[k] = jl[i]; jl[i] = k;
13810c27e3dSHong Zhang     }
13910c27e3dSHong Zhang     iu[k+1] = iu[k] + nzk;
14010c27e3dSHong Zhang 
14110c27e3dSHong Zhang     /* allocate more space to ju if needed */
14210c27e3dSHong Zhang     if (iu[k+1] > umax) {
14310c27e3dSHong Zhang       /* estimate how much additional space we will need */
14410c27e3dSHong Zhang       /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */
14510c27e3dSHong Zhang       /* just double the memory each time */
14610c27e3dSHong Zhang       maxadd = umax;
14710c27e3dSHong Zhang       if (maxadd < nzk) maxadd = (mbs-k)*(nzk+1)/2;
14810c27e3dSHong Zhang       umax += maxadd;
14910c27e3dSHong Zhang 
15010c27e3dSHong Zhang       /* allocate a longer ju */
15110c27e3dSHong Zhang       ierr = PetscMalloc(umax*sizeof(PetscInt),&jutmp);CHKERRQ(ierr);
15210c27e3dSHong Zhang       ierr = PetscMemcpy(jutmp,ju,iu[k]*sizeof(PetscInt));CHKERRQ(ierr);
15310c27e3dSHong Zhang       ierr = PetscFree(ju);CHKERRQ(ierr);
15410c27e3dSHong Zhang       ju   = jutmp;
15510c27e3dSHong Zhang       reallocs++; /* count how many times we realloc */
15610c27e3dSHong Zhang     }
15710c27e3dSHong Zhang 
15810c27e3dSHong Zhang     /* save nonzero structure of k-th row in ju */
15910c27e3dSHong Zhang     i=k;
16010c27e3dSHong Zhang     while (nzk --) {
16110c27e3dSHong Zhang       i           = q[i];
16210c27e3dSHong Zhang       ju[juidx++] = i;
16310c27e3dSHong Zhang     }
16410c27e3dSHong Zhang   }
16510c27e3dSHong Zhang 
1666cf91177SBarry Smith #if defined(PETSC_USE_INFO)
16710c27e3dSHong Zhang   if (ai[mbs] != 0) {
16810c27e3dSHong Zhang     PetscReal af = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]);
169ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr);
170ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
171ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr);
172ae15b995SBarry Smith     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
17310c27e3dSHong Zhang   } else {
174ae15b995SBarry Smith     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
17510c27e3dSHong Zhang   }
17663ba0a88SBarry Smith #endif
17710c27e3dSHong Zhang 
17810c27e3dSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
179d8c74875SBarry Smith   ierr = PetscFree2(jl,q);CHKERRQ(ierr);
18010c27e3dSHong Zhang 
18110c27e3dSHong Zhang   /* put together the new matrix */
182719d5645SBarry Smith   ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(F,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
18310c27e3dSHong Zhang 
184719d5645SBarry Smith   /* ierr = PetscLogObjectParent(B,iperm);CHKERRQ(ierr); */
185719d5645SBarry Smith   b = (Mat_SeqSBAIJ*)(F)->data;
18610c27e3dSHong Zhang   b->singlemalloc = PETSC_FALSE;
187e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
188e6b907acSBarry Smith   b->free_ij       = PETSC_TRUE;
18910c27e3dSHong Zhang   ierr = PetscMalloc((iu[mbs]+1)*sizeof(MatScalar)*bs2,&b->a);CHKERRQ(ierr);
19010c27e3dSHong Zhang   b->j    = ju;
19110c27e3dSHong Zhang   b->i    = iu;
19210c27e3dSHong Zhang   b->diag = 0;
19310c27e3dSHong Zhang   b->ilen = 0;
19410c27e3dSHong Zhang   b->imax = 0;
19510c27e3dSHong Zhang   b->row  = perm;
19610c27e3dSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
19710c27e3dSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
19810c27e3dSHong Zhang   b->icol = perm;
19910c27e3dSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
20010c27e3dSHong Zhang   ierr    = PetscMalloc((bs*mbs+bs)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
20110c27e3dSHong Zhang   /* In b structure:  Free imax, ilen, old a, old j.
20210c27e3dSHong Zhang      Allocate idnew, solve_work, new a, new j */
203719d5645SBarry Smith   ierr = PetscLogObjectMemory(F,(iu[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
20410c27e3dSHong Zhang   b->maxnz = b->nz = iu[mbs];
20510c27e3dSHong Zhang 
206719d5645SBarry Smith   (F)->info.factor_mallocs    = reallocs;
207719d5645SBarry Smith   (F)->info.fill_ratio_given  = f;
20810c27e3dSHong Zhang   if (ai[mbs] != 0) {
209719d5645SBarry Smith     (F)->info.fill_ratio_needed = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]);
21010c27e3dSHong Zhang   } else {
211719d5645SBarry Smith     (F)->info.fill_ratio_needed = 0.0;
21210c27e3dSHong Zhang   }
213719d5645SBarry Smith   ierr = MatSeqSBAIJSetNumericFactorization(F,perm_identity);CHKERRQ(ierr);
21410c27e3dSHong Zhang   PetscFunctionReturn(0);
21510c27e3dSHong Zhang }
21610c27e3dSHong Zhang /*
21710c27e3dSHong Zhang     Symbolic U^T*D*U factorization for SBAIJ format.
21810c27e3dSHong Zhang */
21998a8e78dSHong Zhang #include "petscbt.h"
2207c4f633dSBarry Smith #include "../src/mat/utils/freespace.h"
2214a2ae208SSatish Balay #undef __FUNCT__
2224a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqSBAIJ"
2230481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
22449b5e25fSSatish Balay {
22598a8e78dSHong Zhang   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data;
22698a8e78dSHong Zhang   Mat_SeqSBAIJ       *b;
2276849ba73SBarry Smith   PetscErrorCode     ierr;
22858ebbce7SBarry Smith   PetscTruth         perm_identity,missing;
22998a8e78dSHong Zhang   PetscReal          fill = info->fill;
2305d0c19d7SBarry Smith   const PetscInt     *rip,*ai,*aj;
2315d0c19d7SBarry Smith   PetscInt           i,mbs=a->mbs,bs=A->rmap->bs,reallocs=0,prow,d;
23298a8e78dSHong Zhang   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
2337625dc9aSHong Zhang   PetscInt           nlnk,*lnk,ncols,*cols,*uj,**ui_ptr,*uj_ptr;
234a1a86e44SBarry Smith   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
23598a8e78dSHong Zhang   PetscBT            lnkbt;
23649b5e25fSSatish Balay 
23749b5e25fSSatish Balay   PetscFunctionBegin;
23858ebbce7SBarry Smith   ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr);
23958ebbce7SBarry Smith   if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
24058ebbce7SBarry Smith 
24110c27e3dSHong Zhang   /*
24210c27e3dSHong Zhang    This code originally uses Modified Sparse Row (MSR) storage
24310c27e3dSHong Zhang    (see page 85, "Iterative Methods ..." by Saad) for the output matrix B - bad choise!
24410c27e3dSHong Zhang    Then it is rewritten so the factor B takes seqsbaij format. However the associated
24510c27e3dSHong Zhang    MatCholeskyFactorNumeric_() have not been modified for the cases of bs>1 or !perm_identity,
24610c27e3dSHong Zhang    thus the original code in MSR format is still used for these cases.
24710c27e3dSHong Zhang    The code below should replace MatCholeskyFactorSymbolic_SeqSBAIJ_MSR() whenever
24810c27e3dSHong Zhang    MatCholeskyFactorNumeric_() is modified for using sbaij symbolic factor.
24910c27e3dSHong Zhang   */
250fff829cfSHong Zhang   if (bs > 1){
251719d5645SBarry Smith     ierr = MatCholeskyFactorSymbolic_SeqSBAIJ_MSR(fact,A,perm,info);CHKERRQ(ierr);
252fff829cfSHong Zhang     PetscFunctionReturn(0);
253fff829cfSHong Zhang   }
25410c27e3dSHong Zhang 
25598a8e78dSHong Zhang   /* check whether perm is the identity mapping */
25698a8e78dSHong Zhang   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
25798a8e78dSHong Zhang 
258fff829cfSHong Zhang   if (perm_identity){
259fff829cfSHong Zhang     a->permute = PETSC_FALSE;
260fff829cfSHong Zhang     ai = a->i; aj = a->j;
261fff829cfSHong Zhang   } else {
262d33dd558SHong Zhang     SETERRQ(PETSC_ERR_SUP,"Matrix reordering is not supported for sbaij matrix. Use aij format");
263d33dd558SHong Zhang     /* There are bugs for reordeing. Needs further work.
264d33dd558SHong Zhang        MatReordering for sbaij cannot be efficient. User should use aij formt! */
265fff829cfSHong Zhang     a->permute = PETSC_TRUE;
266fff829cfSHong Zhang     ierr = MatReorderingSeqSBAIJ(A,perm);CHKERRQ(ierr);
267fff829cfSHong Zhang     ai = a->inew; aj = a->jnew;
268fff829cfSHong Zhang   }
269fff829cfSHong Zhang   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
270fff829cfSHong Zhang 
271fff829cfSHong Zhang   /* initialization */
2727625dc9aSHong Zhang   ierr  = PetscMalloc((mbs+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
27398a8e78dSHong Zhang   ui[0] = 0;
27498a8e78dSHong Zhang 
27598a8e78dSHong Zhang   /* jl: linked list for storing indices of the pivot rows
2767625dc9aSHong Zhang      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
277d8c74875SBarry Smith   ierr = PetscMalloc4(mbs,PetscInt*,&ui_ptr,mbs,PetscInt,&il,mbs,PetscInt,&jl,mbs,PetscInt,&cols);CHKERRQ(ierr);
2787625dc9aSHong Zhang   for (i=0; i<mbs; i++){
2797625dc9aSHong Zhang     jl[i] = mbs; il[i] = 0;
280fff829cfSHong Zhang   }
281fff829cfSHong Zhang 
28298a8e78dSHong Zhang   /* create and initialize a linked list for storing column indices of the active row k */
2837625dc9aSHong Zhang   nlnk = mbs + 1;
2847625dc9aSHong Zhang   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
285fff829cfSHong Zhang 
2867625dc9aSHong Zhang   /* initial FreeSpace size is fill*(ai[mbs]+1) */
287a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[mbs]+1)),&free_space);CHKERRQ(ierr);
28898a8e78dSHong Zhang   current_space = free_space;
28998a8e78dSHong Zhang 
2907625dc9aSHong Zhang   for (k=0; k<mbs; k++){  /* for each active row k */
29198a8e78dSHong Zhang     /* initialize lnk by the column indices of row rip[k] of A */
29298a8e78dSHong Zhang     nzk   = 0;
29398a8e78dSHong Zhang     ncols = ai[rip[k]+1] - ai[rip[k]];
2947625dc9aSHong Zhang     for (j=0; j<ncols; j++){
2957625dc9aSHong Zhang       i = *(aj + ai[rip[k]] + j);
2967625dc9aSHong Zhang       cols[j] = rip[i];
2977625dc9aSHong Zhang     }
2987625dc9aSHong Zhang     ierr = PetscLLAdd(ncols,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
29998a8e78dSHong Zhang     nzk += nlnk;
30098a8e78dSHong Zhang 
30198a8e78dSHong Zhang     /* update lnk by computing fill-in for each pivot row to be merged in */
30298a8e78dSHong Zhang     prow = jl[k]; /* 1st pivot row */
303fff829cfSHong Zhang 
304fff829cfSHong Zhang     while (prow < k){
305fff829cfSHong Zhang       nextprow = jl[prow];
30698a8e78dSHong Zhang       /* merge prow into k-th row */
3077625dc9aSHong Zhang       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
30898a8e78dSHong Zhang       jmax = ui[prow+1];
30998a8e78dSHong Zhang       ncols = jmax-jmin;
3107625dc9aSHong Zhang       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
3115a8e39fbSHong Zhang       ierr = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
31298a8e78dSHong Zhang       nzk += nlnk;
313fff829cfSHong Zhang 
31498a8e78dSHong Zhang       /* update il and jl for prow */
315fff829cfSHong Zhang       if (jmin < jmax){
316fff829cfSHong Zhang         il[prow] = jmin;
3177625dc9aSHong Zhang         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
318fff829cfSHong Zhang       }
319fff829cfSHong Zhang       prow = nextprow;
320fff829cfSHong Zhang     }
321fff829cfSHong Zhang 
32298a8e78dSHong Zhang     /* if free space is not available, make more free space */
32398a8e78dSHong Zhang     if (current_space->local_remaining<nzk) {
3247625dc9aSHong Zhang       i = mbs - k + 1; /* num of unfactored rows */
32598a8e78dSHong Zhang       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
326a1a86e44SBarry Smith       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
32798a8e78dSHong Zhang       reallocs++;
32898a8e78dSHong Zhang     }
32998a8e78dSHong Zhang 
33098a8e78dSHong Zhang     /* copy data into free space, then initialize lnk */
3317625dc9aSHong Zhang     ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
33298a8e78dSHong Zhang 
33398a8e78dSHong Zhang     /* add the k-th row into il and jl */
33498a8e78dSHong Zhang     if (nzk-1 > 0){
3357625dc9aSHong Zhang       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
336fff829cfSHong Zhang       jl[k] = jl[i]; jl[i] = k;
33798a8e78dSHong Zhang       il[k] = ui[k] + 1;
338fff829cfSHong Zhang     }
3397625dc9aSHong Zhang     ui_ptr[k] = current_space->array;
34098a8e78dSHong Zhang     current_space->array           += nzk;
34198a8e78dSHong Zhang     current_space->local_used      += nzk;
34298a8e78dSHong Zhang     current_space->local_remaining -= nzk;
343fff829cfSHong Zhang 
34498a8e78dSHong Zhang     ui[k+1] = ui[k] + nzk;
345fff829cfSHong Zhang   }
346fff829cfSHong Zhang 
3476cf91177SBarry Smith #if defined(PETSC_USE_INFO)
3487625dc9aSHong Zhang   if (ai[mbs] != 0) {
3497625dc9aSHong Zhang     PetscReal af = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
350ae15b995SBarry Smith     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
351ae15b995SBarry Smith     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
352ae15b995SBarry Smith     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
353fff829cfSHong Zhang   } else {
354ae15b995SBarry Smith     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
355fff829cfSHong Zhang   }
35663ba0a88SBarry Smith #endif
357fff829cfSHong Zhang 
358fff829cfSHong Zhang   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
359d8c74875SBarry Smith   ierr = PetscFree4(ui_ptr,il,jl,cols);CHKERRQ(ierr);
360fff829cfSHong Zhang 
36198a8e78dSHong Zhang   /* destroy list of free space and other temporary array(s) */
3627625dc9aSHong Zhang   ierr = PetscMalloc((ui[mbs]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
363a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
36498a8e78dSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
365fff829cfSHong Zhang 
36698a8e78dSHong Zhang   /* put together the new matrix in MATSEQSBAIJ format */
367719d5645SBarry Smith   ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(fact,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
36898a8e78dSHong Zhang 
369719d5645SBarry Smith   b = (Mat_SeqSBAIJ*)(fact)->data;
370fff829cfSHong Zhang   b->singlemalloc = PETSC_FALSE;
371e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
372e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
3737625dc9aSHong Zhang   ierr = PetscMalloc((ui[mbs]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
37498a8e78dSHong Zhang   b->j    = uj;
37598a8e78dSHong Zhang   b->i    = ui;
376fff829cfSHong Zhang   b->diag = 0;
377fff829cfSHong Zhang   b->ilen = 0;
378fff829cfSHong Zhang   b->imax = 0;
379fff829cfSHong Zhang   b->row  = perm;
380fff829cfSHong Zhang   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
381fff829cfSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
382fff829cfSHong Zhang   b->icol = perm;
383fff829cfSHong Zhang   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
3847625dc9aSHong Zhang   ierr    = PetscMalloc((mbs+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
385719d5645SBarry Smith   ierr    = PetscLogObjectMemory(fact,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
3867625dc9aSHong Zhang   b->maxnz = b->nz = ui[mbs];
387fff829cfSHong Zhang 
388719d5645SBarry Smith   (fact)->info.factor_mallocs    = reallocs;
389719d5645SBarry Smith   (fact)->info.fill_ratio_given  = fill;
3907625dc9aSHong Zhang   if (ai[mbs] != 0) {
391719d5645SBarry Smith     (fact)->info.fill_ratio_needed = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
392fff829cfSHong Zhang   } else {
393719d5645SBarry Smith     (fact)->info.fill_ratio_needed = 0.0;
394fff829cfSHong Zhang   }
395719d5645SBarry Smith   ierr = MatSeqSBAIJSetNumericFactorization(fact,perm_identity);CHKERRQ(ierr);
396064503c5SHong Zhang   PetscFunctionReturn(0);
397064503c5SHong Zhang }
3984a2ae208SSatish Balay #undef __FUNCT__
3994a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N"
4000481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat C,Mat A,const MatFactorInfo *info)
40149b5e25fSSatish Balay {
4024c16a6a6SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
4034c16a6a6SHong Zhang   IS             perm = b->row;
4046849ba73SBarry Smith   PetscErrorCode ierr;
4055d0c19d7SBarry Smith   const PetscInt *ai,*aj,*perm_ptr,mbs=a->mbs,*bi=b->i,*bj=b->j;
4065d0c19d7SBarry Smith   PetscInt       i,j;
4075d0c19d7SBarry Smith   PetscInt       *a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
408d0f46423SBarry Smith   PetscInt       bs=A->rmap->bs,bs2 = a->bs2,bslog = 0;
4094c16a6a6SHong Zhang   MatScalar      *ba = b->a,*aa,*ap,*dk,*uik;
4104c16a6a6SHong Zhang   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
41128de702eSHong Zhang   MatScalar      *work;
41213f74950SBarry Smith   PetscInt       *pivots;
4134c16a6a6SHong Zhang 
4144c16a6a6SHong Zhang   PetscFunctionBegin;
4154c16a6a6SHong Zhang   /* initialization */
41682502324SSatish Balay   ierr = PetscMalloc(bs2*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
4174c16a6a6SHong Zhang   ierr = PetscMemzero(rtmp,bs2*mbs*sizeof(MatScalar));CHKERRQ(ierr);
418d8c74875SBarry Smith   ierr = PetscMalloc2(mbs,PetscInt,&il,mbs,PetscInt,&jl);CHKERRQ(ierr);
4194c16a6a6SHong Zhang   for (i=0; i<mbs; i++) {
4204c16a6a6SHong Zhang     jl[i] = mbs; il[0] = 0;
4214c16a6a6SHong Zhang   }
422d8c74875SBarry Smith   ierr = PetscMalloc3(bs2,MatScalar,&dk,bs2,MatScalar,&uik,bs,MatScalar,&work);CHKERRQ(ierr);
42313f74950SBarry Smith   ierr = PetscMalloc(bs*sizeof(PetscInt),&pivots);CHKERRQ(ierr);
4244c16a6a6SHong Zhang 
4254c16a6a6SHong Zhang   ierr  = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr);
4264c16a6a6SHong Zhang 
4274c16a6a6SHong Zhang   /* check permutation */
4284c16a6a6SHong Zhang   if (!a->permute){
4294c16a6a6SHong Zhang     ai = a->i; aj = a->j; aa = a->a;
4304c16a6a6SHong Zhang   } else {
4314c16a6a6SHong Zhang     ai   = a->inew; aj = a->jnew;
43282502324SSatish Balay     ierr = PetscMalloc(bs2*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr);
4334c16a6a6SHong Zhang     ierr = PetscMemcpy(aa,a->a,bs2*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr);
43413f74950SBarry Smith     ierr = PetscMalloc(ai[mbs]*sizeof(PetscInt),&a2anew);CHKERRQ(ierr);
43513f74950SBarry Smith     ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(PetscInt));CHKERRQ(ierr);
4364c16a6a6SHong Zhang 
437187a9f4bSHong Zhang     /* flops in while loop */
438187a9f4bSHong Zhang     bslog = 2*bs*bs2;
439187a9f4bSHong Zhang 
4404c16a6a6SHong Zhang     for (i=0; i<mbs; i++){
4414c16a6a6SHong Zhang       jmin = ai[i]; jmax = ai[i+1];
4424c16a6a6SHong Zhang       for (j=jmin; j<jmax; j++){
4434c16a6a6SHong Zhang         while (a2anew[j] != j){
4444c16a6a6SHong Zhang           k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
4454c16a6a6SHong Zhang           for (k1=0; k1<bs2; k1++){
4464c16a6a6SHong Zhang             dk[k1]       = aa[k*bs2+k1];
4474c16a6a6SHong Zhang             aa[k*bs2+k1] = aa[j*bs2+k1];
4484c16a6a6SHong Zhang             aa[j*bs2+k1] = dk[k1];
4494c16a6a6SHong Zhang           }
4504c16a6a6SHong Zhang         }
4514c16a6a6SHong Zhang         /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
4524c16a6a6SHong Zhang         if (i > aj[j]){
4534c16a6a6SHong Zhang           /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */
4544c16a6a6SHong Zhang           ap = aa + j*bs2;                     /* ptr to the beginning of j-th block of aa */
4554c16a6a6SHong Zhang           for (k=0; k<bs2; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */
4564c16a6a6SHong Zhang           for (k=0; k<bs; k++){               /* j-th block of aa <- dk^T */
4574c16a6a6SHong Zhang             for (k1=0; k1<bs; k1++) *ap++ = dk[k + bs*k1];
4584c16a6a6SHong Zhang           }
4594c16a6a6SHong Zhang         }
4604c16a6a6SHong Zhang       }
4614c16a6a6SHong Zhang     }
462323b833fSBarry Smith     ierr = PetscFree(a2anew);CHKERRQ(ierr);
4634c16a6a6SHong Zhang   }
4644c16a6a6SHong Zhang 
4654c16a6a6SHong Zhang   /* for each row k */
4664c16a6a6SHong Zhang   for (k = 0; k<mbs; k++){
4674c16a6a6SHong Zhang 
4684c16a6a6SHong Zhang     /*initialize k-th row with elements nonzero in row perm(k) of A */
4694c16a6a6SHong Zhang     jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];
470057f5ba7SHong Zhang 
4714c16a6a6SHong Zhang     ap = aa + jmin*bs2;
4724c16a6a6SHong Zhang     for (j = jmin; j < jmax; j++){
4734c16a6a6SHong Zhang       vj = perm_ptr[aj[j]];         /* block col. index */
4744c16a6a6SHong Zhang       rtmp_ptr = rtmp + vj*bs2;
4754c16a6a6SHong Zhang       for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++;
4764c16a6a6SHong Zhang     }
4774c16a6a6SHong Zhang 
4784c16a6a6SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
4794c16a6a6SHong Zhang     ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
4804c16a6a6SHong Zhang     i = jl[k]; /* first row to be added to k_th row  */
4814c16a6a6SHong Zhang 
482057f5ba7SHong Zhang     while (i < k){
4834c16a6a6SHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
4844c16a6a6SHong Zhang 
4854c16a6a6SHong Zhang       /* compute multiplier */
4864c16a6a6SHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
4874c16a6a6SHong Zhang 
4884c16a6a6SHong Zhang       /* uik = -inv(Di)*U_bar(i,k) */
4894c16a6a6SHong Zhang       diag = ba + i*bs2;
4904c16a6a6SHong Zhang       u    = ba + ili*bs2;
4914c16a6a6SHong Zhang       ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
4924c16a6a6SHong Zhang       Kernel_A_gets_A_minus_B_times_C(bs,uik,diag,u);
4934c16a6a6SHong Zhang 
4944c16a6a6SHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
4954c16a6a6SHong Zhang       Kernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u);
496dc0b31edSSatish Balay       ierr = PetscLogFlops(bslog*2.0);CHKERRQ(ierr);
4974c16a6a6SHong Zhang 
4984c16a6a6SHong Zhang       /* update -U(i,k) */
4994c16a6a6SHong Zhang       ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
5004c16a6a6SHong Zhang 
5014c16a6a6SHong Zhang       /* add multiple of row i to k-th row ... */
5024c16a6a6SHong Zhang       jmin = ili + 1; jmax = bi[i+1];
5034c16a6a6SHong Zhang       if (jmin < jmax){
5044c16a6a6SHong Zhang         for (j=jmin; j<jmax; j++) {
5054c16a6a6SHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j) */
5064c16a6a6SHong Zhang           rtmp_ptr = rtmp + bj[j]*bs2;
5074c16a6a6SHong Zhang           u = ba + j*bs2;
5084c16a6a6SHong Zhang           Kernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u);
5094c16a6a6SHong Zhang         }
510187a9f4bSHong Zhang         ierr = PetscLogFlops(bslog*(jmax-jmin));CHKERRQ(ierr);
5114c16a6a6SHong Zhang 
5124c16a6a6SHong Zhang         /* ... add i to row list for next nonzero entry */
5134c16a6a6SHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
5144c16a6a6SHong Zhang         j     = bj[jmin];
5154c16a6a6SHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
5164c16a6a6SHong Zhang       }
5174c16a6a6SHong Zhang       i = nexti;
5184c16a6a6SHong Zhang     }
5194c16a6a6SHong Zhang 
5204c16a6a6SHong Zhang     /* save nonzero entries in k-th row of U ... */
5214c16a6a6SHong Zhang 
5224c16a6a6SHong Zhang     /* invert diagonal block */
5234c16a6a6SHong Zhang     diag = ba+k*bs2;
5244c16a6a6SHong Zhang     ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr);
525d230e6fdSBarry Smith     ierr = Kernel_A_gets_inverse_A(bs,diag,pivots,work);CHKERRQ(ierr);
5264c16a6a6SHong Zhang 
5274c16a6a6SHong Zhang     jmin = bi[k]; jmax = bi[k+1];
5284c16a6a6SHong Zhang     if (jmin < jmax) {
5294c16a6a6SHong Zhang       for (j=jmin; j<jmax; j++){
5304c16a6a6SHong Zhang          vj = bj[j];           /* block col. index of U */
5314c16a6a6SHong Zhang          u   = ba + j*bs2;
5324c16a6a6SHong Zhang          rtmp_ptr = rtmp + vj*bs2;
5334c16a6a6SHong Zhang          for (k1=0; k1<bs2; k1++){
5344c16a6a6SHong Zhang            *u++        = *rtmp_ptr;
5354c16a6a6SHong Zhang            *rtmp_ptr++ = 0.0;
5364c16a6a6SHong Zhang          }
5374c16a6a6SHong Zhang       }
5384c16a6a6SHong Zhang 
5394c16a6a6SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
5404c16a6a6SHong Zhang       il[k] = jmin;
5414c16a6a6SHong Zhang       i     = bj[jmin];
5424c16a6a6SHong Zhang       jl[k] = jl[i]; jl[i] = k;
5434c16a6a6SHong Zhang     }
5444c16a6a6SHong Zhang   }
5454c16a6a6SHong Zhang 
5464c16a6a6SHong Zhang   ierr = PetscFree(rtmp);CHKERRQ(ierr);
547d8c74875SBarry Smith   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
548d8c74875SBarry Smith   ierr = PetscFree3(dk,uik,work);CHKERRQ(ierr);
5494c16a6a6SHong Zhang   ierr = PetscFree(pivots);CHKERRQ(ierr);
5504c16a6a6SHong Zhang   if (a->permute){
5514c16a6a6SHong Zhang     ierr = PetscFree(aa);CHKERRQ(ierr);
5524c16a6a6SHong Zhang   }
5534c16a6a6SHong Zhang 
5544c16a6a6SHong Zhang   ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr);
555db4efbfdSBarry Smith   C->ops->solve          = MatSolve_SeqSBAIJ_N;
556db4efbfdSBarry Smith   C->ops->solvetranspose = MatSolve_SeqSBAIJ_N;
557db4efbfdSBarry Smith   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_N;
558db4efbfdSBarry Smith   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_N;
559db4efbfdSBarry Smith 
5604c16a6a6SHong Zhang   C->assembled    = PETSC_TRUE;
5614c16a6a6SHong Zhang   C->preallocated = PETSC_TRUE;
562efee365bSSatish Balay   ierr = PetscLogFlops(1.3333*bs*bs2*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
5634c16a6a6SHong Zhang   PetscFunctionReturn(0);
5644c16a6a6SHong Zhang }
565d4edadd4SHong Zhang 
5664a2ae208SSatish Balay #undef __FUNCT__
5674a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering"
5680481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info)
569671cb588SHong Zhang {
570671cb588SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
571dfbe8321SBarry Smith   PetscErrorCode ierr;
57213f74950SBarry Smith   PetscInt       i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
57313f74950SBarry Smith   PetscInt       *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
574d0f46423SBarry Smith   PetscInt       bs=A->rmap->bs,bs2 = a->bs2,bslog;
575671cb588SHong Zhang   MatScalar      *ba = b->a,*aa,*ap,*dk,*uik;
576671cb588SHong Zhang   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
57728de702eSHong Zhang   MatScalar      *work;
57813f74950SBarry Smith   PetscInt       *pivots;
579671cb588SHong Zhang 
580671cb588SHong Zhang   PetscFunctionBegin;
58182502324SSatish Balay   ierr = PetscMalloc(bs2*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
582671cb588SHong Zhang   ierr = PetscMemzero(rtmp,bs2*mbs*sizeof(MatScalar));CHKERRQ(ierr);
583d8c74875SBarry Smith   ierr = PetscMalloc2(mbs,PetscInt,&il,mbs,PetscInt,&jl);CHKERRQ(ierr);
584671cb588SHong Zhang   for (i=0; i<mbs; i++) {
585671cb588SHong Zhang     jl[i] = mbs; il[0] = 0;
586671cb588SHong Zhang   }
587d8c74875SBarry Smith   ierr = PetscMalloc3(bs2,MatScalar,&dk,bs2,MatScalar,&uik,bs,MatScalar,&work);CHKERRQ(ierr);
58813f74950SBarry Smith   ierr = PetscMalloc(bs*sizeof(PetscInt),&pivots);CHKERRQ(ierr);
589671cb588SHong Zhang 
590671cb588SHong Zhang   ai = a->i; aj = a->j; aa = a->a;
591671cb588SHong Zhang 
592187a9f4bSHong Zhang   /* flops in while loop */
593187a9f4bSHong Zhang   bslog = 2*bs*bs2;
594187a9f4bSHong Zhang 
595671cb588SHong Zhang   /* for each row k */
596671cb588SHong Zhang   for (k = 0; k<mbs; k++){
597671cb588SHong Zhang 
598671cb588SHong Zhang     /*initialize k-th row with elements nonzero in row k of A */
599671cb588SHong Zhang     jmin = ai[k]; jmax = ai[k+1];
600671cb588SHong Zhang     ap = aa + jmin*bs2;
601671cb588SHong Zhang     for (j = jmin; j < jmax; j++){
602671cb588SHong Zhang       vj = aj[j];         /* block col. index */
603671cb588SHong Zhang       rtmp_ptr = rtmp + vj*bs2;
604671cb588SHong Zhang       for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++;
605671cb588SHong Zhang     }
606671cb588SHong Zhang 
607671cb588SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
608671cb588SHong Zhang     ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
609671cb588SHong Zhang     i = jl[k]; /* first row to be added to k_th row  */
610671cb588SHong Zhang 
611057f5ba7SHong Zhang     while (i < k){
612671cb588SHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
613671cb588SHong Zhang 
614671cb588SHong Zhang       /* compute multiplier */
615671cb588SHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
616671cb588SHong Zhang 
617671cb588SHong Zhang       /* uik = -inv(Di)*U_bar(i,k) */
618671cb588SHong Zhang       diag = ba + i*bs2;
619671cb588SHong Zhang       u    = ba + ili*bs2;
620671cb588SHong Zhang       ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
621671cb588SHong Zhang       Kernel_A_gets_A_minus_B_times_C(bs,uik,diag,u);
622671cb588SHong Zhang 
623671cb588SHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
624671cb588SHong Zhang       Kernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u);
625dc0b31edSSatish Balay       ierr = PetscLogFlops(bslog*2.0);CHKERRQ(ierr);
626671cb588SHong Zhang 
627671cb588SHong Zhang       /* update -U(i,k) */
628671cb588SHong Zhang       ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
629671cb588SHong Zhang 
630671cb588SHong Zhang       /* add multiple of row i to k-th row ... */
631671cb588SHong Zhang       jmin = ili + 1; jmax = bi[i+1];
632671cb588SHong Zhang       if (jmin < jmax){
633671cb588SHong Zhang         for (j=jmin; j<jmax; j++) {
634671cb588SHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j) */
635671cb588SHong Zhang           rtmp_ptr = rtmp + bj[j]*bs2;
636671cb588SHong Zhang           u = ba + j*bs2;
637671cb588SHong Zhang           Kernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u);
638671cb588SHong Zhang         }
639187a9f4bSHong Zhang         ierr = PetscLogFlops(bslog*(jmax-jmin));CHKERRQ(ierr);
640671cb588SHong Zhang 
641671cb588SHong Zhang         /* ... add i to row list for next nonzero entry */
642671cb588SHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
643671cb588SHong Zhang         j     = bj[jmin];
644671cb588SHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
645671cb588SHong Zhang       }
646671cb588SHong Zhang       i = nexti;
647671cb588SHong Zhang     }
648671cb588SHong Zhang 
649671cb588SHong Zhang     /* save nonzero entries in k-th row of U ... */
650671cb588SHong Zhang 
651671cb588SHong Zhang     /* invert diagonal block */
652671cb588SHong Zhang     diag = ba+k*bs2;
653671cb588SHong Zhang     ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr);
654d230e6fdSBarry Smith     ierr = Kernel_A_gets_inverse_A(bs,diag,pivots,work);CHKERRQ(ierr);
655671cb588SHong Zhang 
656671cb588SHong Zhang     jmin = bi[k]; jmax = bi[k+1];
657671cb588SHong Zhang     if (jmin < jmax) {
658671cb588SHong Zhang       for (j=jmin; j<jmax; j++){
659671cb588SHong Zhang          vj = bj[j];           /* block col. index of U */
660671cb588SHong Zhang          u   = ba + j*bs2;
661671cb588SHong Zhang          rtmp_ptr = rtmp + vj*bs2;
662671cb588SHong Zhang          for (k1=0; k1<bs2; k1++){
663671cb588SHong Zhang            *u++        = *rtmp_ptr;
664671cb588SHong Zhang            *rtmp_ptr++ = 0.0;
665671cb588SHong Zhang          }
666671cb588SHong Zhang       }
667671cb588SHong Zhang 
668671cb588SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
669671cb588SHong Zhang       il[k] = jmin;
670671cb588SHong Zhang       i     = bj[jmin];
671671cb588SHong Zhang       jl[k] = jl[i]; jl[i] = k;
672671cb588SHong Zhang     }
673671cb588SHong Zhang   }
674671cb588SHong Zhang 
675671cb588SHong Zhang   ierr = PetscFree(rtmp);CHKERRQ(ierr);
676d8c74875SBarry Smith   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
677d8c74875SBarry Smith   ierr = PetscFree3(dk,uik,work);CHKERRQ(ierr);
678671cb588SHong Zhang   ierr = PetscFree(pivots);CHKERRQ(ierr);
679671cb588SHong Zhang 
680db4efbfdSBarry Smith   C->ops->solve          = MatSolve_SeqSBAIJ_N_NaturalOrdering;
681db4efbfdSBarry Smith   C->ops->solvetranspose = MatSolve_SeqSBAIJ_N_NaturalOrdering;
682db4efbfdSBarry Smith   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_N_NaturalOrdering;
683db4efbfdSBarry Smith   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_N_NaturalOrdering;
684671cb588SHong Zhang   C->assembled = PETSC_TRUE;
685671cb588SHong Zhang   C->preallocated = PETSC_TRUE;
686efee365bSSatish Balay   ierr = PetscLogFlops(1.3333*bs*bs2*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
687671cb588SHong Zhang   PetscFunctionReturn(0);
688671cb588SHong Zhang }
689671cb588SHong Zhang 
69049b5e25fSSatish Balay /*
691fcf159c0SHong Zhang     Numeric U^T*D*U factorization for SBAIJ format. Modified from SNF of YSMP.
692cc0c071aSHong Zhang     Version for blocks 2 by 2.
69349b5e25fSSatish Balay */
6944a2ae208SSatish Balay #undef __FUNCT__
6954a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2"
6960481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat C,Mat A,const MatFactorInfo *info)
69749b5e25fSSatish Balay {
69891602c66SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
699cc0c071aSHong Zhang   IS             perm = b->row;
7006849ba73SBarry Smith   PetscErrorCode ierr;
7015d0c19d7SBarry Smith   const PetscInt *ai,*aj,*perm_ptr;
7025d0c19d7SBarry Smith   PetscInt       i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
7035d0c19d7SBarry Smith   PetscInt       *a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
704d8c74875SBarry Smith   MatScalar      *ba = b->a,*aa,*ap;
705d8c74875SBarry Smith   MatScalar      *u,*diag,*rtmp,*rtmp_ptr,dk[4],uik[4];
70662bba022SBarry Smith   PetscReal      shift = info->shiftinblocks;
70749b5e25fSSatish Balay 
70849b5e25fSSatish Balay   PetscFunctionBegin;
70991602c66SHong Zhang   /* initialization */
71091602c66SHong Zhang   /* il and jl record the first nonzero element in each row of the accessing
71191602c66SHong Zhang      window U(0:k, k:mbs-1).
71291602c66SHong Zhang      jl:    list of rows to be added to uneliminated rows
71391602c66SHong Zhang             i>= k: jl(i) is the first row to be added to row i
71491602c66SHong Zhang             i<  k: jl(i) is the row following row i in some list of rows
71591602c66SHong Zhang             jl(i) = mbs indicates the end of a list
71691602c66SHong Zhang      il(i): points to the first nonzero element in columns k,...,mbs-1 of
71791602c66SHong Zhang             row i of U */
71882502324SSatish Balay   ierr = PetscMalloc(4*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
719cc0c071aSHong Zhang   ierr = PetscMemzero(rtmp,4*mbs*sizeof(MatScalar));CHKERRQ(ierr);
720d8c74875SBarry Smith   ierr = PetscMalloc2(mbs,PetscInt,&il,mbs,PetscInt,&jl);CHKERRQ(ierr);
72191602c66SHong Zhang   for (i=0; i<mbs; i++) {
7223845f261SHong Zhang     jl[i] = mbs; il[0] = 0;
72391602c66SHong Zhang   }
724cc0c071aSHong Zhang   ierr = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr);
725cc0c071aSHong Zhang 
726cc0c071aSHong Zhang   /* check permutation */
727cc0c071aSHong Zhang   if (!a->permute){
728cc0c071aSHong Zhang     ai = a->i; aj = a->j; aa = a->a;
729cc0c071aSHong Zhang   } else {
730cc0c071aSHong Zhang     ai   = a->inew; aj = a->jnew;
73182502324SSatish Balay     ierr = PetscMalloc(4*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr);
732cc0c071aSHong Zhang     ierr = PetscMemcpy(aa,a->a,4*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr);
73313f74950SBarry Smith     ierr = PetscMalloc(ai[mbs]*sizeof(PetscInt),&a2anew);CHKERRQ(ierr);
73413f74950SBarry Smith     ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(PetscInt));CHKERRQ(ierr);
735cc0c071aSHong Zhang 
736cc0c071aSHong Zhang     for (i=0; i<mbs; i++){
737cc0c071aSHong Zhang       jmin = ai[i]; jmax = ai[i+1];
738cc0c071aSHong Zhang       for (j=jmin; j<jmax; j++){
739cc0c071aSHong Zhang         while (a2anew[j] != j){
740cc0c071aSHong Zhang           k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
741cc0c071aSHong Zhang           for (k1=0; k1<4; k1++){
742cc0c071aSHong Zhang             dk[k1]       = aa[k*4+k1];
743cc0c071aSHong Zhang             aa[k*4+k1] = aa[j*4+k1];
744cc0c071aSHong Zhang             aa[j*4+k1] = dk[k1];
745cc0c071aSHong Zhang           }
746cc0c071aSHong Zhang         }
747cc0c071aSHong Zhang         /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
748cc0c071aSHong Zhang         if (i > aj[j]){
749a1723e09SHong Zhang           /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */
750cc0c071aSHong Zhang           ap = aa + j*4;     /* ptr to the beginning of the block */
751cc0c071aSHong Zhang           dk[1] = ap[1];     /* swap ap[1] and ap[2] */
752cc0c071aSHong Zhang           ap[1] = ap[2];
753cc0c071aSHong Zhang           ap[2] = dk[1];
754cc0c071aSHong Zhang         }
755cc0c071aSHong Zhang       }
756cc0c071aSHong Zhang     }
757ac355199SBarry Smith     ierr = PetscFree(a2anew);CHKERRQ(ierr);
758cc0c071aSHong Zhang   }
7593845f261SHong Zhang 
76091602c66SHong Zhang   /* for each row k */
76191602c66SHong Zhang   for (k = 0; k<mbs; k++){
76291602c66SHong Zhang 
76391602c66SHong Zhang     /*initialize k-th row with elements nonzero in row perm(k) of A */
764cc0c071aSHong Zhang     jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];
765cc0c071aSHong Zhang     ap = aa + jmin*4;
76691602c66SHong Zhang     for (j = jmin; j < jmax; j++){
767cc0c071aSHong Zhang       vj = perm_ptr[aj[j]];         /* block col. index */
768cc0c071aSHong Zhang       rtmp_ptr = rtmp + vj*4;
769cc0c071aSHong Zhang       for (i=0; i<4; i++) *rtmp_ptr++ = *ap++;
77091602c66SHong Zhang     }
77191602c66SHong Zhang 
77291602c66SHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
773cc0c071aSHong Zhang     ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr);
77491602c66SHong Zhang     i = jl[k]; /* first row to be added to k_th row  */
77591602c66SHong Zhang 
776057f5ba7SHong Zhang     while (i < k){
77791602c66SHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
77891602c66SHong Zhang 
7793845f261SHong Zhang       /* compute multiplier */
78091602c66SHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
7813845f261SHong Zhang 
7823845f261SHong Zhang       /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */
783cc0c071aSHong Zhang       diag = ba + i*4;
784cc0c071aSHong Zhang       u    = ba + ili*4;
785cc0c071aSHong Zhang       uik[0] = -(diag[0]*u[0] + diag[2]*u[1]);
786cc0c071aSHong Zhang       uik[1] = -(diag[1]*u[0] + diag[3]*u[1]);
787cc0c071aSHong Zhang       uik[2] = -(diag[0]*u[2] + diag[2]*u[3]);
788cc0c071aSHong Zhang       uik[3] = -(diag[1]*u[2] + diag[3]*u[3]);
7893845f261SHong Zhang 
7903845f261SHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */
791cc0c071aSHong Zhang       dk[0] += uik[0]*u[0] + uik[1]*u[1];
792cc0c071aSHong Zhang       dk[1] += uik[2]*u[0] + uik[3]*u[1];
793cc0c071aSHong Zhang       dk[2] += uik[0]*u[2] + uik[1]*u[3];
794cc0c071aSHong Zhang       dk[3] += uik[2]*u[2] + uik[3]*u[3];
7953845f261SHong Zhang 
796dc0b31edSSatish Balay       ierr = PetscLogFlops(16.0*2.0);CHKERRQ(ierr);
797187a9f4bSHong Zhang 
7983845f261SHong Zhang       /* update -U(i,k): ba[ili] = uik */
799cc0c071aSHong Zhang       ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr);
80091602c66SHong Zhang 
80191602c66SHong Zhang       /* add multiple of row i to k-th row ... */
80291602c66SHong Zhang       jmin = ili + 1; jmax = bi[i+1];
80391602c66SHong Zhang       if (jmin < jmax){
8043845f261SHong Zhang         for (j=jmin; j<jmax; j++) {
8053845f261SHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */
806cc0c071aSHong Zhang           rtmp_ptr = rtmp + bj[j]*4;
807cc0c071aSHong Zhang           u = ba + j*4;
808cc0c071aSHong Zhang           rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1];
809cc0c071aSHong Zhang           rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1];
810cc0c071aSHong Zhang           rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3];
811cc0c071aSHong Zhang           rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3];
8123845f261SHong Zhang         }
813dc0b31edSSatish Balay         ierr = PetscLogFlops(16.0*(jmax-jmin));CHKERRQ(ierr);
8143845f261SHong Zhang 
81591602c66SHong Zhang         /* ... add i to row list for next nonzero entry */
81691602c66SHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
81791602c66SHong Zhang         j     = bj[jmin];
81891602c66SHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
81991602c66SHong Zhang       }
820a1723e09SHong Zhang       i = nexti;
82191602c66SHong Zhang     }
822cc0c071aSHong Zhang 
82391602c66SHong Zhang     /* save nonzero entries in k-th row of U ... */
8243845f261SHong Zhang 
825cc0c071aSHong Zhang     /* invert diagonal block */
826cc0c071aSHong Zhang     diag = ba+k*4;
827cc0c071aSHong Zhang     ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr);
82862bba022SBarry Smith     ierr = Kernel_A_gets_inverse_A_2(diag,shift);CHKERRQ(ierr);
8293845f261SHong Zhang 
83091602c66SHong Zhang     jmin = bi[k]; jmax = bi[k+1];
83191602c66SHong Zhang     if (jmin < jmax) {
83291602c66SHong Zhang       for (j=jmin; j<jmax; j++){
833cc0c071aSHong Zhang          vj = bj[j];           /* block col. index of U */
834cc0c071aSHong Zhang          u   = ba + j*4;
835cc0c071aSHong Zhang          rtmp_ptr = rtmp + vj*4;
836cc0c071aSHong Zhang          for (k1=0; k1<4; k1++){
837cc0c071aSHong Zhang            *u++        = *rtmp_ptr;
838cc0c071aSHong Zhang            *rtmp_ptr++ = 0.0;
839cc0c071aSHong Zhang          }
840cc0c071aSHong Zhang       }
8413845f261SHong Zhang 
84291602c66SHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
84391602c66SHong Zhang       il[k] = jmin;
84491602c66SHong Zhang       i     = bj[jmin];
84591602c66SHong Zhang       jl[k] = jl[i]; jl[i] = k;
84691602c66SHong Zhang     }
84791602c66SHong Zhang   }
8483845f261SHong Zhang 
84949b5e25fSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
850d8c74875SBarry Smith   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
85191602c66SHong Zhang   if (a->permute) {
85291602c66SHong Zhang     ierr = PetscFree(aa);CHKERRQ(ierr);
85391602c66SHong Zhang   }
854cc0c071aSHong Zhang   ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr);
855db4efbfdSBarry Smith   C->ops->solve          = MatSolve_SeqSBAIJ_2;
856db4efbfdSBarry Smith   C->ops->solvetranspose = MatSolve_SeqSBAIJ_2;
85749b5e25fSSatish Balay   C->assembled = PETSC_TRUE;
8585c0bcdfcSHong Zhang   C->preallocated = PETSC_TRUE;
859efee365bSSatish Balay   ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
86049b5e25fSSatish Balay   PetscFunctionReturn(0);
86149b5e25fSSatish Balay }
86291602c66SHong Zhang 
86349b5e25fSSatish Balay /*
86449b5e25fSSatish Balay       Version for when blocks are 2 by 2 Using natural ordering
86549b5e25fSSatish Balay */
8664a2ae208SSatish Balay #undef __FUNCT__
8674a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering"
8680481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info)
86949b5e25fSSatish Balay {
870ab27746eSHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
871dfbe8321SBarry Smith   PetscErrorCode ierr;
87213f74950SBarry Smith   PetscInt       i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
87313f74950SBarry Smith   PetscInt       *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
874d8c74875SBarry Smith   MatScalar      *ba = b->a,*aa,*ap,dk[8],uik[8];
875ab27746eSHong Zhang   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
87662bba022SBarry Smith   PetscReal      shift = info->shiftinblocks;
87749b5e25fSSatish Balay 
87849b5e25fSSatish Balay   PetscFunctionBegin;
879ab27746eSHong Zhang   /* initialization */
880ab27746eSHong Zhang   /* il and jl record the first nonzero element in each row of the accessing
881ab27746eSHong Zhang      window U(0:k, k:mbs-1).
882ab27746eSHong Zhang      jl:    list of rows to be added to uneliminated rows
883ab27746eSHong Zhang             i>= k: jl(i) is the first row to be added to row i
884ab27746eSHong Zhang             i<  k: jl(i) is the row following row i in some list of rows
885ab27746eSHong Zhang             jl(i) = mbs indicates the end of a list
886ab27746eSHong Zhang      il(i): points to the first nonzero element in columns k,...,mbs-1 of
887ab27746eSHong Zhang             row i of U */
88882502324SSatish Balay   ierr = PetscMalloc(4*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
889ab27746eSHong Zhang   ierr = PetscMemzero(rtmp,4*mbs*sizeof(MatScalar));CHKERRQ(ierr);
890d8c74875SBarry Smith   ierr = PetscMalloc2(mbs,PetscInt,&il,mbs,PetscInt,&jl);CHKERRQ(ierr);
891ab27746eSHong Zhang   for (i=0; i<mbs; i++) {
892ab27746eSHong Zhang     jl[i] = mbs; il[0] = 0;
89349b5e25fSSatish Balay   }
894ab27746eSHong Zhang   ai = a->i; aj = a->j; aa = a->a;
895ab27746eSHong Zhang 
896ab27746eSHong Zhang   /* for each row k */
897ab27746eSHong Zhang   for (k = 0; k<mbs; k++){
898ab27746eSHong Zhang 
899ab27746eSHong Zhang     /*initialize k-th row with elements nonzero in row k of A */
900ab27746eSHong Zhang     jmin = ai[k]; jmax = ai[k+1];
901ab27746eSHong Zhang     ap = aa + jmin*4;
902ab27746eSHong Zhang     for (j = jmin; j < jmax; j++){
903ab27746eSHong Zhang       vj = aj[j];         /* block col. index */
904ab27746eSHong Zhang       rtmp_ptr = rtmp + vj*4;
905ab27746eSHong Zhang       for (i=0; i<4; i++) *rtmp_ptr++ = *ap++;
90649b5e25fSSatish Balay     }
907ab27746eSHong Zhang 
908ab27746eSHong Zhang     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
909ab27746eSHong Zhang     ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr);
910ab27746eSHong Zhang     i = jl[k]; /* first row to be added to k_th row  */
911ab27746eSHong Zhang 
912057f5ba7SHong Zhang     while (i < k){
913ab27746eSHong Zhang       nexti = jl[i]; /* next row to be added to k_th row */
914ab27746eSHong Zhang 
915ab27746eSHong Zhang       /* compute multiplier */
916ab27746eSHong Zhang       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
917ab27746eSHong Zhang 
918ab27746eSHong Zhang       /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */
919ab27746eSHong Zhang       diag = ba + i*4;
920ab27746eSHong Zhang       u    = ba + ili*4;
921ab27746eSHong Zhang       uik[0] = -(diag[0]*u[0] + diag[2]*u[1]);
922ab27746eSHong Zhang       uik[1] = -(diag[1]*u[0] + diag[3]*u[1]);
923ab27746eSHong Zhang       uik[2] = -(diag[0]*u[2] + diag[2]*u[3]);
924ab27746eSHong Zhang       uik[3] = -(diag[1]*u[2] + diag[3]*u[3]);
925ab27746eSHong Zhang 
926ab27746eSHong Zhang       /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */
927ab27746eSHong Zhang       dk[0] += uik[0]*u[0] + uik[1]*u[1];
928ab27746eSHong Zhang       dk[1] += uik[2]*u[0] + uik[3]*u[1];
929ab27746eSHong Zhang       dk[2] += uik[0]*u[2] + uik[1]*u[3];
930ab27746eSHong Zhang       dk[3] += uik[2]*u[2] + uik[3]*u[3];
931ab27746eSHong Zhang 
932dc0b31edSSatish Balay       ierr = PetscLogFlops(16.0*2.0);CHKERRQ(ierr);
933187a9f4bSHong Zhang 
934ab27746eSHong Zhang       /* update -U(i,k): ba[ili] = uik */
935ab27746eSHong Zhang       ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr);
936ab27746eSHong Zhang 
937ab27746eSHong Zhang       /* add multiple of row i to k-th row ... */
938ab27746eSHong Zhang       jmin = ili + 1; jmax = bi[i+1];
939ab27746eSHong Zhang       if (jmin < jmax){
940ab27746eSHong Zhang         for (j=jmin; j<jmax; j++) {
941ab27746eSHong Zhang           /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */
942ab27746eSHong Zhang           rtmp_ptr = rtmp + bj[j]*4;
943ab27746eSHong Zhang           u = ba + j*4;
944ab27746eSHong Zhang           rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1];
945ab27746eSHong Zhang           rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1];
946ab27746eSHong Zhang           rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3];
947ab27746eSHong Zhang           rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3];
94849b5e25fSSatish Balay         }
949dc0b31edSSatish Balay         ierr = PetscLogFlops(16.0*(jmax-jmin));CHKERRQ(ierr);
950ab27746eSHong Zhang 
951ab27746eSHong Zhang         /* ... add i to row list for next nonzero entry */
952ab27746eSHong Zhang         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
953ab27746eSHong Zhang         j     = bj[jmin];
954ab27746eSHong Zhang         jl[i] = jl[j]; jl[j] = i; /* update jl */
95549b5e25fSSatish Balay       }
956ab27746eSHong Zhang       i = nexti;
95749b5e25fSSatish Balay     }
958ab27746eSHong Zhang 
959ab27746eSHong Zhang     /* save nonzero entries in k-th row of U ... */
960ab27746eSHong Zhang 
96149b5e25fSSatish Balay     /* invert diagonal block */
962ab27746eSHong Zhang     diag = ba+k*4;
963ab27746eSHong Zhang     ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr);
96462bba022SBarry Smith     ierr = Kernel_A_gets_inverse_A_2(diag,shift);CHKERRQ(ierr);
965ab27746eSHong Zhang 
966ab27746eSHong Zhang     jmin = bi[k]; jmax = bi[k+1];
967ab27746eSHong Zhang     if (jmin < jmax) {
968ab27746eSHong Zhang       for (j=jmin; j<jmax; j++){
969ab27746eSHong Zhang          vj = bj[j];           /* block col. index of U */
970ab27746eSHong Zhang          u   = ba + j*4;
971ab27746eSHong Zhang          rtmp_ptr = rtmp + vj*4;
972ab27746eSHong Zhang          for (k1=0; k1<4; k1++){
973ab27746eSHong Zhang            *u++        = *rtmp_ptr;
974ab27746eSHong Zhang            *rtmp_ptr++ = 0.0;
975ab27746eSHong Zhang          }
976ab27746eSHong Zhang       }
977ab27746eSHong Zhang 
978ab27746eSHong Zhang       /* ... add k to row list for first nonzero entry in k-th row */
979ab27746eSHong Zhang       il[k] = jmin;
980ab27746eSHong Zhang       i     = bj[jmin];
981ab27746eSHong Zhang       jl[k] = jl[i]; jl[i] = k;
982ab27746eSHong Zhang     }
98349b5e25fSSatish Balay   }
98449b5e25fSSatish Balay 
98549b5e25fSSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
986d8c74875SBarry Smith   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
987ab27746eSHong Zhang 
988db4efbfdSBarry Smith   C->ops->solve          = MatSolve_SeqSBAIJ_2_NaturalOrdering;
989db4efbfdSBarry Smith   C->ops->solvetranspose = MatSolve_SeqSBAIJ_2_NaturalOrdering;
990db4efbfdSBarry Smith   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_2_NaturalOrdering;
991db4efbfdSBarry Smith   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_2_NaturalOrdering;
99249b5e25fSSatish Balay   C->assembled = PETSC_TRUE;
9935c0bcdfcSHong Zhang   C->preallocated = PETSC_TRUE;
994efee365bSSatish Balay   ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
99549b5e25fSSatish Balay   PetscFunctionReturn(0);
99649b5e25fSSatish Balay }
99749b5e25fSSatish Balay 
99849b5e25fSSatish Balay /*
99998a8e78dSHong Zhang     Numeric U^T*D*U factorization for SBAIJ format.
100091602c66SHong Zhang     Version for blocks are 1 by 1.
100149b5e25fSSatish Balay */
10024a2ae208SSatish Balay #undef __FUNCT__
10034a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1"
10040481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat C,Mat A,const MatFactorInfo *info)
100549b5e25fSSatish Balay {
100649b5e25fSSatish Balay   Mat_SeqSBAIJ   *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ *)C->data;
100749b5e25fSSatish Balay   IS             ip=b->row;
10086849ba73SBarry Smith   PetscErrorCode ierr;
10095d0c19d7SBarry Smith   const PetscInt *ai,*aj,*rip;
10105d0c19d7SBarry Smith   PetscInt       *a2anew,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j,*bcol;
1011997a0949SHong Zhang   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
1012997a0949SHong Zhang   MatScalar      *rtmp,*ba=b->a,*bval,*aa,dk,uikdi;
10133cea4cbeSHong Zhang   PetscReal      zeropivot,rs,shiftnz;
1014fbf22428SSatish Balay   PetscReal      shiftpd;
10153cea4cbeSHong Zhang   ChShift_Ctx    sctx;
10163cea4cbeSHong Zhang   PetscInt       newshift;
101749b5e25fSSatish Balay 
101849b5e25fSSatish Balay   PetscFunctionBegin;
10193cea4cbeSHong Zhang   /* initialization */
10203cea4cbeSHong Zhang   shiftnz   = info->shiftnz;
10213cea4cbeSHong Zhang   shiftpd   = info->shiftpd;
10223cea4cbeSHong Zhang   zeropivot = info->zeropivot;
10233cea4cbeSHong Zhang 
102449b5e25fSSatish Balay   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
1025cb718733SHong Zhang   if (!a->permute){
10262d007305SHong Zhang     ai = a->i; aj = a->j; aa = a->a;
10272d007305SHong Zhang   } else {
10282d007305SHong Zhang     ai = a->inew; aj = a->jnew;
1029fff829cfSHong Zhang     nz = ai[mbs];
1030fff829cfSHong Zhang     ierr = PetscMalloc(nz*sizeof(MatScalar),&aa);CHKERRQ(ierr);
1031fff829cfSHong Zhang     a2anew = a->a2anew;
1032997a0949SHong Zhang     bval   = a->a;
1033fff829cfSHong Zhang     for (j=0; j<nz; j++){
1034997a0949SHong Zhang       aa[a2anew[j]] = *(bval++);
10352d007305SHong Zhang     }
10362d007305SHong Zhang   }
10372d007305SHong Zhang 
103891602c66SHong Zhang   /* initialization */
103949b5e25fSSatish Balay   /* il and jl record the first nonzero element in each row of the accessing
104049b5e25fSSatish Balay      window U(0:k, k:mbs-1).
104149b5e25fSSatish Balay      jl:    list of rows to be added to uneliminated rows
104249b5e25fSSatish Balay             i>= k: jl(i) is the first row to be added to row i
104349b5e25fSSatish Balay             i<  k: jl(i) is the row following row i in some list of rows
104449b5e25fSSatish Balay             jl(i) = mbs indicates the end of a list
104549b5e25fSSatish Balay      il(i): points to the first nonzero element in columns k,...,mbs-1 of
104649b5e25fSSatish Balay             row i of U */
1047d8c74875SBarry Smith   ierr = PetscMalloc3(mbs,MatScalar,&rtmp,mbs,PetscInt,&il,mbs,PetscInt,&jl);CHKERRQ(ierr);
1048997a0949SHong Zhang 
10493cea4cbeSHong Zhang   sctx.shift_amount = 0;
10503cea4cbeSHong Zhang   sctx.nshift       = 0;
1051997a0949SHong Zhang   do {
10523cea4cbeSHong Zhang     sctx.chshift = PETSC_FALSE;
105349b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
105449b5e25fSSatish Balay       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
105549b5e25fSSatish Balay     }
1056997a0949SHong Zhang 
105749b5e25fSSatish Balay     for (k = 0; k<mbs; k++){
1058997a0949SHong Zhang       /*initialize k-th row by the perm[k]-th row of A */
105949b5e25fSSatish Balay       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
10607701de02SHong Zhang       bval = ba + bi[k];
106149b5e25fSSatish Balay       for (j = jmin; j < jmax; j++){
1062997a0949SHong Zhang         col = rip[aj[j]];
1063997a0949SHong Zhang         rtmp[col] = aa[j];
10647701de02SHong Zhang         *bval++  = 0.0; /* for in-place factorization */
106549b5e25fSSatish Balay       }
10663cea4cbeSHong Zhang 
10673cea4cbeSHong Zhang       /* shift the diagonal of the matrix */
10683cea4cbeSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
106949b5e25fSSatish Balay 
107091602c66SHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
107149b5e25fSSatish Balay       dk = rtmp[k];
107249b5e25fSSatish Balay       i = jl[k]; /* first row to be added to k_th row  */
107349b5e25fSSatish Balay 
1074057f5ba7SHong Zhang       while (i < k){
107549b5e25fSSatish Balay         nexti = jl[i]; /* next row to be added to k_th row */
107649b5e25fSSatish Balay 
1077fff829cfSHong Zhang         /* compute multiplier, update diag(k) and U(i,k) */
107849b5e25fSSatish Balay         ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
1079997a0949SHong Zhang         uikdi = - ba[ili]*ba[bi[i]];  /* diagonal(k) */
108049b5e25fSSatish Balay         dk += uikdi*ba[ili];
1081658e7b3eSHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
108249b5e25fSSatish Balay 
1083997a0949SHong Zhang         /* add multiple of row i to k-th row */
108449b5e25fSSatish Balay         jmin = ili + 1; jmax = bi[i+1];
108549b5e25fSSatish Balay         if (jmin < jmax){
108649b5e25fSSatish Balay           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
1087dc0b31edSSatish Balay           ierr = PetscLogFlops(2.0*(jmax-jmin));CHKERRQ(ierr);
1088187a9f4bSHong Zhang 
1089fff829cfSHong Zhang           /* update il and jl for row i */
1090fff829cfSHong Zhang           il[i] = jmin;
1091fff829cfSHong Zhang           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
109249b5e25fSSatish Balay         }
1093ab27746eSHong Zhang         i = nexti;
109449b5e25fSSatish Balay       }
109549b5e25fSSatish Balay 
10963cea4cbeSHong Zhang       /* shift the diagonals when zero pivot is detected */
10973cea4cbeSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
10983cea4cbeSHong Zhang       rs   = 0.0;
1099997a0949SHong Zhang       jmin = bi[k]+1;
1100997a0949SHong Zhang       nz   = bi[k+1] - jmin;
1101997a0949SHong Zhang       if (nz){
1102997a0949SHong Zhang         bcol = bj + jmin;
1103997a0949SHong Zhang         while (nz--){
110489f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
110589f845c8SHong Zhang           bcol++;
1106997a0949SHong Zhang         }
1107997a0949SHong Zhang       }
11083cea4cbeSHong Zhang 
11093cea4cbeSHong Zhang       sctx.rs = rs;
11103cea4cbeSHong Zhang       sctx.pv = dk;
111145938b79SHong Zhang       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
111245938b79SHong Zhang       if (newshift == 1) break;    /* sctx.shift_amount is updated */
111349b5e25fSSatish Balay 
1114997a0949SHong Zhang       /* copy data into U(k,:) */
111598a8e78dSHong Zhang       ba[bi[k]] = 1.0/dk; /* U(k,k) */
1116fff829cfSHong Zhang       jmin = bi[k]+1; jmax = bi[k+1];
111749b5e25fSSatish Balay       if (jmin < jmax) {
111849b5e25fSSatish Balay         for (j=jmin; j<jmax; j++){
1119997a0949SHong Zhang           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
112049b5e25fSSatish Balay         }
112198a8e78dSHong Zhang         /* add the k-th row into il and jl */
112249b5e25fSSatish Balay         il[k] = jmin;
112398a8e78dSHong Zhang         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
112449b5e25fSSatish Balay       }
112549b5e25fSSatish Balay     }
11263cea4cbeSHong Zhang   } while (sctx.chshift);
1127d8c74875SBarry Smith   ierr = PetscFree3(rtmp,il,jl);CHKERRQ(ierr);
112898a8e78dSHong Zhang   if (a->permute){ierr = PetscFree(aa);CHKERRQ(ierr);}
112949b5e25fSSatish Balay 
113049b5e25fSSatish Balay   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
1131*0a3c351aSHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_1_inplace;
1132db4efbfdSBarry Smith   C->ops->solves         = MatSolves_SeqSBAIJ_1;
1133*0a3c351aSHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_1_inplace;
1134*0a3c351aSHong Zhang   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_1_inplace;
1135*0a3c351aSHong Zhang   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_1_inplace;
113649b5e25fSSatish Balay   C->assembled    = PETSC_TRUE;
11375c0bcdfcSHong Zhang   C->preallocated = PETSC_TRUE;
1138d0f46423SBarry Smith   ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr);
11393cea4cbeSHong Zhang   if (sctx.nshift){
11403cea4cbeSHong Zhang     if (shiftnz) {
11411e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
11423cea4cbeSHong Zhang     } else if (shiftpd) {
11431e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
1144997a0949SHong Zhang     }
1145997a0949SHong Zhang   }
114649b5e25fSSatish Balay   PetscFunctionReturn(0);
114749b5e25fSSatish Balay }
114849b5e25fSSatish Balay 
1149671cb588SHong Zhang /*
1150671cb588SHong Zhang   Version for when blocks are 1 by 1 Using natural ordering
1151671cb588SHong Zhang */
11524a2ae208SSatish Balay #undef __FUNCT__
11534a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering"
11540481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info)
1155671cb588SHong Zhang {
1156671cb588SHong Zhang   Mat_SeqSBAIJ   *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ *)C->data;
1157dfbe8321SBarry Smith   PetscErrorCode ierr;
115813f74950SBarry Smith   PetscInt       i,j,mbs = a->mbs;
115913f74950SBarry Smith   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
11603cea4cbeSHong Zhang   PetscInt       k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz;
1161653a6975SHong Zhang   MatScalar      *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval;
11623cea4cbeSHong Zhang   PetscReal      zeropivot,rs,shiftnz;
1163fbf22428SSatish Balay   PetscReal      shiftpd;
11643cea4cbeSHong Zhang   ChShift_Ctx    sctx;
11653cea4cbeSHong Zhang   PetscInt       newshift;
1166653a6975SHong Zhang 
1167653a6975SHong Zhang   PetscFunctionBegin;
1168653a6975SHong Zhang   /* initialization */
11693cea4cbeSHong Zhang   shiftnz   = info->shiftnz;
11703cea4cbeSHong Zhang   shiftpd   = info->shiftpd;
11713cea4cbeSHong Zhang   zeropivot = info->zeropivot;
11723cea4cbeSHong Zhang 
1173653a6975SHong Zhang   /* il and jl record the first nonzero element in each row of the accessing
1174653a6975SHong Zhang      window U(0:k, k:mbs-1).
1175653a6975SHong Zhang      jl:    list of rows to be added to uneliminated rows
1176653a6975SHong Zhang             i>= k: jl(i) is the first row to be added to row i
1177653a6975SHong Zhang             i<  k: jl(i) is the row following row i in some list of rows
1178653a6975SHong Zhang             jl(i) = mbs indicates the end of a list
1179653a6975SHong Zhang      il(i): points to the first nonzero element in U(i,k:mbs-1)
1180653a6975SHong Zhang   */
1181653a6975SHong Zhang   ierr = PetscMalloc(mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
1182d8c74875SBarry Smith   ierr = PetscMalloc2(mbs,PetscInt,&il,mbs,PetscInt,&jl);CHKERRQ(ierr);
1183b00f7748SHong Zhang 
11843cea4cbeSHong Zhang   sctx.shift_amount = 0;
11853cea4cbeSHong Zhang   sctx.nshift       = 0;
1186b00f7748SHong Zhang   do {
11873cea4cbeSHong Zhang     sctx.chshift = PETSC_FALSE;
1188653a6975SHong Zhang     for (i=0; i<mbs; i++) {
1189653a6975SHong Zhang       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
1190653a6975SHong Zhang     }
1191653a6975SHong Zhang 
1192997a0949SHong Zhang     for (k = 0; k<mbs; k++){
1193653a6975SHong Zhang       /*initialize k-th row with elements nonzero in row perm(k) of A */
1194653a6975SHong Zhang       nz   = ai[k+1] - ai[k];
1195653a6975SHong Zhang       acol = aj + ai[k];
1196653a6975SHong Zhang       aval = aa + ai[k];
1197653a6975SHong Zhang       bval = ba + bi[k];
1198653a6975SHong Zhang       while (nz -- ){
1199653a6975SHong Zhang         rtmp[*acol++] = *aval++;
1200653a6975SHong Zhang         *bval++       = 0.0; /* for in-place factorization */
1201653a6975SHong Zhang       }
12023cea4cbeSHong Zhang 
12033cea4cbeSHong Zhang       /* shift the diagonal of the matrix */
12043cea4cbeSHong Zhang       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
1205653a6975SHong Zhang 
1206653a6975SHong Zhang       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
1207653a6975SHong Zhang       dk = rtmp[k];
1208653a6975SHong Zhang       i  = jl[k]; /* first row to be added to k_th row  */
1209653a6975SHong Zhang 
1210653a6975SHong Zhang       while (i < k){
1211653a6975SHong Zhang         nexti = jl[i]; /* next row to be added to k_th row */
1212653a6975SHong Zhang         /* compute multiplier, update D(k) and U(i,k) */
1213653a6975SHong Zhang         ili   = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
1214653a6975SHong Zhang         uikdi = - ba[ili]*ba[bi[i]];
1215653a6975SHong Zhang         dk   += uikdi*ba[ili];
1216653a6975SHong Zhang         ba[ili] = uikdi; /* -U(i,k) */
1217653a6975SHong Zhang 
1218653a6975SHong Zhang         /* add multiple of row i to k-th row ... */
1219653a6975SHong Zhang         jmin = ili + 1;
1220653a6975SHong Zhang         nz   = bi[i+1] - jmin;
1221653a6975SHong Zhang         if (nz > 0){
1222653a6975SHong Zhang           bcol = bj + jmin;
1223653a6975SHong Zhang           bval = ba + jmin;
1224dc0b31edSSatish Balay           ierr = PetscLogFlops(2.0*nz);CHKERRQ(ierr);
1225653a6975SHong Zhang           while (nz --) rtmp[*bcol++] += uikdi*(*bval++);
1226187a9f4bSHong Zhang 
1227997a0949SHong Zhang           /* update il and jl for i-th row */
1228997a0949SHong Zhang           il[i] = jmin;
1229997a0949SHong Zhang           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
1230653a6975SHong Zhang         }
1231653a6975SHong Zhang         i = nexti;
1232653a6975SHong Zhang       }
1233653a6975SHong Zhang 
12343cea4cbeSHong Zhang       /* shift the diagonals when zero pivot is detected */
12353cea4cbeSHong Zhang       /* compute rs=sum of abs(off-diagonal) */
12363cea4cbeSHong Zhang       rs   = 0.0;
123721c26570Svictorle       jmin = bi[k]+1;
123821c26570Svictorle       nz   = bi[k+1] - jmin;
123921c26570Svictorle       if (nz){
124021c26570Svictorle         bcol = bj + jmin;
124121c26570Svictorle         while (nz--){
124289f845c8SHong Zhang           rs += PetscAbsScalar(rtmp[*bcol]);
124389f845c8SHong Zhang           bcol++;
124421c26570Svictorle         }
124521c26570Svictorle       }
12463cea4cbeSHong Zhang 
12473cea4cbeSHong Zhang       sctx.rs = rs;
12483cea4cbeSHong Zhang       sctx.pv = dk;
124945938b79SHong Zhang       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
125045938b79SHong Zhang       if (newshift == 1) break;    /* sctx.shift_amount is updated */
1251653a6975SHong Zhang 
1252997a0949SHong Zhang       /* copy data into U(k,:) */
1253653a6975SHong Zhang       ba[bi[k]] = 1.0/dk;
1254653a6975SHong Zhang       jmin      = bi[k]+1;
1255653a6975SHong Zhang       nz        = bi[k+1] - jmin;
1256653a6975SHong Zhang       if (nz){
1257653a6975SHong Zhang         bcol = bj + jmin;
1258653a6975SHong Zhang         bval = ba + jmin;
1259653a6975SHong Zhang         while (nz--){
1260653a6975SHong Zhang           *bval++       = rtmp[*bcol];
1261653a6975SHong Zhang           rtmp[*bcol++] = 0.0;
1262653a6975SHong Zhang         }
1263997a0949SHong Zhang         /* add k-th row into il and jl */
1264653a6975SHong Zhang         il[k] = jmin;
1265997a0949SHong Zhang         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
1266653a6975SHong Zhang       }
1267b00f7748SHong Zhang     } /* end of for (k = 0; k<mbs; k++) */
12683cea4cbeSHong Zhang   } while (sctx.chshift);
1269653a6975SHong Zhang   ierr = PetscFree(rtmp);CHKERRQ(ierr);
1270d8c74875SBarry Smith   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
1271653a6975SHong Zhang 
1272*0a3c351aSHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
1273db4efbfdSBarry Smith   C->ops->solves         = MatSolves_SeqSBAIJ_1;
1274*0a3c351aSHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
1275*0a3c351aSHong Zhang   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
1276*0a3c351aSHong Zhang   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
1277db4efbfdSBarry Smith 
1278653a6975SHong Zhang   C->assembled    = PETSC_TRUE;
1279653a6975SHong Zhang   C->preallocated = PETSC_TRUE;
1280d0f46423SBarry Smith   ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr);
12813cea4cbeSHong Zhang   if (sctx.nshift){
12823cea4cbeSHong Zhang     if (shiftnz) {
12831e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
12843cea4cbeSHong Zhang     } else if (shiftpd) {
12851e2582c4SBarry Smith       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
1286b00f7748SHong Zhang     }
1287fb10cecfSBarry Smith   }
1288653a6975SHong Zhang   PetscFunctionReturn(0);
1289653a6975SHong Zhang }
1290653a6975SHong Zhang 
1291653a6975SHong Zhang #undef __FUNCT__
12924a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactor_SeqSBAIJ"
12930481f469SBarry Smith PetscErrorCode MatCholeskyFactor_SeqSBAIJ(Mat A,IS perm,const MatFactorInfo *info)
129449b5e25fSSatish Balay {
1295dfbe8321SBarry Smith   PetscErrorCode ierr;
129649b5e25fSSatish Balay   Mat            C;
129749b5e25fSSatish Balay 
129849b5e25fSSatish Balay   PetscFunctionBegin;
1299719d5645SBarry Smith   ierr = MatGetFactor(A,"petsc",MAT_FACTOR_CHOLESKY,&C);CHKERRQ(ierr);
1300719d5645SBarry Smith   ierr = MatCholeskyFactorSymbolic(C,A,perm,info);CHKERRQ(ierr);
1301719d5645SBarry Smith   ierr = MatCholeskyFactorNumeric(C,A,info);CHKERRQ(ierr);
1302db4efbfdSBarry Smith   A->ops->solve            = C->ops->solve;
1303db4efbfdSBarry Smith   A->ops->solvetranspose   = C->ops->solvetranspose;
13044445ddedSHong Zhang   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
130549b5e25fSSatish Balay   PetscFunctionReturn(0);
130649b5e25fSSatish Balay }
130749b5e25fSSatish Balay 
130849b5e25fSSatish Balay 
1309